QFORM Subroutine

public subroutine QFORM(m, n, q, ldq, wa)

Arguments

Type IntentOptional AttributesName
integer :: m
integer :: n
real(kind=dp) :: q(ldq,m)
integer :: ldq
real(kind=dp) :: wa(m)

Contents

Source Code


Source Code

  SUBROUTINE QFORM(m,n,q,ldq,wa)

    !  This subroutine proceeds from the computed QR factorization of
    !  an M by N matrix A to accumulate the M by M orthogonal matrix
    !  Q from its factored form.
    !
    !  The subroutine statement is
    !
    !  subroutine qform(m,n,q,ldq,wa)
    !
    !  where
    !
    !  M is a positive integer input variable set to the number
    !  of rows of A and the order of Q.
    !
    !  N is a positive integer input variable set to the number
    !  of columns of A.
    !
    !  Q is an M by M array. On input the full lower trapezoid in
    !  the first min(M,N) columns of Q contains the factored form.
    !  On output Q has been accumulated into a square matrix.
    !
    !  LDQ is a positive integer input variable not less than M
    !  which specifies the leading dimension of the array Q.
    !
    !  WA is a work array of length M.

    !  Argonne National Laboratory. Minpack project. March 1980.
    !  Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More

    IMPLICIT NONE

    INTEGER m,n,ldq,i,j,jm1,k,l,minmn,np1

    real(dp) q(ldq,m),wa(m)
    real(dp) one,sum,temp,zero

    one = 1.0D0
    zero = 0.0D0

    !  Zero out upper triangle of q in the first min(m,n) columns.

    minmn = min(m,n)
    if (minmn < 2) goto 30
    do j = 2, minmn
       jm1 = j - 1
       do i = 1, jm1
          q(i,j) = zero
       end do
    end do

30  continue

    !  Initialize remaining columns to those of the identity matrix.

    np1 = n + 1
    if (m < np1) goto 60
    do j = np1, m
       do i = 1, m
          q(i,j) = zero
       end do
       q(j,j) = one
    end do

60  continue

    !  Accumulate q from its factored form.

    do l = 1, minmn
       k = minmn - l + 1
       do i = k, m
          wa(i) = q(i,k)
          q(i,k) = zero
       end do
       q(k,k) = one
       if (wa(k) == zero) goto 110
       do j = k, m
          sum = zero
          do i = k, m
             sum = sum + q(i,j)*wa(i)
          end do
          temp = sum/wa(k)
          do i = k, m
             q(i,j) = q(i,j) - temp*wa(i)
          end do
       end do

110    continue
    end do

    return
  end SUBROUTINE QFORM