enorm Function

public pure function enorm(n, x)

given an n-vector x, this function calculates the euclidean norm of x.

the euclidean norm is computed by accumulating the sum of squares in three different sums. the sums of squares for the small and large components are scaled so that no overflows occur. non-destructive underflows are permitted. underflows and overflows do not occur in the computation of the unscaled sum of squares for the intermediate components. the definitions of small, intermediate and large components depend on two constants, rdwarf and rgiant. the main restrictions on these constants are that rdwarf2 not underflow and rgiant2 not overflow. the constants given here are suitable for every known computer.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable.

real(kind=wp), intent(in) :: x(n)

an input array of length n.

Return Value real(kind=wp)


Called by

proc~~enorm~~CalledByGraph proc~enorm enorm proc~dogleg dogleg proc~dogleg->proc~enorm proc~lmstr lmstr proc~lmstr->proc~enorm proc~qrfac qrfac proc~lmstr->proc~qrfac proc~lmpar lmpar proc~lmstr->proc~lmpar proc~hybrj hybrj proc~hybrj->proc~enorm proc~hybrj->proc~dogleg proc~hybrj->proc~qrfac proc~qrfac->proc~enorm proc~lmpar->proc~enorm proc~lmdif lmdif proc~lmdif->proc~enorm proc~lmdif->proc~qrfac proc~lmdif->proc~lmpar proc~hybrd hybrd proc~hybrd->proc~enorm proc~hybrd->proc~dogleg proc~hybrd->proc~qrfac proc~lmder lmder proc~lmder->proc~enorm proc~lmder->proc~qrfac proc~lmder->proc~lmpar proc~minpack_lmdif minpack_lmdif proc~minpack_lmdif->proc~lmdif proc~minpack_hybrd minpack_hybrd proc~minpack_hybrd->proc~hybrd proc~lmstr1 lmstr1 proc~lmstr1->proc~lmstr proc~minpack_lmstr minpack_lmstr proc~minpack_lmstr->proc~lmstr proc~hybrj1 hybrj1 proc~hybrj1->proc~hybrj proc~minpack_lmder minpack_lmder proc~minpack_lmder->proc~lmder proc~lmdif1 lmdif1 proc~lmdif1->proc~lmdif proc~lmder1 lmder1 proc~lmder1->proc~lmder proc~hybrd1 hybrd1 proc~hybrd1->proc~hybrd proc~minpack_hybrj minpack_hybrj proc~minpack_hybrj->proc~hybrj proc~minpack_hybrj1 minpack_hybrj1 proc~minpack_hybrj1->proc~hybrj1 proc~minpack_lmstr1 minpack_lmstr1 proc~minpack_lmstr1->proc~lmstr1 proc~minpack_lmdif1 minpack_lmdif1 proc~minpack_lmdif1->proc~lmdif1 proc~minpack_lmder1 minpack_lmder1 proc~minpack_lmder1->proc~lmder1 proc~minpack_hybrd1 minpack_hybrd1 proc~minpack_hybrd1->proc~hybrd1

Contents

Source Code


Source Code

    pure real(wp) function enorm(n, x)

        implicit none

        integer, intent(in) :: n !! a positive integer input variable.
        real(wp), intent(in) :: x(n) !! an input array of length n.

        integer :: i
        real(wp) :: agiant, s1, s2, s3, xabs, x1max, x3max

        real(wp), parameter :: rdwarf = 3.834e-20_wp
        real(wp), parameter :: rgiant = 1.304e19_wp

        s1 = zero
        s2 = zero
        s3 = zero
        x1max = zero
        x3max = zero
        agiant = rgiant/real(n, wp)
        do i = 1, n
            xabs = abs(x(i))
            if (xabs > rdwarf .and. xabs < agiant) then
                ! sum for intermediate components.
                s2 = s2 + xabs**2
            elseif (xabs <= rdwarf) then
                ! sum for small components.
                if (xabs <= x3max) then
                    if (xabs /= zero) s3 = s3 + (xabs/x3max)**2
                else
                    s3 = one + s3*(x3max/xabs)**2
                    x3max = xabs
                end if
                ! sum for large components.
            elseif (xabs <= x1max) then
                s1 = s1 + (xabs/x1max)**2
            else
                s1 = one + s1*(x1max/xabs)**2
                x1max = xabs
            end if
        end do

        ! calculation of norm.

        if (s1 /= zero) then
            enorm = x1max*sqrt(s1 + (s2/x1max)/x1max)
        elseif (s2 == zero) then
            enorm = x3max*sqrt(s3)
        else
            if (s2 >= x3max) enorm = sqrt(s2*(one + (x3max/s2)*(x3max*s3)))
            if (s2 < x3max) enorm = sqrt(x3max*((s2/x3max) + (x3max*s3)))
        end if

    end function enorm