F08ZFF (DGGRQF) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F08ZFF (DGGRQF)

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

F08ZFF (DGGRQF) computes a generalized RQ factorization of a real matrix pair A,B, where A is an m by n matrix and B is a p by n matrix.

2  Specification

SUBROUTINE F08ZFF ( M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
INTEGER  M, P, N, LDA, LDB, LWORK, INFO
REAL (KIND=nag_wp)  A(LDA,*), TAUA(min(M,N)), B(LDB,*), TAUB(min(P,N)), WORK(max(1,LWORK))
The routine may be called by its LAPACK name dggrqf.

3  Description

F08ZFF (DGGRQF) forms the generalized RQ factorization of an m by n matrix A and a p by n matrix B 
A = RQ ,   B= ZTQ ,
where Q is an n by n orthogonal matrix, Z is a p by p orthogonal matrix and R and T are of the form
R = n-mmm(0R12) ;   if ​ mn , nm-n(R11) n R21 ;   if ​ m>n ,
with R12 or R21 upper triangular,
T = nn(T11) p-n 0 ;   if ​ pn , pn-pp(T11T12) ;   if ​ p<n ,
with T11 upper triangular.
In particular, if B is square and nonsingular, the generalized RQ factorization of A and B implicitly gives the RQ factorization of AB-1 as
AB-1= R T-1 ZT .

4  References

Anderson E, Bai Z, Bischof C, Blackford S, Demmel J, Dongarra J J, Du Croz J J, Greenbaum A, Hammarling S, McKenney A and Sorensen D (1999) LAPACK Users' Guide (3rd Edition) SIAM, Philadelphia http://www.netlib.org/lapack/lug
Anderson E, Bai Z and Dongarra J (1992) Generalized QR factorization and its applications Linear Algebra Appl. (Volume 162–164) 243–271
Hammarling S (1987) The numerical solution of the general Gauss-Markov linear model Mathematics in Signal Processing (eds T S Durrani, J B Abbiss, J E Hudson, R N Madan, J G McWhirter and T A Moore) 441–456 Oxford University Press
Paige C C (1990) Some aspects of generalized QR factorizations . In Reliable Numerical Computation (eds M G Cox and S Hammarling) 73–91 Oxford University Press

5  Parameters

1:     M – INTEGERInput
On entry: m, the number of rows of the matrix A.
Constraint: M0.
2:     P – INTEGERInput
On entry: p, the number of rows of the matrix B.
Constraint: P0.
3:     N – INTEGERInput
On entry: n, the number of columns of the matrices A and B.
Constraint: N0.
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: the m by n matrix A.
On exit: if mn, the upper triangle of the subarray A1:mn-m+1:n contains the m by m upper triangular matrix R12.
If mn, the elements on and above the m-nth subdiagonal contain the m by n upper trapezoidal matrix R; the remaining elements, with the array TAUA, represent the orthogonal matrix Q as a product of minm,n elementary reflectors (see Section 3.3.6 in the F08 Chapter Introduction).
5:     LDA – INTEGERInput
On entry: the first dimension of the array A as declared in the (sub)program from which F08ZFF (DGGRQF) is called.
Constraint: LDAmax1,M.
6:     TAUA(minM,N) – REAL (KIND=nag_wp) arrayOutput
On exit: the scalar factors of the elementary reflectors which represent the orthogonal matrix Q.
7:     B(LDB,*) – REAL (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array B must be at least max1,N.
On entry: the p by n matrix B.
On exit: the elements on and above the diagonal of the array contain the minp,n by n upper trapezoidal matrix T (T is upper triangular if pn); the elements below the diagonal, with the array TAUB, represent the orthogonal matrix Z as a product of elementary reflectors (see Section 3.3.6 in the F08 Chapter Introduction).
8:     LDB – INTEGERInput
On entry: the first dimension of the array B as declared in the (sub)program from which F08ZFF (DGGRQF) is called.
Constraint: LDBmax1,P.
9:     TAUB(minP,N) – REAL (KIND=nag_wp) arrayOutput
On exit: the scalar factors of the elementary reflectors which represent the orthogonal matrix Z.
10:   WORK(max1,LWORK) – REAL (KIND=nag_wp) arrayWorkspace
On exit: if INFO=0, WORK1 contains the minimum value of LWORK required for optimal performance.
11:   LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F08ZFF (DGGRQF) 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, LWORKmaxN,M,P×maxnb1,nb2,nb3, where nb1 is the optimal block size for the RQ factorization of an m by n matrix by F08CHF (DGERQF), nb2 is the optimal block size for the QR factorization of a p by n matrix by F08AEF (DGEQRF), and nb3 is the optimal block size for a call of F08CKF (DORMRQ).
Constraint: LWORKmax1,N,M,P or LWORK=-1.
12:   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 generalized RQ factorization is the exact factorization for nearby matrices A+E and B+F, where
E2 = Oε A2   and   F2= Oε B2 ,
and ε is the machine precision.

8  Further Comments

The orthogonal matrices Q and Z may be formed explicitly by calls to F08CJF (DORGRQ) and F08AFF (DORGQR) respectively. F08CKF (DORMRQ) may be used to multiply Q by another matrix and F08AGF (DORMQR) may be used to multiply Z by another matrix.
The complex analogue of this routine is F08ZTF (ZGGRQF).

9  Example

This example solves the linear equality constrained least squares problem
minx c-Ax2   subject to   Bx= d
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 ,   B= 1 0 -1 0 0 1 0 -1 ,
c = -1.50 -2.14 1.23 -0.54 -1.68 0.82   and   d= 0 0 .
The constraints Bx=d correspond to x1=x3 and x2=x4.
The solution is obtained by first computing a generalized RQ factorization of the matrix pair B,A. The example illustrates the general solution process.
Note that the block size (NB) of 64 assumed in this example is not realistic for such a small problem, but should be suitable for large problems.

9.1  Program Text

Program Text (f08zffe.f90)

9.2  Program Data

Program Data (f08zffe.d)

9.3  Program Results

Program Results (f08zffe.r)


F08ZFF (DGGRQF) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

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