gelsd#

Functions

void sgelsd(
    const INT           m,
    const INT           n,
    const INT           nrhs,
          f32* restrict A,
    const INT           lda,
          f32* restrict B,
    const INT           ldb,
          f32* restrict S,
    const f32           rcond,
          INT*          rank,
          f32* restrict work,
    const INT           lwork,
          INT* restrict iwork,
          INT*          info
);
void sgelsd(const INT m, const INT n, const INT nrhs, f32 *restrict A, const INT lda, f32 *restrict B, const INT ldb, f32 *restrict S, const f32 rcond, INT *rank, f32 *restrict work, const INT lwork, INT *restrict iwork, INT *info)#

SGELSD computes the minimum-norm solution to a real linear least squares problem: minimize 2-norm(| b - A*x |) using the singular value decomposition (SVD) of A.

A is an M-by-N matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a “bidiagonal least squares problem” (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than RCOND times the largest singular value.

Parameters

in
m

The number of rows of A. m >= 0.

in
n

The number of columns of A. n >= 0.

in
nrhs

The number of right hand sides. nrhs >= 0.

inout
A

Double precision array, dimension (lda, n). On entry, the M-by-N matrix A. On exit, A has been destroyed.

in
lda

The leading dimension of A. lda >= max(1, m).

inout
B

Double precision array, dimension (ldb, nrhs). On entry, the M-by-NRHS right hand side matrix B. On exit, B is overwritten by the N-by-NRHS solution matrix X. If m >= n and RANK = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of elements n+1:m in that column.

in
ldb

The leading dimension of B. ldb >= max(1, max(m, n)).

out
S

Double precision array, dimension (min(m, n)). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).

in
rcond

RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

out
rank

The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).

out
work

Double precision array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The dimension of work. lwork must be at least 1. The exact minimum amount of workspace needed depends on M, N and NRHS. As long as LWORK is at least 12*N + 2*N*SMLSIZ + 8*N*NLVL + N*NRHS + (SMLSIZ+1)^2, if M is greater than or equal to N or 12*M + 2*M*SMLSIZ + 8*M*NLVL + M*NRHS + (SMLSIZ+1)^2, if M is less than N, the code will execute correctly. SMLSIZ is returned by ILAENV and is equal to the maximum size of the subproblems at the bottom of the computation tree (usually about 25), and NLVL = MAX( 0, INT( LOG_2( MIN( M,N )/(SMLSIZ+1) ) ) + 1 ) For good performance, LWORK should generally be larger. If LWORK = -1, a workspace query is assumed.

out
iwork

Integer array, dimension (max(1, liwork)). liwork >= max(1, 3 * MINMN * NLVL + 11 * MINMN), where MINMN = MIN( M,N ). On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: the algorithm for computing the SVD failed to converge; if info = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.

Functions

void dgelsd(
    const INT           m,
    const INT           n,
    const INT           nrhs,
          f64* restrict A,
    const INT           lda,
          f64* restrict B,
    const INT           ldb,
          f64* restrict S,
    const f64           rcond,
          INT*          rank,
          f64* restrict work,
    const INT           lwork,
          INT* restrict iwork,
          INT*          info
);
void dgelsd(const INT m, const INT n, const INT nrhs, f64 *restrict A, const INT lda, f64 *restrict B, const INT ldb, f64 *restrict S, const f64 rcond, INT *rank, f64 *restrict work, const INT lwork, INT *restrict iwork, INT *info)#

DGELSD computes the minimum-norm solution to a real linear least squares problem: minimize 2-norm(| b - A*x |) using the singular value decomposition (SVD) of A.

A is an M-by-N matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a “bidiagonal least squares problem” (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than RCOND times the largest singular value.

Parameters

in
m

The number of rows of A. m >= 0.

in
n

The number of columns of A. n >= 0.

in
nrhs

The number of right hand sides. nrhs >= 0.

inout
A

Double precision array, dimension (lda, n). On entry, the M-by-N matrix A. On exit, A has been destroyed.

in
lda

The leading dimension of A. lda >= max(1, m).

inout
B

Double precision array, dimension (ldb, nrhs). On entry, the M-by-NRHS right hand side matrix B. On exit, B is overwritten by the N-by-NRHS solution matrix X. If m >= n and RANK = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of elements n+1:m in that column.

in
ldb

The leading dimension of B. ldb >= max(1, max(m, n)).

out
S

Double precision array, dimension (min(m, n)). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).

in
rcond

RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

out
rank

The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).

out
work

Double precision array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The dimension of work. lwork must be at least 1. The exact minimum amount of workspace needed depends on M, N and NRHS. As long as LWORK is at least 12*N + 2*N*SMLSIZ + 8*N*NLVL + N*NRHS + (SMLSIZ+1)^2, if M is greater than or equal to N or 12*M + 2*M*SMLSIZ + 8*M*NLVL + M*NRHS + (SMLSIZ+1)^2, if M is less than N, the code will execute correctly. SMLSIZ is returned by ILAENV and is equal to the maximum size of the subproblems at the bottom of the computation tree (usually about 25), and NLVL = MAX( 0, INT( LOG_2( MIN( M,N )/(SMLSIZ+1) ) ) + 1 ) For good performance, LWORK should generally be larger. If LWORK = -1, a workspace query is assumed.

