Select scaling values for delta
according to the algorithm given in the ODRPACK95
reference guide.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
Number of observations. |
||
integer, | intent(in) | :: | m |
Number of columns of data in the independent variable. |
||
real(kind=wp), | intent(in) | :: | x(n,m) |
Independent variable. |
||
real(kind=wp), | intent(out) | :: | tt(ldtt,m) |
Scaling values for |
||
integer, | intent(in) | :: | ldtt |
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 scale_delta(n, m, x, tt, ldtt) !! Select scaling values for `delta` according to the algorithm given in the ODRPACK95 !! reference guide. use odrpack_kinds, only: zero, one integer, intent(in) :: n !! Number of observations. integer, intent(in) :: m !! Number of columns of data in the independent variable. real(wp), intent(in) :: x(n, m) !! Independent variable. real(wp), intent(out) :: tt(ldtt, m) !! Scaling values for `delta`. integer, intent(in) :: ldtt !! 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 (TRUE) or not (FALSE). ! I: An indexing variable. ! J: An indexing 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 tt(1:n, j) = one 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) = 1/abs(x(i, j)) else tt(i, j) = 1/xmax end if else tt(i, j) = 10/xmin end if end do end if end do end subroutine scale_delta