dscld Subroutine

public pure subroutine dscld(n, m, x, ldx, tt, ldtt)

Uses

  • proc~~dscld~~UsesGraph proc~dscld dscld module~odrpack_kinds odrpack_kinds proc~dscld->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack_kinds->iso_fortran_env

Select scaling values for delta according to the algorithm given in the ODRPACK95 reference guide.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

The number of observations.

integer, intent(in) :: m

The number of columns of data in the independent variable.

real(kind=wp), intent(in) :: x(ldx,m)

The independent variable.

integer, intent(in) :: ldx

The leading dimension of array x.

real(kind=wp), intent(out) :: tt(ldtt,m)

The scaling values for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.


Called by

proc~~dscld~~CalledByGraph proc~dscld dscld proc~diniwk diniwk proc~diniwk->proc~dscld proc~doddrv doddrv proc~doddrv->proc~diniwk proc~dodcnt dodcnt proc~dodcnt->proc~doddrv proc~odr odr proc~odr->proc~dodcnt proc~odr_long_c odr_long_c proc~odr_long_c->proc~odr proc~odr_medium_c odr_medium_c proc~odr_medium_c->proc~odr proc~odr_short_c odr_short_c proc~odr_short_c->proc~odr program~example1 example1 program~example1->proc~odr program~example2 example2 program~example2->proc~odr program~example3 example3 program~example3->proc~odr program~example4 example4 program~example4->proc~odr program~example5 example5 program~example5->proc~odr

Variables

Type Visibility Attributes Name Initial
real(kind=wp), public :: xmax
real(kind=wp), public :: xmin
integer, public :: i
integer, public :: j
logical, public :: bigdif

Source Code

   pure subroutine dscld(n, m, x, ldx, tt, ldtt)
   !! Select scaling values for `delta` according to the algorithm given in the ODRPACK95
   !! reference guide.
      ! Routines Called (NONE)
      ! Date Written   860529   (YYMMDD)
      ! Revision Date  920304   (YYMMDD)

      use odrpack_kinds, only: zero, one, ten

      integer, intent(in) :: n
         !! The number of observations.
      integer, intent(in) :: m
         !! The number of columns of data in the independent variable.
      real(wp), intent(in) :: x(ldx, m)
         !! The independent variable.
      integer, intent(in) :: ldx
         !! The leading dimension of array `x`.
      real(wp), intent(out) :: tt(ldtt, m)
         !! The scaling values for `delta`.
      integer, intent(in) :: ldtt
         !! The leading dimension of array `tt`.

      ! Local scalars
      real(wp) :: xmax, xmin
      integer :: i, j
      logical :: bigdif

      ! Variable Definitions (alphabetically)
      !  BIGDIF:  The variable designating whether there is a significant
      !           difference in the magnitudes of the nonzero elements of
      !           X (BIGDIF=.TRUE.) or not (BIGDIF=.FALSE.).
      !  I:       An indexing variable.
      !  J:       An indexing variable.
      !  LDTT:    The leading dimension of array TT.
      !  LDX:     The leading dimension of array X.
      !  M:       The number of columns of data in the independent variable.
      !  N:       The number of observations.
      !  TT:      THE SCALING VALUES FOR DELTA.
      !  X:       The independent variable.
      !  XMAX:    The largest nonzero magnitude.
      !  XMIN:    THE SMALLEST NONZERO MAGNITUDE.

      do j = 1, m
         xmax = abs(x(1, j))
         do i = 2, n
            xmax = max(xmax, abs(x(i, j)))
         end do
         if (xmax == zero) then
            !  All input values of X(I,J), I=1,...,N, are zero
            do i = 1, n
               tt(i, j) = one
            end do
         else
            !  Some of the input values are nonzero
            xmin = xmax
            do i = 1, n
               if (x(i, j) /= zero) then
                  xmin = min(xmin, abs(x(i, j)))
               end if
            end do
            bigdif = log10(xmax) - log10(xmin) >= one
            do i = 1, n
               if (x(i, j) /= zero) then
                  if (bigdif) then
                     tt(i, j) = one/abs(x(i, j))
                  else
                     tt(i, j) = one/xmax
                  end if
               else
                  tt(i, j) = ten/xmin
               end if
            end do
         end if
      end do

   end subroutine dscld