hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_lapack_dgebrd (f08ke)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_dgebrd (f08ke) reduces a real m by n matrix to bidiagonal form.

Syntax

[a, d, e, tauq, taup, info] = f08ke(a, 'm', m, 'n', n)
[a, d, e, tauq, taup, info] = nag_lapack_dgebrd(a, 'm', m, 'n', n)

Description

nag_lapack_dgebrd (f08ke) reduces a real m by n matrix A to bidiagonal form B by an orthogonal transformation: A=QBPT, where Q and PT are orthogonal matrices of order m and n respectively.
If mn, the reduction is given by:
A =Q B1 0 PT = Q1 B1 PT ,  
where B1 is an n by n upper bidiagonal matrix and Q1 consists of the first n columns of Q.
If m<n, the reduction is given by
A =Q B1 0 PT = Q B1 P1T ,  
where B1 is an m by m lower bidiagonal matrix and P1T consists of the first m rows of PT.
The orthogonal matrices Q and P are not formed explicitly but are represented as products of elementary reflectors (see the F08 Chapter Introduction for details). Functions are provided to work with Q and P in this representation (see Further Comments).

References

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

Parameters

Compulsory Input Parameters

1:     alda: – double array
The first dimension of the array a must be at least max1,m.
The second dimension of the array a must be at least max1,n.
The m by n matrix A.

Optional Input Parameters

1:     m int64int32nag_int scalar
Default: the first dimension of the array a.
m, the number of rows of the matrix A.
Constraint: m0.
2:     n int64int32nag_int scalar
Default: the second dimension of the array a.
n, the number of columns of the matrix A.
Constraint: n0.

Output Parameters

1:     alda: – double array
The first dimension of the array a will be max1,m.
The second dimension of the array a will be max1,n.
If mn, the diagonal and first superdiagonal store the upper bidiagonal matrix B, elements below the diagonal store details of the orthogonal matrix Q and elements above the first superdiagonal store details of the orthogonal matrix P.
If m<n, the diagonal and first subdiagonal store the lower bidiagonal matrix B, elements below the first subdiagonal store details of the orthogonal matrix Q and elements above the diagonal store details of the orthogonal matrix P.
2:     d: – double array
The dimension of the array d will be max1,minm,n
The diagonal elements of the bidiagonal matrix B.
3:     e: – double array
The dimension of the array e will be max1,minm,n-1
The off-diagonal elements of the bidiagonal matrix B.
4:     tauq: – double array
The dimension of the array tauq will be max1,minm,n
Further details of the orthogonal matrix Q.
5:     taup: – double array
The dimension of the array taup will be max1,minm,n
Further details of the orthogonal matrix P.
6:     info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

   info=-i
If info=-i, parameter i had an illegal value on entry. The parameters are numbered as follows:
1: m, 2: n, 3: a, 4: lda, 5: d, 6: e, 7: tauq, 8: taup, 9: work, 10: lwork, 11: info.
It is possible that info refers to a parameter that is omitted from the MATLAB interface. This usually indicates that an error in one of the other input parameters has caused an incorrect value to be inferred.

Accuracy

The computed bidiagonal form B satisfies QBPT=A+E, where
E2 c n ε A2 ,  
cn is a modestly increasing function of n, and ε is the machine precision.
The elements of B themselves may be sensitive to small perturbations in A or to rounding errors in the computation, but this does not affect the stability of the singular values and vectors.

Further Comments

The total number of floating-point operations is approximately 43n23m-n if mn or 43m23n-m if m<n.
If mn, it can be more efficient to first call nag_lapack_dgeqrf (f08ae) to perform a QR factorization of A, and then to call nag_lapack_dgebrd (f08ke) to reduce the factor R to bidiagonal form. This requires approximately 2n2m+n floating-point operations.
If mn, it can be more efficient to first call nag_lapack_dgelqf (f08ah) to perform an LQ factorization of A, and then to call nag_lapack_dgebrd (f08ke) to reduce the factor L to bidiagonal form. This requires approximately 2m2m+n operations.
To form the orthogonal matrices PT and/or Q nag_lapack_dgebrd (f08ke) may be followed by calls to nag_lapack_dorgbr (f08kf):
to form the m by m orthogonal matrix Q 
[a, info] = f08kf('Q', k, a, tauq);
but note that the second dimension of the array a must be at least m, which may be larger than was required by nag_lapack_dgebrd (f08ke);
to form the n by n orthogonal matrix PT 
[a, info] = f08kf('P', k, a, taup);
but note that the first dimension of the array a, specified by the argument lda, must be at least n, which may be larger than was required by nag_lapack_dgebrd (f08ke).
To apply Q or P to a real rectangular matrix C, nag_lapack_dgebrd (f08ke) may be followed by a call to nag_lapack_dormbr (f08kg).
The complex analogue of this function is nag_lapack_zgebrd (f08ks).

Example

This example reduces the matrix A to bidiagonal form, 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 .  
function f08ke_example


fprintf('f08ke example results\n\n');

m = 6;
n = int64(4);
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];

% Reduce A to bidiagonal form
[~, d, e, tauq, taup, info] = f08ke(a);

fprintf(' Bidiagonal matrix B\n   Main diagonal  ');
fprintf(' %7.3f',d);
fprintf('\n   super-diagonal ');
fprintf(' %7.3f',e);
fprintf('\n');


f08ke example results

 Bidiagonal matrix B
   Main diagonal     3.618   2.416  -1.921  -1.427
   super-diagonal    1.259   1.526  -1.189

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015