Select scaling values for delta
according to the algorithm given in the ODRPACK95
reference guide.
Type | Intent | Optional | 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 |
||
real(kind=wp), | intent(out) | :: | tt(ldtt,m) |
The scaling values for |
||
integer, | intent(in) | :: | ldtt |
The leading dimension of array |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=wp), | public | :: | xmax | ||||
real(kind=wp), | public | :: | xmin | ||||
integer, | public | :: | i | ||||
integer, | public | :: | j | ||||
logical, | public | :: | bigdif |
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