nag_zgeqrf (f08asc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_zgeqrf (f08asc)

 Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_zgeqrf (f08asc) computes the QR factorization of a complex m by n matrix.

2  Specification

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

3  Description

nag_zgeqrf (f08asc) forms the QR factorization of an arbitrary rectangular complex m by n matrix. No pivoting is performed.
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 = Q1R ,  
where Q1 consists of the first n columns of Q, and Q2 the remaining m-n columns.
If m<n, R is 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 in the first k columns of the array a represents a QR factorization of the first k columns of the original matrix A.

4  References

Golub G H and Van Loan C F (1996) Matrix Computations (3rd 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 2.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 A.
Constraint: m0.
3:     n IntegerInput
On entry: n, the number of columns of the matrix A.
Constraint: n0.
4:     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.
5:     pda IntegerInput
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.
6:     tau[dim] ComplexOutput
Note: the dimension, dim, of the array tau must be at least max1,minm,n.
On exit: further details of the unitary matrix Q.
7:     fail NagError *Input/Output
The NAG error argument (see Section 2.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.
On entry, pda=value.
Constraint: pda>0.
NE_INT_2
On entry, pda=value and m=value.
Constraint: pdamax1,m.
On entry, pda=value and n=value.
Constraint: pdamax1,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.
An unexpected error has been triggered by this function. Please contact NAG.
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_zgeqrf (f08asc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_zgeqrf (f08asc) 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 real floating-point operations is approximately 83 n2 3m-n  if mn or 83 m2 3n-m  if m<n.
To form the unitary matrix Q nag_zgeqrf (f08asc) may be followed by a call to nag_zungqr (f08atc):
nag_zungqr(order,m,m,MIN(m,n),&a,pda,tau,&fail)
but note that the second dimension of the array a must be at least m, which may be larger than was required by nag_zgeqrf (f08asc).
When mn, it is often only the first n columns of Q that are required, and they may be formed by the call:
nag_zungqr(order,m,n,n,&a,pda,tau,&fail)
To apply Q to an arbitrary complex rectangular matrix C, nag_zgeqrf (f08asc) may be followed by a call to nag_zunmqr (f08auc). For example,
nag_zunmqr(order,Nag_LeftSide,Nag_ConjTrans,m,p,MIN(m,n),&a,pda,
  tau,&c,pdc,&fail)
forms C=QHC, where C is m by p.
To compute a QR factorization with column pivoting, use nag_zgeqpf (f08bsc).
The real analogue of this function is nag_dgeqrf (f08aec).

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 = -1.54+0.76i 3.17-2.09i 0.12-1.92i -6.53+4.18i -9.08-4.31i 7.28+0.73i 7.49+3.65i 0.91-3.97i -5.63-2.12i -5.46-1.64i 2.37+8.03i -2.84-5.86i .  

10.1  Program Text

Program Text (f08asce.c)

10.2  Program Data

Program Data (f08asce.d)

10.3  Program Results

Program Results (f08asce.r)


nag_zgeqrf (f08asc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

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