Computes the CS decomposition of an orthogonal/unitary matrix in bidiagonal-block form.
FORTRAN 77:
call sbbcsd( jobu1, jobu2, jobv1t, jobv2t, trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2, v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e, b21d, b21e, b21e, b22e, work, lwork, info )
call dbbcsd( jobu1, jobu2, jobv1t, jobv2t, trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2, v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e, b21d, b21e, b21e, b22e, work, lwork, info )
call cbbcsd( jobu1, jobu2, jobv1t, jobv2t, trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2, v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e, b21d, b21e, b21e, b22e, rwork, rlwork, info )
call zbbcsd( jobu1, jobu2, jobv1t, jobv2t, trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2, v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e, b21d, b21e, b21e, b22e, rwork, rlwork, info )
Fortran 95:
call bbcsd( theta,phi,u1,u2,v1t,v2t[,b11d][,b11e][,b12d][,b12e][,b21d][,b21e][,b22d][,b22e][,jobu1][,jobu2][,jobv1t][,jobv2t][,trans][,info] )
C:
lapack_int LAPACKE_sbbcsd( int matrix_order, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, float* theta, float* phi, float* u1, lapack_int ldu1, float* u2, lapack_int ldu2, float* v1t, lapack_int ldv1t, float* v2t, lapack_int ldv2t, float* b11d, float* b11e, float* b12d, float* b12e, float* b21d, float* b21e, float* b22d, float* b22e );
lapack_int LAPACKE_dbbcsd( int matrix_order, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, double* theta, double* phi, double* u1, lapack_int ldu1, double* u2, lapack_int ldu2, double* v1t, lapack_int ldv1t, double* v2t, lapack_int ldv2t, double* b11d, double* b11e, double* b12d, double* b12e, double* b21d, double* b21e, double* b22d, double* b22e );
lapack_int LAPACKE_cbbcsd( int matrix_order, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, float* theta, float* phi, lapack_complex_float* u1, lapack_int ldu1, lapack_complex_float* u2, lapack_int ldu2, lapack_complex_float* v1t, lapack_int ldv1t, lapack_complex_float* v2t, lapack_int ldv2t, float* b11d, float* b11e, float* b12d, float* b12e, float* b21d, float* b21e, float* b22d, float* b22e );
lapack_int LAPACKE_zbbcsd( int matrix_order, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, double* theta, double* phi, lapack_complex_double* u1, lapack_int ldu1, lapack_complex_double* u2, lapack_int ldu2, lapack_complex_double* v1t, lapack_int ldv1t, lapack_complex_double* v2t, lapack_int ldv2t, double* b11d, double* b11e, double* b12d, double* b12e, double* b21d, double* b21e, double* b22d, double* b22e );
The FORTRAN 77 interfaces are specified in the mkl_lapack.fi and mkl_lapack.h include files, the Fortran 95 interfaces are specified in the lapack.f90 include file, and the C interfaces are specified in the mkl_lapacke.h include file.
mkl_lapack.fiThe routine ?bbcsd computes the CS decomposition of an orthogonal or unitary matrix in bidiagonal-block form:
or
respectively.
x is m-by-m with the top-left block p-by-q. Note that q must not be larger than p, m-p, or m-q. If q is not the smallest index, x must be transposed and/or permuted in constant time using the trans option. See ?orcsd/?uncsd for details.
The bidiagonal matrices b11, b12, b21, and b22 are represented implicitly by angles theta(1:q) and phi(1:q-1).
The orthogonal/unitary matrices u1, u2, v1t, and v2t are input/output. The input matrices are pre- or post-multiplied by the appropriate singular vector matrices.
The data types are given for the Fortran interface. A <datatype> placeholder, if present, is used for the C interface data types in the C interface section above. See C Interface Conventions for the C interface principal conventions and type definitions.
CHARACTER. If equals Y, then u1 is updated. Otherwise, u1 is not updated.
CHARACTER. If equals Y, then u2 is updated. Otherwise, u2 is not updated.
CHARACTER. If equals Y, then v1t is updated. Otherwise, v1t is not updated.
CHARACTER. If equals Y, then v2t is updated. Otherwise, v2t is not updated.
CHARACTER
INTEGER. The number of rows and columns of the orthogonal/unitary matrix X in bidiagonal-block form.
INTEGER. The number of rows in the top-left block of x. 0 ? p ? m.
INTEGER. The number of columns in the top-left block of x. 0 ? q ? min(p,m-p,m-q).
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q).
On entry, the angles theta(1), ..., theta(q) that, along with phi(1), ..., phi(q-1), define the matrix in bidiagonal-block form.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q-1).
The angles phi(1), ..., phi(q-1) that, along with theta(1), ..., theta(q), define the matrix in bidiagonal-block form.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (ldu1,p).
On entry, an ldu1-by-p matrix.
INTEGER. The leading dimension of the array u1.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (ldu2,m-p).
On entry, an ldu2-by-(m-p) matrix.
INTEGER. The leading dimension of the array u2.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (ldv1t,q).
On entry, an ldv1t-by-q matrix.
INTEGER. The leading dimension of the array v1t.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (ldv2t,m-q).
On entry, an ldv2t-by-(m-q) matrix.
INTEGER. The leading dimension of the array v2t.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Workspace array, DIMENSION (max(1,lwork)).
INTEGER. The size of the work array. lwork ? max(1,8*q)
If lwork = -1, then 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 by xerbla.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
On exit, the angles whose cosines and sines define th edaigonal blocks in the CS decomposition.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
On exit, u1 is postmultiplied by the left singular vector matrix common to [ b11 ; 0 ] and [ b12 0 0 ; 0 -I 0 0 ].
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
On exit, u2 is postmultiplied by the left singular vector matrix common to [ b21 ; 0 ] and [ b22 0 0 ; 0 0 I ].
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q).
On exit, v1t is premultiplied by the transpose of the right singular vector matrix common to [ b11 ; 0 ] and [ b21 ; 0 ].
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
On exit, v2t is premultiplied by the transpose of the right singular vector matrix common to [ b12 0 0 ; 0 -I 0 ] and [ b22 0 0 ; 0 0 I ].
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q).
When ?bbcsd converges, b11d contains the cosines of theta(1), ..., theta(q). If ?bbcsd fails to converge, b11d contains the diagonal of the partially reduced top left block.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q-1).
When ?bbcsd converges, b11e contains zeros. If ?bbcsd fails to converge, b11e contains the superdiagonal of the partially reduced top left block.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q).
When ?bbcsd converges, b12d contains the negative sines of theta(1), ..., theta(q). If ?bbcsd fails to converge, b12d contains the diagonal of the partially reduced top right block.
REAL for sbbcsd
DOUBLE PRECISION for dbbcsd
COMPLEX for cbbcsd
DOUBLE COMPLEX for zbbcsd
Array, DIMENSION (q-1).
When ?bbcsd converges, b12e contains zeros. If ?bbcsd fails to converge, b11e contains the superdiagonal of the partially reduced top right block.
INTEGER.
= 0: successful exit
< 0: if info = -i, the i-th argument has an illegal value
> 0: if ?bbcsd did not converge, info specifies the number of nonzero entries in phi, and b11d, b11e, etc. and contains the partially reduced matrix.
Routines in Fortran 95 interface have fewer arguments in the calling sequence than their FORTRAN 77 counterparts. For general conventions applied to skip redundant or reconstructible arguments, see Fortran 95 Interface Conventions.
Specific details for the routine ?bbcsd interface are as follows:
Holds the vector of length q.
Holds the vector of length q-1.
Holds the matrix of size (p,p).
Holds the matrix of size (m-p,m-p).
Holds the matrix of size (q,q).
Holds the matrix of size (m-q,m-q).
Holds the vector of length q.
Holds the vector of length q-1.
Holds the vector of length q.
Holds the vector of length q-1.
Holds the vector of length q.
Holds the vector of length q-1.
Holds the vector of length q.
Holds the vector of length q-1.
Indicates whether u1 is computed. Must be 'Y' or 'O'.
Indicates whether u2 is computed. Must be 'Y' or 'O'.
Indicates whether v1t is computed. Must be 'Y' or 'O'.
Indicates whether v2t is computed. Must be 'Y' or 'O'.
Must be 'N' or 'T'.
Copyright © 1994 - 2011, Intel Corporation. All rights reserved.