NAG Library Routine Document
f08vqf (zggsvd3)
1
Purpose
f08vqf (zggsvd3) computes the generalized singular value decomposition (GSVD) of an by complex matrix and a by complex matrix .
2
Specification
Fortran Interface
Subroutine f08vqf ( |
jobu, jobv, jobq, m, n, p, k, l, a, lda, b, ldb, alpha, beta, u, ldu, v, ldv, q, ldq, work, lwork, rwork, iwork, info) |
Integer, Intent (In) | :: | m, n, p, lda, ldb, ldu, ldv, ldq, lwork | Integer, Intent (Out) | :: | k, l, iwork(n), info | Real (Kind=nag_wp), Intent (Out) | :: | alpha(n), beta(n), rwork(2*n) | Complex (Kind=nag_wp), Intent (Inout) | :: | a(lda,*), b(ldb,*), u(ldu,*), v(ldv,*), q(ldq,*) | Complex (Kind=nag_wp), Intent (Out) | :: | work(max(1,lwork)) | Character (1), Intent (In) | :: | jobu, jobv, jobq |
|
C Header Interface
#include <nagmk26.h>
void |
f08vqf_ (const char *jobu, const char *jobv, const char *jobq, const Integer *m, const Integer *n, const Integer *p, Integer *k, Integer *l, Complex a[], const Integer *lda, Complex b[], const Integer *ldb, double alpha[], double beta[], Complex u[], const Integer *ldu, Complex v[], const Integer *ldv, Complex q[], const Integer *ldq, Complex work[], const Integer *lwork, double rwork[], Integer iwork[], Integer *info, const Charlen length_jobu, const Charlen length_jobv, const Charlen length_jobq) |
|
The routine may be called by its
LAPACK
name zggsvd3.
3
Description
Given an
by
complex matrix
and a
by
complex matrix
, the generalized singular value decomposition is given by
where
,
and
are unitary matrices. Let
be the effective numerical rank of
and
be the effective numerical rank of the matrix
, then the first
generalized singular values are infinite and the remaining
are finite.
is a
by
nonsingular upper triangular matrix,
and
are
by
and
by
‘diagonal’ matrices structured as follows:
if
,
where
and
is stored as a submatrix of
with elements
stored as
on exit.
If
,
where
and
is stored as a submatrix of
with
stored as
, and
is stored as a submatrix of
with
stored as
.
The routine computes , , and, optionally, the unitary transformation matrices , and .
In particular, if
is an
by
nonsingular matrix, then the GSVD of
and
implicitly gives the SVD of
:
If
has orthonormal columns, then the GSVD of
and
is also equal to the CS decomposition of
and
. Furthermore, the GSVD can be used to derive the solution of the eigenvalue problem:
In some literature, the GSVD of
and
is presented in the form
where
and
are orthogonal and
is nonsingular, and
and
are ‘diagonal’. The former GSVD form can be converted to the latter form by setting
A two stage process is used to compute the GSVD of the matrix pair
. The pair is first reduced to upper triangular form by unitary transformations using
f08vuf (zggsvp3). The GSVD of the resulting upper triangular matrix pair is then performed by
f08ysf (ztgsja) which uses a variant of the Kogbetliantz algorithm (a cyclic Jacobi method).
4
References
Anderson E, Bai Z, Bischof C, Blackford S, Demmel J, Dongarra J J, Du Croz J J, Greenbaum A, Hammarling S, McKenney A and Sorensen D (1999)
LAPACK Users' Guide (3rd Edition) SIAM, Philadelphia
http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (2012) Matrix Computations (4th Edition) Johns Hopkins University Press, Baltimore
5
Arguments
- 1: – Character(1)Input
-
On entry: if
, the unitary matrix
is computed.
If , is not computed.
Constraint:
or .
- 2: – Character(1)Input
-
On entry: if
, the unitary matrix
is computed.
If , is not computed.
Constraint:
or .
- 3: – Character(1)Input
-
On entry: if
, the unitary matrix
is computed.
If , is not computed.
Constraint:
or .
- 4: – IntegerInput
-
On entry: , the number of rows of the matrix .
Constraint:
.
- 5: – IntegerInput
-
On entry: , the number of columns of the matrices and .
Constraint:
.
- 6: – IntegerInput
-
On entry: , the number of rows of the matrix .
Constraint:
.
- 7: – IntegerOutput
- 8: – IntegerOutput
-
On exit:
k and
l specify the dimension of the subblocks
and
as described in
Section 3;
is the effective numerical rank of
.
- 9: – Complex (Kind=nag_wp) arrayInput/Output
-
Note: the second dimension of the array
a
must be at least
.
On entry: the by matrix .
On exit: contains the triangular matrix
, or part of
. See
Section 3 for details.
- 10: – IntegerInput
-
On entry: the first dimension of the array
a as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
Constraint:
.
- 11: – Complex (Kind=nag_wp) arrayInput/Output
-
Note: the second dimension of the array
b
must be at least
.
On entry: the by matrix .
On exit: contains the triangular matrix
if
. See
Section 3 for details.
- 12: – IntegerInput
-
On entry: the first dimension of the array
b as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
Constraint:
.
- 13: – Real (Kind=nag_wp) arrayOutput
-
On exit: see the description of
beta.
- 14: – Real (Kind=nag_wp) arrayOutput
-
On exit:
alpha and
beta contain the generalized singular value pairs of
and
,
and
;
- ,
- ,
and if
,
- ,
- ,
or if
,
- ,
- ,
- ,
- , and
- ,
- .
The notation above refers to consecutive elements
, for .
- 15: – Complex (Kind=nag_wp) arrayOutput
-
Note: the second dimension of the array
u
must be at least
if
, and at least
otherwise.
On exit: if
,
u contains the
by
unitary matrix
.
If
,
u is not referenced.
- 16: – IntegerInput
-
On entry: the first dimension of the array
u as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
Constraints:
- if , ;
- otherwise .
- 17: – Complex (Kind=nag_wp) arrayOutput
-
Note: the second dimension of the array
v
must be at least
if
, and at least
otherwise.
On exit: if
,
v contains the
by
unitary matrix
.
If
,
v is not referenced.
- 18: – IntegerInput
-
On entry: the first dimension of the array
v as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
Constraints:
- if , ;
- otherwise .
- 19: – Complex (Kind=nag_wp) arrayOutput
-
Note: the second dimension of the array
q
must be at least
if
, and at least
otherwise.
On exit: if
,
q contains the
by
unitary matrix
.
If
,
q is not referenced.
- 20: – IntegerInput
-
On entry: the first dimension of the array
q as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
Constraints:
- if , ;
- otherwise .
- 21: – Complex (Kind=nag_wp) arrayWorkspace
-
On exit: if
, the real part of
contains the minimum value of
lwork required for optimal performance.
- 22: – IntegerInput
-
On entry: the dimension of the array
work as declared in the (sub)program from which
f08vqf (zggsvd3) is called.
If
, a workspace query is assumed; the routine only calculates the optimal size of the
work array, returns this value as the first entry of the
work array, and no error message related to
lwork is issued.
Suggested value:
for optimal performance,
lwork must generally be larger than the minimum; increase workspace by, say,
, where
is the optimal
block size.
Constraints:
- if , ;
- if , .
- 23: – Real (Kind=nag_wp) arrayWorkspace
-
- 24: – Integer arrayOutput
-
On exit: stores the sorting information. More precisely, if
is the ordered set of indices of
alpha containing
(denote as
, see
beta), then the corresponding elements
contain the swap pivots,
, that sorts
such that
is in descending numerical order.
The following pseudocode sorts the set :
- 25: – IntegerOutput
On exit:
unless the routine detects an error (see
Section 6).
6
Error Indicators and Warnings
If , argument had an illegal value. An explanatory message is output, and execution of the program is terminated.
-
The Jacobi-type procedure failed to converge.
7
Accuracy
The computed generalized singular value decomposition is nearly the exact generalized singular value decomposition for nearby matrices
and
, where
and
is the
machine precision. See Section 4.12 of
Anderson et al. (1999) for further details.
8
Parallelism and Performance
f08vqf (zggsvd3) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08vqf (zggsvd3) 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 routine. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
This routine replaces the deprecated routine
f08vnf (zggsvd) which used an unblocked algorithm and therefore did not make best use of level
BLAS routines.
The diagonal elements of the matrix are real.
The real analogue of this routine is
f08vcf (dggsvd3).
10
Example
This example finds the generalized singular value decomposition
where
and
together with estimates for the condition number of
and the error bound for the computed generalized singular values.
The example program assumes that , and would need slight modification if this is not the case.
10.1
Program Text
Program Text (f08vqfe.f90)
10.2
Program Data
Program Data (f08vqfe.d)
10.3
Program Results
Program Results (f08vqfe.r)