out
iwork

Integer array, dimension (max(1, liwork)). liwork >= max(1, 3 * MINMN * NLVL + 11 * MINMN), where MINMN = MIN( M,N ). On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: the algorithm for computing the SVD failed to converge; if info = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.

Functions

void cgelsd(
    const INT           m,
    const INT           n,
    const INT           nrhs,
          c64* restrict A,
    const INT           lda,
          c64* restrict B,
    const INT           ldb,
          f32* restrict S,
    const f32           rcond,
          INT*          rank,
          c64* restrict work,
    const INT           lwork,
          f32* restrict rwork,
          INT* restrict iwork,
          INT*          info
);
void cgelsd(const INT m, const INT n, const INT nrhs, c64 *restrict A, const INT lda, c64 *restrict B, const INT ldb, f32 *restrict S, const f32 rcond, INT *rank, c64 *restrict work, const INT lwork, f32 *restrict rwork, INT *restrict iwork, INT *info)#

CGELSD computes the minimum-norm solution to a complex linear least squares problem: minimize 2-norm(| b - A*x |) using the singular value decomposition (SVD) of A.

A is an M-by-N matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a “bidiagonal least squares problem” (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than RCOND times the largest singular value.

Parameters

in
m

The number of rows of A. m >= 0.

in
n

The number of columns of A. n >= 0.

in
nrhs

The number of right hand sides. nrhs >= 0.

inout
A

Complex*16 array, dimension (lda, n). On entry, the M-by-N matrix A. On exit, A has been destroyed.

in
lda

The leading dimension of A. lda >= max(1, m).

inout
B

Complex*16 array, dimension (ldb, nrhs). On entry, the M-by-NRHS right hand side matrix B. On exit, B is overwritten by the N-by-NRHS solution matrix X. If m >= n and RANK = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of the modulus of elements n+1:m in that column.

in
ldb

The leading dimension of B. ldb >= max(1, max(m, n)).

out
S

Single precision array, dimension (min(m, n)). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).

in
rcond

RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

out
rank

The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).

out
work

Complex*16 array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The dimension of work. lwork must be at least 1. If LWORK = -1, a workspace query is assumed.

out
rwork

Single precision array, dimension (max(1, lrwork)). On exit, if info = 0, rwork[0] returns the minimum lrwork.

out
iwork

Integer array, dimension (max(1, liwork)). On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: the algorithm for computing the SVD failed to converge; if info = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.

Functions

void zgelsd(
    const INT            m,
    const INT            n,
    const INT            nrhs,
          c128* restrict A,
    const INT            lda,
          c128* restrict B,
    const INT            ldb,
          f64*  restrict S,
    const f64            rcond,
          INT*           rank,
          c128* restrict work,
    const INT            lwork,
          f64*  restrict rwork,
          INT*  restrict iwork,
          INT*           info
);
void zgelsd(const INT m, const INT n, const INT nrhs, c128 *restrict A, const INT lda, c128 *restrict B, const INT ldb, f64 *restrict S, const f64 rcond, INT *rank, c128 *restrict work, const INT lwork, f64 *restrict rwork, INT *restrict iwork, INT *info)#

ZGELSD computes the minimum-norm solution to a complex linear least squares problem: minimize 2-norm(| b - A*x |) using the singular value decomposition (SVD) of A.

A is an M-by-N matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a “bidiagonal least squares problem” (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than RCOND times the largest singular value.

Parameters

in
m

The number of rows of A. m >= 0.

in
n

The number of columns of A. n >= 0.

in
nrhs

The number of right hand sides. nrhs >= 0.

inout
A

Complex*16 array, dimension (lda, n). On entry, the M-by-N matrix A. On exit, A has been destroyed.

in
lda

The leading dimension of A. lda >= max(1, m).

inout
B

Complex*16 array, dimension (ldb, nrhs). On entry, the M-by-NRHS right hand side matrix B. On exit, B is overwritten by the N-by-NRHS solution matrix X. If m >= n and RANK = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of the modulus of elements n+1:m in that column.

in
ldb

The leading dimension of B. ldb >= max(1, max(m, n)).

out
S

Double precision array, dimension (min(m, n)). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).

in
rcond

RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

out
rank

The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).

out
work

Complex*16 array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The dimension of work. lwork must be at least 1. If LWORK = -1, a workspace query is assumed.

out
rwork

Double precision array, dimension (max(1, lrwork)). On exit, if info = 0, rwork[0] returns the minimum lrwork.

out
iwork

Integer array, dimension (max(1, liwork)). On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: the algorithm for computing the SVD failed to converge; if info = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.