NAG Library Function Document

nag_ztpqrt (f08bpc)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

nag_ztpqrt (f08bpc) computes the QR factorization of a complex m+n by n triangular-pentagonal matrix.

2
Specification

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

3
Description

nag_ztpqrt (f08bpc) forms the QR factorization of a complex m+n by n triangular-pentagonal matrix C,
C= A B  
where A is an upper triangular n by n matrix and B is an m by n pentagonal matrix consisting of an m-l by n rectangular matrix B1 on top of an l by n upper trapezoidal matrix B2:
B= B1 B2 .  
The upper trapezoidal matrix B2 consists of the first l rows of an n by n upper triangular matrix, where 0lminm,n. If l=0, B is m by n rectangular; if l=n and m=n, B is upper triangular.
A recursive, explicitly blocked, QR factorization (see nag_zgeqrt (f08apc)) is performed on the matrix C. The upper triangular matrix R, details of the unitary matrix Q, and further details (the block reflector factors) of Q are returned.
Typically the matrix A or B2 contains the matrix R from the QR factorization of a subproblem and nag_ztpqrt (f08bpc) performs the QR update operation from the inclusion of matrix B1.
For example, consider the QR factorization of an l by n matrix B^ with l<n: B^ = Q^R^ , R^ = R1^ R2^ , where R1^ is l by l upper triangular and R2^ is n-l by n rectangular (this can be performed by nag_zgeqrt (f08apc)). Given an initial least squares problem B^ X^ = Y^  where X and Y are l by nrhs matrices, we have R^ X^ = Q^H Y^ .
Now, adding an additional m-l rows to the original system gives the augmented least squares problem
BX=Y  
where B is an m by n matrix formed by adding m-l rows on top of R^ and Y is an m by nrhs matrix formed by adding m-l rows on top of Q^HY^.
nag_ztpqrt (f08bpc) can then be used to perform the QR factorization of the pentagonal matrix B; the n by n matrix A will be zero on input and contain R on output.
In the case where B^ is r by n, rn, R^ is n by n upper triangular (forming A) on top of r-n rows of zeros (forming first r-n rows of B). Augmentation is then performed by adding rows to the bottom of B with l=0.

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:     order Nag_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.3.1.3 in How to Use the NAG Library and its Documentation for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2:     m IntegerInput
On entry: m, the number of rows of the matrix B.
Constraint: m0.
3:     n IntegerInput
On entry: n, the number of columns of the matrix B and the order of the upper triangular matrix A.
Constraint: n0.
4:     l IntegerInput
On entry: l, the number of rows of the trapezoidal part of B (i.e., B2).
Constraint: 0lminm,n.
5:     nb IntegerInput
On entry: the explicitly chosen block-size to be used in the algorithm for computing the QR factorization. See Section 9 for details.
Constraints:
  • nb1;
  • if n>0, nbn.
6:     a[dim] ComplexInput/Output
Note: the dimension, dim, of the array a must be at least max1,pda×n.
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 n by n upper triangular matrix A.
On exit: the upper triangle is overwritten by the corresponding elements of the n by n upper triangular matrix R.
7:     pda IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array a.
Constraint: pdamax1,n.
8:     b[dim] ComplexInput/Output
Note: the dimension, dim, of the array b must be at least
  • max1,pdb×n when order=Nag_ColMajor;
  • max1,m×pdb when order=Nag_RowMajor.
The i,jth element of the matrix B is stored in
  • b[j-1×pdb+i-1] when order=Nag_ColMajor;
  • b[i-1×pdb+j-1] when order=Nag_RowMajor.
On entry: the m by n pentagonal matrix B composed of an m-l by n rectangular matrix B1 above an l by n upper trapezoidal matrix B2.
On exit: details of the unitary matrix Q.
9:     pdb IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array b.
Constraints:
  • if order=Nag_ColMajor, pdbmax1,m;
  • if order=Nag_RowMajor, pdbmax1,n.
10:   t[dim] ComplexOutput
Note: the dimension, dim, of the array t must be at least
  • max1,pdt×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.
11:   pdt IntegerInput
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, pdtn.
12:   fail NagError *Input/Output
The NAG error argument (see Section 3.7 in How to Use the NAG Library and its Documentation).

6
Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
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, nb=value and n=value.
Constraint: nb1 and
if n>0, nbn.
On entry, pda=value and n=value.
Constraint: pdamax1,n.
On entry, pdb=value and m=value.
Constraint: pdbmax1,m.
On entry, pdb=value and n=value.
Constraint: pdbmax1,n.
On entry, pdt=value and n=value.
Constraint: pdtn.
On entry, pdt=value and nb=value.
Constraint: pdtnb.
NE_INT_3
On entry, l=value, m=value and n=value.
Constraint: 0lminm,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.
See Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.

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_ztpqrt (f08bpc) 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 x06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9
Further Comments

The total number of floating-point operations is approximately 23 n2 3m-n  if mn or 23 m2 3n-m  if m<n.
The block size, nb, used by nag_ztpqrt (f08bpc) 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 apply Q to an arbitrary complex rectangular matrix C, nag_ztpqrt (f08bpc) may be followed by a call to nag_ztpmqrt (f08bqc). For example,
nag_ztpmqrt(Nag_ColMajor,Nag_LeftSide,Nag_Trans,m,p,n,l,nb,b,pdb,
t,pdt,c,pdc,&c(n+1,1),ldc,&fail)
forms C=QHC, where C is m+n by p.
To form the unitary matrix Q explicitly set p=m+n, initialize C to the identity matrix and make a call to nag_ztpmqrt (f08bqc) as above.

10
Example

This example finds the basic solutions for 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 .  
A QR factorization is performed on the first 4 rows of A using nag_zgeqrt (f08apc) after which the first 4 rows of B are updated by applying QT using nag_zgemqrt (f08aqc). The remaining row is added by performing a QR update using nag_ztpqrt (f08bpc); B is updated by applying the new QT using nag_ztpmqrt (f08bqc); the solution is finally obtained by triangular solve using R from the updated QR.

10.1
Program Text

Program Text (f08bpce.c)

10.2
Program Data

Program Data (f08bpce.d)

10.3
Program Results

Program Results (f08bpce.r)

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