Select scaling values for beta
according to the algorithm given in the ODRPACK95
reference guide.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | np |
Number of function parameters. |
||
real(kind=wp), | intent(in) | :: | beta(np) |
Function parameters. |
||
real(kind=wp), | intent(out) | :: | ssf(np) |
Scaling values for |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=wp), | public | :: | bmax | ||||
real(kind=wp), | public | :: | bmin | ||||
integer, | public | :: | k | ||||
logical, | public | :: | bigdif |
pure subroutine scale_beta(np, beta, ssf) !! Select scaling values for `beta` according to the algorithm given in the ODRPACK95 !! reference guide. use odrpack_kinds, only: zero, one integer, intent(in) :: np !! Number of function parameters. real(wp), intent(in) :: beta(np) !! Function parameters. real(wp), intent(out) :: ssf(np) !! Scaling values for `beta`. ! Local scalars real(wp) :: bmax, bmin integer :: k logical ::bigdif ! Variable Definitions (alphabetically) ! BIGDIF: The variable designating whether there is a significant difference in the ! magnitudes of the nonzero elements of BETA (TRUE) or not (FALSE). ! BMAX: The largest nonzero magnitude. ! BMIN: The smallest nonzero magnitude. ! K: An indexing variable. bmax = abs(beta(1)) do k = 2, np bmax = max(bmax, abs(beta(k))) end do if (bmax == zero) then ! All input values of BETA are zero ssf(1:np) = one else ! Some of the input values are nonzero bmin = bmax do k = 1, np if (beta(k) /= zero) then bmin = min(bmin, abs(beta(k))) end if end do bigdif = log10(bmax) - log10(bmin) >= one do k = 1, np if (beta(k) == zero) then ssf(k) = 10/bmin else if (bigdif) then ssf(k) = 1/abs(beta(k)) else ssf(k) = 1/bmax end if end if end do end if end subroutine scale_beta