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.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer | :: | m | ||||
integer | :: | n | ||||
real(kind=wp) | :: | q(ldq,m) | ||||
integer | :: | ldq | ||||
real(kind=wp) | :: | wa(m) |
subroutine qform(m,n,q,ldq,wa) implicit none integer m , n , ldq real(wp) q(ldq,m) , wa(m) integer i , j , jm1 , k , l , minmn , np1 real(wp) sum , temp ! ! zero out upper triangle of q in the first min(m,n) columns. ! minmn = min(m,n) if ( minmn>=2 ) then do j = 2 , minmn jm1 = j - 1 do i = 1 , jm1 q(i,j) = zero enddo enddo endif ! ! initialize remaining columns to those of the identity matrix. ! np1 = n + 1 if ( m>=np1 ) then do j = np1 , m do i = 1 , m q(i,j) = zero enddo q(j,j) = one enddo endif ! ! 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 enddo q(k,k) = one if ( wa(k)/=zero ) then do j = k , m sum = zero do i = k , m sum = sum + q(i,j)*wa(i) enddo temp = sum/wa(k) do i = k , m q(i,j) = q(i,j) - temp*wa(i) enddo enddo endif enddo end subroutine qform