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_dormbr (f08kg)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_dormbr (f08kg) multiplies an arbitrary real m by n matrix C by one of the real orthogonal matrices Q or P which were determined by nag_lapack_dgebrd (f08ke) when reducing a real matrix to bidiagonal form.

Syntax

[c, info] = f08kg(vect, side, trans, k, a, tau, c, 'm', m, 'n', n)
[c, info] = nag_lapack_dormbr(vect, side, trans, k, a, tau, c, 'm', m, 'n', n)

Description

nag_lapack_dormbr (f08kg) is intended to be used after a call to nag_lapack_dgebrd (f08ke), which reduces a real rectangular matrix A to bidiagonal form B by an orthogonal transformation: A=QBPT. nag_lapack_dgebrd (f08ke) represents the matrices Q and PT as products of elementary reflectors.
This function may be used to form one of the matrix products
QC , QTC , CQ , CQT , PC , PTC , CP ​ or ​ CPT ,  
overwriting the result on C (which may be any real rectangular matrix).

References

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

Parameters

Note: in the descriptions below, r denotes the order of Q or PT: if side='L', r=m and if side='R', r=n.

Compulsory Input Parameters

1:     vect – string (length ≥ 1)
Indicates whether Q or QT or P or PT is to be applied to C.
vect='Q'
Q or QT is applied to C.
vect='P'
P or PT is applied to C.
Constraint: vect='Q' or 'P'.
2:     side – string (length ≥ 1)
Indicates how Q or QT or P or PT is to be applied to C.
side='L'
Q or QT or P or PT is applied to C from the left.
side='R'
Q or QT or P or PT is applied to C from the right.
Constraint: side='L' or 'R'.
3:     trans – string (length ≥ 1)
Indicates whether Q or P or QT or PT is to be applied to C.
trans='N'
Q or P is applied to C.
trans='T'
QT or PT is applied to C.
Constraint: trans='N' or 'T'.
4:     k int64int32nag_int scalar
If vect='Q', the number of columns in the original matrix A.
If vect='P', the number of rows in the original matrix A.
Constraint: k0.
5:     alda: – double array
The first dimension, lda, of the array a must satisfy
  • if vect='Q', lda max1,r ;
  • if vect='P', lda max1,minr,k .
The second dimension of the array a must be at least max1,minr,k  if vect='Q' and at least max1,r if vect='P'.
Details of the vectors which define the elementary reflectors, as returned by nag_lapack_dgebrd (f08ke).
6:     tau: – double array
The dimension of the array tau must be at least max1,minr,k
Further details of the elementary reflectors, as returned by nag_lapack_dgebrd (f08ke) in its argument tauq if vect='Q', or in its argument taup if vect='P'.
7:     cldc: – double array
The first dimension of the array c must be at least max1,m.
The second dimension of the array c must be at least max1,n.
The matrix C.

Optional Input Parameters

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

Output Parameters

1:     cldc: – double array
The first dimension of the array c will be max1,m.
The second dimension of the array c will be max1,n.
c stores QC or QTC or CQ or CTQ or PC or PTC or CP or CTP as specified by vect, side and trans.
2:     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: vect, 2: side, 3: trans, 4: m, 5: n, 6: k, 7: a, 8: lda, 9: tau, 10: c, 11: ldc, 12: work, 13: lwork, 14: 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 result differs from the exact result by a matrix E such that
E2 = Oε C2 ,  
where ε is the machine precision.

Further Comments

The total number of floating-point operations is approximately where k is the value of the argument k
The complex analogue of this function is nag_lapack_zunmbr (f08ku).

Example

