scale_beta Subroutine

public pure subroutine scale_beta(np, beta, ssf)

Uses

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

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

Arguments

Type IntentOptional 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 beta.


Called by

proc~~scale_beta~~CalledByGraph proc~scale_beta scale_beta proc~init_work init_work proc~init_work->proc~scale_beta proc~odr odr proc~odr->proc~init_work 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 :: bmax
real(kind=wp), public :: bmin
integer, public :: k
logical, public :: bigdif

Source Code

   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