nag_zgeqrt (f08apc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_zgeqrt (f08apc)

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_zgeqrt (f08apc) recursively computes, with explicit blocking, the QR factorization of a complex m by n matrix.

2  Specification

#include <nag.h>
#include <nagf08.h>
void  nag_zgeqrt (Nag_OrderType order, Integer m, Integer n, Integer nb, Complex a[], Integer pda, Complex t[], Integer pdt, NagError *fail)

3  Description

nag_zgeqrt (f08apc) forms the QR factorization of an arbitrary rectangular complex m by n matrix. No pivoting is performed.
It differs from nag_zgeqrf (f08asc) in that it: requires an explicit block size; stores reflector factors that are upper triangular matrices of the chosen block size (rather than scalars); and recursively computes the QR factorization based on the algorithm of Elmroth and Gustavson (2000).
If mn, the factorization is given by:
A = Q R 0 ,
where R is an n by n upper triangular matrix (with real diagonal elements) and Q is an m by m unitary matrix. It is sometimes more convenient to write the factorization as
A = Q1 Q2 R 0 ,
which reduces to
A = Q1 R ,
where Q1 consists of the first n columns of Q, and Q2 the remaining m-n columns.
If m<n, R is upper trapezoidal, and the factorization can be written
A = Q R1 R2 ,
where R1 is upper triangular and R2 is rectangular.
The matrix Q is not formed explicitly but is represented as a product of minm,n elementary reflectors (see the f08 Chapter Introduction for details). Functions are provided to work with Q in this representation (see Section 9).
Note also that for any k<n, the information returned represents a QR factorization of the first k columns of the original matrix A.

4  References

Elmroth E and Gustavson F (2000) Applying Recursion to Serial and Parallel QR Factorization Leads to Better Performance IBM Journal of Research and Development. (Volume 44) 4 605–624
Golub G H and Van Loan C F (2012) Matrix Computations (4th Edition) Johns Hopkins University Press, Baltimore

5  Arguments

1:     orderNag_OrderTypeInput
On entry: the order argument specifies the two-dimensional storage scheme being used, i.e., row-major ordering or column-major ordering. C language defined storage is specified by order=Nag_RowMajor. See Section 3.2.1.3 in the Essential Introduction for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2:     mIntegerInput
On entry: m, the number of rows of the matrix A.
Constraint: m0.
3:     nIntegerInput
On entry: n, the number of columns of the matrix A.
Constraint: n0.
4:     nbIntegerInput
On entry: the explicitly chosen block size to be used in computing the QR factorization. See Section 9 for details.
Constraints:
  • nb1;
  • if minm,n>0, nbminm,n.
5:     a[dim]ComplexInput/Output
Note: the dimension, dim, of the array a must be at least
  • max1,pda×n when order=Nag_ColMajor;
  • max1,m×pda when order=Nag_RowMajor.
The i,jth element of the matrix A is stored in
  • a[j-1×pda+i-1] when order=Nag_ColMajor;
  • a[i-1×pda+j-1] when order=Nag_RowMajor.
On entry: the m by n matrix A.
On exit: if mn, the elements below the diagonal are overwritten by details of the unitary matrix Q and the upper triangle is overwritten by the corresponding elements of the n by n upper triangular matrix R.
If m<n, the strictly lower triangular part is overwritten by details of the unitary matrix Q and the remaining elements are overwritten by the corresponding elements of the m by n upper trapezoidal matrix R.
The diagonal elements of R are real.
6:     pdaIntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array a.
Constraints:
  • if order=Nag_ColMajor, pdamax1,m;
  • if order=Nag_RowMajor, pdamax1,n.
7:     t[dim]ComplexOutput
Note: the dimension, dim, of the array t must be at least
  • max1,pdt×minm,n when order=Nag_ColMajor;
  • max1,nb×pdt when order=Nag_RowMajor.
The i,jth element of the matrix T is stored in
  • t[j-1×pdt+i-1] when order=Nag_ColMajor;
  • t[i-1×pdt+j-1] when order=Nag_RowMajor.
On exit: further details of the unitary matrix Q. The number of blocks is b=knb, where k=minm,n and each block is of order nb except for the last block, which is of order k-b-1×nb. For each of the blocks, an upper triangular block reflector factor is computed: T1,T2,,Tb. These are stored in the nb by n matrix T as T=T1|T2||Tb.
8:     pdtIntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array t.
Constraints:
  • if order=Nag_ColMajor, pdtnb;
  • if order=Nag_RowMajor, pdtmax1,minm,n.
9:     failNagError *Input/Output
The NAG error argument (see Section 3.6 in the Essential Introduction).

6  Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_INT
On entry, m=value.
Constraint: m0.
On entry, n=value.
Constraint: n0.
NE_INT_2
On entry, pda=value and m=value.
Constraint: pdamax1,m.
On entry, pda=value and n=value.
Constraint: pdamax1,n.
On entry, pdt=value and nb=value.
Constraint: pdtnb.
NE_INT_3
On entry, nb=value, m=value and n=value.
Constraint: nb1 and
if minm,n>0, nbminm,n.
On entry, pdt=value, m=value and n=value.
Constraint: pdtmax1,minm,n.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.

7  Accuracy

The computed factorization is the exact factorization of a nearby matrix A+E, where
E2 = Oε A2 ,
and ε is the machine precision.

8  Parallelism and Performance

nag_zgeqrt (f08apc) is not threaded by NAG in any implementation.
nag_zgeqrt (f08apc) makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

The total number of real floating-point operations is approximately 83 n2 3m-n  if mn or 83 m2 3n-m  if m<n.
To apply Q to an arbitrary complex rectangular matrix C, nag_zgeqrt (f08apc) may be followed by a call to nag_zgemqrt (f08aqc). For example,
nag_zgemqrt(order,Nag_LeftSide,Nag_ConjTrans,m,p,MIN(m,n),nb,a,pda,
t,pdt,c,pdc,&fail)
forms C=QHC, where C is m by p.
To form the unitary matrix Q explicitly, simply initialize the m by m matrix C to the identity matrix and form C=QC using nag_zgemqrt (f08aqc) as above.
The block size, nb, used by nag_zgeqrt (f08apc) is supplied explicitly through the interface. For moderate and large sizes of matrix, the block size can have a marked effect on the efficiency of the algorithm with the optimal value being dependent on problem size and platform. A value of nb=64minm,n is likely to achieve good efficiency and it is unlikely that an optimal value would exceed 340.
To compute a QR factorization with column pivoting, use nag_ztpqrt (f08bpc) or nag_zgeqpf (f08bsc).
The real analogue of this function is nag_dgeqrt (f08abc).

10  Example

This example solves the linear least squares problems
minimize Axi - bi 2 ,   i=1,2
where b1 and b2 are the columns of the matrix B,
A = 0.96-0.81i -0.03+0.96i -0.91+2.06i -0.05+0.41i -0.98+1.98i -1.20+0.19i -0.66+0.42i -0.81+0.56i 0.62-0.46i 1.01+0.02i 0.63-0.17i -1.11+0.60i -0.37+0.38i 0.19-0.54i -0.98-0.36i 0.22-0.20i 0.83+0.51i 0.20+0.01i -0.17-0.46i 1.47+1.59i 1.08-0.28i 0.20-0.12i -0.07+1.23i 0.26+0.26i
and
B = -2.09+1.93i 3.26-2.70i 3.34-3.53i -6.22+1.16i -4.94-2.04i 7.94-3.13i 0.17+4.23i 1.04-4.26i -5.19+3.63i -2.31-2.12i 0.98+2.53i -1.39-4.05i .

10.1  Program Text

Program Text (f08apce.c)

10.2  Program Data

Program Data (f08apce.d)

10.3  Program Results

Program Results (f08apce.r)


nag_zgeqrt (f08apc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

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