For this function two examples are presented. Both illustrate how the reduction to bidiagonal form of a matrix A may be preceded by a QR or LQ factorization of A.
In the first example, m>n, and
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 .  
The function first performs a QR factorization of A as A=QaR and then reduces the factor R to bidiagonal form B: R=QbBPT. Finally it forms Qa and calls nag_lapack_dormbr (f08kg) to form Q=QaQb.
In the second example, m<n, and
A = -5.42 3.28 -3.68 0.27 2.06 0.46 -1.65 -3.40 -3.20 -1.03 -4.06 -0.01 -0.37 2.35 1.90 4.31 -1.76 1.13 -3.15 -0.11 1.99 -2.70 0.26 4.50 .  
The function first performs an LQ factorization of A as A=LPaT and then reduces the factor L to bidiagonal form B: L=QBPbT. Finally it forms PbT and calls nag_lapack_dormbr (f08kg) to form PT=PbTPaT.
function f08kg_example


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

% Two cases of preceding reduction to bidiagonal form by QR or LQ
% Case 1: m > n, precede by QR
ex1;
% Case 2: m < n, precede by LQ
ex2;



function ex1
  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];
  % Factorize A = QR
  [QR, tau, info] = f08ae(a);

  % Generate Q from QR
  [Q, info] = f08af(QR, tau);

  % Extract R from QR
  R = triu(QR(1:n,1:n));

  % Bidiagonalize R = Q1 B P^T
  [B, d, e, tauq, taup, info] = ...
  f08ke(R);

  % Update Q: Q2 = Q*Q1 (so A = QR = Q2 B P^T)
  vect = 'Q';
  side = 'Right';
  trans = 'No transpose';
  [Q2, info] = f08kg( ...
                      vect, side, trans, n, B, tauq, Q);

  fprintf('Example 1: bidiagonal matrix B\n   Main diagonal  ');
  fprintf(' %7.3f',d);
  fprintf('\n   super-diagonal ');
  fprintf(' %7.3f',e);
  fprintf('\n\n');
  disp('Example 1: Orthogonal matrix Q');
  disp(Q2);

function ex2
  m = int64(4);
  n = int64(6);
  a = [ -5.42   3.28  -3.68   0.27   2.06   0.46;
        -1.65  -3.40  -3.20  -1.03  -4.06  -0.01;
        -0.37   2.35   1.90   4.31  -1.76   1.13;
        -3.15  -0.11   1.99  -2.70   0.26   4.50];

  % Factorize A = LQ
  [LQ, tau, info] = f08ah(a);

  % Generate Q from LQ
  [Q, info] = f08aj(LQ, tau);

  % Extract L from LQ
  L = tril(LQ(1:m,1:m));

  % Bidiagonalize L = Q1 B P^T
  [B, d, e, tauq, taup, info] = ...
  f08ke(L);

  % Update Q: P2 = P^T*Q (so A = LQ = Q1 B P2)
  vect = 'P';
  side = 'Left';
  trans = 'Transpose';
  [P2, info] = f08kg( ...
                      vect, side, trans, n, B, taup, Q);

  fprintf('Example 2: bidiagonal matrix B\n   Main diagonal  ');
  fprintf(' %7.3f',d);
  fprintf('\n   super-diagonal ');
  fprintf(' %7.3f',e);
  fprintf('\n\n');
  disp('Example 2: Orthogonal matrix P^T');
  disp(P2);
f08kg example results

Example 1: bidiagonal matrix B
   Main diagonal     3.618  -2.416   1.921  -1.427
   super-diagonal    1.259  -1.526   1.189

Example 1: Orthogonal matrix Q
   -0.1576   -0.2690    0.2612    0.8513
   -0.5335    0.5311   -0.2922    0.0184
    0.6358    0.3495   -0.0250   -0.0210
   -0.5335    0.0035    0.1537   -0.2592
    0.0415    0.5572   -0.2917    0.4523
   -0.0055    0.4614    0.8585   -0.0532

Example 2: bidiagonal matrix B
   Main diagonal    -7.772   6.157  -6.058   5.793
   super-diagonal    1.193   0.573  -1.914

Example 2: Orthogonal matrix P^T
   -0.7104    0.4299   -0.4824    0.0354    0.2700    0.0603
    0.3583    0.1382   -0.4110    0.4044    0.0951   -0.7148
   -0.0507    0.4244    0.3795    0.7402   -0.2773    0.2203
    0.2442    0.4016    0.4158   -0.1354    0.7666   -0.0137


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