F08AFF (DORGQR) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F08AFF (DORGQR)

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

F08AFF (DORGQR) generates all or part of the real orthogonal matrix Q from a QR factorization computed by F08AEF (DGEQRF), F08BEF (DGEQPF) or F08BFF (DGEQP3).

2  Specification

SUBROUTINE F08AFF ( M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
INTEGER  M, N, K, LDA, LWORK, INFO
REAL (KIND=nag_wp)  A(LDA,*), TAU(*), WORK(max(1,LWORK))
The routine may be called by its LAPACK name dorgqr.

3  Description

F08AFF (DORGQR) is intended to be used after a call to F08AEF (DGEQRF), F08BEF (DGEQPF) or F08BFF (DGEQP3). which perform a QR factorization of a real matrix A. The orthogonal matrix Q is represented as a product of elementary reflectors.
This routine may be used to generate Q explicitly as a square matrix, or to form only its leading columns.
Usually Q is determined from the QR factorization of an m by p matrix A with mp. The whole of Q may be computed by:
CALL DORGQR(M,M,P,A,LDA,TAU,WORK,LWORK,INFO)
(note that the array A must have at least m columns) or its leading p columns by:
CALL DORGQR(M,P,P,A,LDA,TAU,WORK,LWORK,INFO)
The columns of Q returned by the last call form an orthonormal basis for the space spanned by the columns of A; thus F08AEF (DGEQRF) followed by F08AFF (DORGQR) can be used to orthogonalize the columns of A.
The information returned by the QR factorization routines also yields the QR factorization of the leading k columns of A, where k<p. The orthogonal matrix arising from this factorization can be computed by:
CALL DORGQR(M,M,K,A,LDA,TAU,WORK,LWORK,INFO)
or its leading k columns by:
CALL DORGQR(M,K,K,A,LDA,TAU,WORK,LWORK,INFO)

4  References

Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

5  Parameters

1:     M – INTEGERInput
On entry: m, the order of the orthogonal matrix Q.
Constraint: M0.
2:     N – INTEGERInput
On entry: n, the number of columns of the matrix Q.
Constraint: MN0.
3:     K – INTEGERInput
On entry: k, the number of elementary reflectors whose product defines the matrix Q.
Constraint: NK0.
4:     A(LDA,*) – REAL (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array A must be at least max1,N.
On entry: details of the vectors which define the elementary reflectors, as returned by F08AEF (DGEQRF), F08BEF (DGEQPF) or F08BFF (DGEQP3).
On exit: the m by n matrix Q.
5:     LDA – INTEGERInput
On entry: the first dimension of the array A as declared in the (sub)program from which F08AFF (DORGQR) is called.
Constraint: LDAmax1,M.
6:     TAU(*) – REAL (KIND=nag_wp) arrayInput
Note: the dimension of the array TAU must be at least max1,K.
On entry: further details of the elementary reflectors, as returned by F08AEF (DGEQRF), F08BEF (DGEQPF) or F08BFF (DGEQP3).
7:     WORK(max1,LWORK) – REAL (KIND=nag_wp) arrayWorkspace
On exit: if INFO=0, WORK1 contains the minimum value of LWORK required for optimal performance.
8:     LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F08AFF (DORGQR) is called.
If LWORK=-1, a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued.
Suggested value: for optimal performance, LWORKN×nb, where nb is the optimal block size.
Constraint: LWORKmax1,N or LWORK=-1.
9:     INFO – INTEGEROutput
On exit: INFO=0 unless the routine detects an error (see Section 6).

6  Error Indicators and Warnings

Errors or warnings detected by the routine:
INFO<0
If INFO=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.

7  Accuracy

The computed matrix Q differs from an exactly orthogonal matrix by a matrix E such that
E2 = Oε ,
where ε is the machine precision.

8  Further Comments

The total number of floating point operations is approximately 4mnk-2 m+n k2 + 43 k3 ; when n=k, the number is approximately 23 n2 3m-n .
The complex analogue of this routine is F08ATF (ZUNGQR).

9  Example

This example forms the leading 4 columns of the orthogonal matrix Q from the QR factorization of the matrix A, where
A = -0.57 -1.28 -0.39 0.25 -1.93 1.08 -0.31 -2.14 2.30 0.24 0.40 -0.35 -1.93 0.64 -0.66 0.08 0.15 0.30 0.15 -2.13 -0.02 1.03 -1.43 0.50 .
The columns of Q form an orthonormal basis for the space spanned by the columns of A.

9.1  Program Text

Program Text (f08affe.f90)

9.2  Program Data

Program Data (f08affe.d)

9.3  Program Results

Program Results (f08affe.r)


F08AFF (DORGQR) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012