Visible to Intel only — GUID: GUID-F52443B0-8063-4E7B-A7ED-4A7905D58A11
Visible to Intel only — GUID: GUID-F52443B0-8063-4E7B-A7ED-4A7905D58A11
?sbgst
Reduces a real symmetric-definite generalized eigenproblem for banded matrices to the standard form using the factorization performed by ?pbstf.
Syntax
lapack_int LAPACKE_ssbgst (int matrix_layout, char vect, char uplo, lapack_int n, lapack_int ka, lapack_int kb, float* ab, lapack_int ldab, const float* bb, lapack_int ldbb, float* x, lapack_int ldx);
lapack_int LAPACKE_dsbgst (int matrix_layout, char vect, char uplo, lapack_int n, lapack_int ka, lapack_int kb, double* ab, lapack_int ldab, const double* bb, lapack_int ldbb, double* x, lapack_int ldx);
Include Files
- mkl.h
Description
To reduce the real symmetric-definite generalized eigenproblem A*z = λ*B*z to the standard form C*y=λ*y, where A, B and C are banded, this routine must be preceded by a call to pbstf, which computes the split Cholesky factorization of the positive-definite matrix B: B=ST*S. The split Cholesky factorization, compared with the ordinary Cholesky factorization, allows the work to be approximately halved.
This routine overwrites A with C = XT*A*X, where X = inv(S)*Q and Q is an orthogonal matrix chosen (implicitly) to preserve the bandwidth of A. The routine also has an option to allow the accumulation of X, and then, if z is an eigenvector of C, X*z is an eigenvector of the original system.
Input Parameters
- matrix_layout
-
Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).
- vect
-
Must be 'N' or 'V'.
If vect = 'N', then matrix X is not returned;
If vect = 'V', then matrix X is returned.
- uplo
-
Must be 'U' or 'L'.
If uplo = 'U', ab stores the upper triangular part of A.
If uplo = 'L', ab stores the lower triangular part of A.
- n
-
The order of the matrices A and B (n≥ 0).
- ka
-
The number of super- or sub-diagonals in A
(ka≥ 0).
- kb
-
The number of super- or sub-diagonals in B
(ka≥kb≥ 0).
- ab, bb
-
ab(size at least max(1, ldab*n) for column major layout and at least max(1, ldab*(ka + 1)) for row major layout) is an array containing either upper or lower triangular part of the symmetric matrix A (as specified by uplo) in band storage format.
bb(size at least max(1, ldbb*n) for column major layout and at least max(1, ldbb*(kb + 1)) for row major layout) is an array containing the banded split Cholesky factor of B as specified by uplo, n and kb and returned by pbstf/pbstf.
- ldab
-
The leading dimension of the array ab; must be at least ka+1 for column major layout and max(1, n) for row major layout.
- ldbb
-
The leading dimension of the array bb; must be at least kb+1 for column major layout and max(1, n) for row major layout.
- ldx
-
The leading dimension of the output array x. Constraints: if vect = 'N', then ldx≥ 1;
if vect = 'V', then ldx≥ max(1, n).
Output Parameters
- ab
-
On exit, this array is overwritten by the upper or lower triangle of C as specified by uplo.
- x
-
Array.
If vect = 'V', then x (size at least max(1, ldx*n)) contains the n-by-n matrix X = inv(S)*Q.
If vect = 'N', then x is not referenced.
Return Values
This function returns a value info.
If info=0, the execution is successful.
If info = -i, the i-th parameter had an illegal value.
Application Notes
Forming the reduced matrix C involves implicit multiplication by inv(B). When the routine is used as a step in the computation of eigenvalues and eigenvectors of the original problem, there may be a significant loss of accuracy if B is ill-conditioned with respect to inversion.
If ka and kb are much less than n then the total number of floating-point operations is approximately 6n2*kb, when vect = 'N'. Additional (3/2)n3*(kb/ka) operations are required when vect = 'V'.