Processing math: 100%

g1 Subroutine

private subroutine g1(a, b, c, s, sig)

Compute orthogonal rotation matrix.

Compute matrix [cssc] so that

[cssc][ab]=[a2+b20]

Compute σ=a2+b2

σ is computed last to allow for the possibility that sig may be in the same location as a or b.

References

  • C.L. Lawson, R.J. Hanson, 'Solving least squares problems' Prentice Hall, 1974. (revised 1995 edition)
  • lawson-hanson from Netlib.

History

  • Jacob Williams, refactored into modern Fortran, Jan. 2016.

Arguments

Type IntentOptional Attributes Name
real(kind=wp) :: a
real(kind=wp) :: b
real(kind=wp), intent(out) :: c
real(kind=wp), intent(out) :: s
real(kind=wp) :: sig

Called by

proc~~g1~~CalledByGraph proc~g1 g1 proc~nnls nnls proc~nnls->proc~g1 proc~ldp ldp proc~ldp->proc~nnls proc~lsi lsi proc~lsi->proc~ldp proc~lsei lsei proc~lsei->proc~lsi proc~lsq lsq proc~lsq->proc~lsei proc~slsqpb slsqpb proc~slsqpb->proc~lsq proc~slsqp slsqp proc~slsqp->proc~slsqpb proc~slsqp_wrapper slsqp_solver%slsqp_wrapper proc~slsqp_wrapper->proc~slsqp

Source Code

    subroutine g1(a,b,c,s,sig)

    implicit none

    real(wp)             :: a
    real(wp)             :: b
    real(wp)             :: sig
    real(wp),intent(out) :: c
    real(wp),intent(out) :: s

    real(wp) :: xr, yr

    if ( abs(a)>abs(b) ) then
        xr = b/a
        yr = sqrt(one+xr**2)
        c = sign(one/yr,a)
        s = c*xr
        sig = abs(a)*yr
    else
        if ( abs(b)>zero ) then
            xr = a/b
            yr = sqrt(one+xr**2)
            s = sign(one/yr,b)
            c = s*xr
            sig = abs(b)*yr
        else
            sig = zero
            c = zero
            s = one
        end if
    end if

    end subroutine g1