MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
gglse: Least squares solves Ax = b subject to Bx = d (driver)

Functions

magma_int_t magma_cgglse (magma_int_t m, magma_int_t n, magma_int_t p, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex *B, magma_int_t ldb, magmaFloatComplex *c, magmaFloatComplex *d, magmaFloatComplex *x, magmaFloatComplex *work, magma_int_t lwork, magma_int_t *info)
 CGGLSE solves the linear equality-constrained least squares (LSE) problem: More...
 
magma_int_t magma_dgglse (magma_int_t m, magma_int_t n, magma_int_t p, double *A, magma_int_t lda, double *B, magma_int_t ldb, double *c, double *d, double *x, double *work, magma_int_t lwork, magma_int_t *info)
 DGGLSE solves the linear equality-constrained least squares (LSE) problem: More...
 
magma_int_t magma_sgglse (magma_int_t m, magma_int_t n, magma_int_t p, float *A, magma_int_t lda, float *B, magma_int_t ldb, float *c, float *d, float *x, float *work, magma_int_t lwork, magma_int_t *info)
 SGGLSE solves the linear equality-constrained least squares (LSE) problem: More...
 
magma_int_t magma_zgglse (magma_int_t m, magma_int_t n, magma_int_t p, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex *B, magma_int_t ldb, magmaDoubleComplex *c, magmaDoubleComplex *d, magmaDoubleComplex *x, magmaDoubleComplex *work, magma_int_t lwork, magma_int_t *info)
 ZGGLSE solves the linear equality-constrained least squares (LSE) problem: More...
 

Detailed Description

Function Documentation

magma_int_t magma_cgglse ( magma_int_t  m,
magma_int_t  n,
magma_int_t  p,
magmaFloatComplex *  A,
magma_int_t  lda,
magmaFloatComplex *  B,
magma_int_t  ldb,
magmaFloatComplex *  c,
magmaFloatComplex *  d,
magmaFloatComplex *  x,
magmaFloatComplex *  work,
magma_int_t  lwork,
magma_int_t *  info 
)

CGGLSE solves the linear equality-constrained least squares (LSE) problem:

    minimize || c - A*x ||_2   subject to   B*x = d

where A is an M-by-N matrix, B is a P-by-N matrix, c is a given M-vector, and d is a given P-vector. It is assumed that P <= N <= M+P, and

     rank(B) = P and  rank( ( A ) ) = N.
                          ( ( B ) )

These conditions ensure that the LSE problem has a unique solution, which is obtained using a GRQ factorization of the matrices B and A.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrices A and B. N >= 0.
[in]pINTEGER The number of rows of the matrix B. 0 <= P <= N <= M+P.
[in,out]ACOMPLEX array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, A is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[in,out]BCOMPLEX array, dimension (LDB,N) On entry, the P-by-N matrix B. On exit, B is destroyed.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,P).
[in,out]cCOMPLEX array, dimension (M) On entry, C contains the right hand side vector for the least squares part of the LSE problem. On exit, the residual sum of squares for the solution is given by the sum of squares of elements N-P+1 to M of vector C.
[in,out]dCOMPLEX array, dimension (P) On entry, D contains the right hand side vector for the constrained equation. On exit, D is destroyed.
[out]xCOMPLEX array, dimension (N) On exit, x is the solution of the LSE problem.
[out]work(workspace) COMPLEX array, dimension (LWORK) On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkINTEGER The dimension of the array WORK. LWORK >= max(1,M+N+P). For optimum performance LWORK >= P+min(M,N)+max(M,N)*NB, where NB is an upper bound for the optimal blocksizes for CGEQRF, CGERQF, CUNMQR and CUNMRQ.

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.

Parameters
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value.
magma_int_t magma_dgglse ( magma_int_t  m,
magma_int_t  n,
magma_int_t  p,
double *  A,
magma_int_t  lda,
double *  B,
magma_int_t  ldb,
double *  c,
double *  d,
double *  x,
double *  work,
magma_int_t  lwork,
magma_int_t *  info 
)

DGGLSE solves the linear equality-constrained least squares (LSE) problem:

    minimize || c - A*x ||_2   subject to   B*x = d

where A is an M-by-N matrix, B is a P-by-N matrix, c is a given M-vector, and d is a given P-vector. It is assumed that P <= N <= M+P, and

     rank(B) = P and  rank( ( A ) ) = N.
                          ( ( B ) )

These conditions ensure that the LSE problem has a unique solution, which is obtained using a GRQ factorization of the matrices B and A.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrices A and B. N >= 0.
[in]pINTEGER The number of rows of the matrix B. 0 <= P <= N <= M+P.
[in,out]ADOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, A is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[in,out]BDOUBLE PRECISION array, dimension (LDB,N) On entry, the P-by-N matrix B. On exit, B is destroyed.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,P).
[in,out]cDOUBLE PRECISION array, dimension (M) On entry, C contains the right hand side vector for the least squares part of the LSE problem. On exit, the residual sum of squares for the solution is given by the sum of squares of elements N-P+1 to M of vector C.
[in,out]dDOUBLE PRECISION array, dimension (P) On entry, D contains the right hand side vector for the constrained equation. On exit, D is destroyed.
[out]xDOUBLE PRECISION array, dimension (N) On exit, x is the solution of the LSE problem.
[out]work(workspace) DOUBLE PRECISION array, dimension (LWORK) On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkINTEGER The dimension of the array WORK. LWORK >= max(1,M+N+P). For optimum performance LWORK >= P+min(M,N)+max(M,N)*NB, where NB is an upper bound for the optimal blocksizes for DGEQRF, CGERQF, DORMQR and CUNMRQ.

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.

Parameters
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value.
magma_int_t magma_sgglse ( magma_int_t  m,
magma_int_t  n,
magma_int_t  p,
float *  A,
magma_int_t  lda,
float *  B,
magma_int_t  ldb,
float *  c,
float *  d,
float *  x,
float *  work,
magma_int_t  lwork,
magma_int_t *  info 
)

SGGLSE solves the linear equality-constrained least squares (LSE) problem:

    minimize || c - A*x ||_2   subject to   B*x = d

where A is an M-by-N matrix, B is a P-by-N matrix, c is a given M-vector, and d is a given P-vector. It is assumed that P <= N <= M+P, and

     rank(B) = P and  rank( ( A ) ) = N.
                          ( ( B ) )

These conditions ensure that the LSE problem has a unique solution, which is obtained using a GRQ factorization of the matrices B and A.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrices A and B. N >= 0.
[in]pINTEGER The number of rows of the matrix B. 0 <= P <= N <= M+P.
[in,out]AREAL array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, A is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[in,out]BREAL array, dimension (LDB,N) On entry, the P-by-N matrix B. On exit, B is destroyed.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,P).
[in,out]cREAL array, dimension (M) On entry, C contains the right hand side vector for the least squares part of the LSE problem. On exit, the residual sum of squares for the solution is given by the sum of squares of elements N-P+1 to M of vector C.
[in,out]dREAL array, dimension (P) On entry, D contains the right hand side vector for the constrained equation. On exit, D is destroyed.
[out]xREAL array, dimension (N) On exit, x is the solution of the LSE problem.
[out]work(workspace) REAL array, dimension (LWORK) On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkINTEGER The dimension of the array WORK. LWORK >= max(1,M+N+P). For optimum performance LWORK >= P+min(M,N)+max(M,N)*NB, where NB is an upper bound for the optimal blocksizes for SGEQRF, CGERQF, SORMQR and CUNMRQ.

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.

Parameters
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value.
magma_int_t magma_zgglse ( magma_int_t  m,
magma_int_t  n,
magma_int_t  p,
magmaDoubleComplex *  A,
magma_int_t  lda,
magmaDoubleComplex *  B,
magma_int_t  ldb,
magmaDoubleComplex *  c,
magmaDoubleComplex *  d,
magmaDoubleComplex *  x,
magmaDoubleComplex *  work,
magma_int_t  lwork,
magma_int_t *  info 
)

ZGGLSE solves the linear equality-constrained least squares (LSE) problem:

    minimize || c - A*x ||_2   subject to   B*x = d

where A is an M-by-N matrix, B is a P-by-N matrix, c is a given M-vector, and d is a given P-vector. It is assumed that P <= N <= M+P, and

     rank(B) = P and  rank( ( A ) ) = N.
                          ( ( B ) )

These conditions ensure that the LSE problem has a unique solution, which is obtained using a GRQ factorization of the matrices B and A.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrices A and B. N >= 0.
[in]pINTEGER The number of rows of the matrix B. 0 <= P <= N <= M+P.
[in,out]ACOMPLEX_16 array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, A is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[in,out]BCOMPLEX_16 array, dimension (LDB,N) On entry, the P-by-N matrix B. On exit, B is destroyed.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,P).
[in,out]cCOMPLEX_16 array, dimension (M) On entry, C contains the right hand side vector for the least squares part of the LSE problem. On exit, the residual sum of squares for the solution is given by the sum of squares of elements N-P+1 to M of vector C.
[in,out]dCOMPLEX_16 array, dimension (P) On entry, D contains the right hand side vector for the constrained equation. On exit, D is destroyed.
[out]xCOMPLEX_16 array, dimension (N) On exit, x is the solution of the LSE problem.
[out]work(workspace) COMPLEX_16 array, dimension (LWORK) On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkINTEGER The dimension of the array WORK. LWORK >= max(1,M+N+P). For optimum performance LWORK >= P+min(M,N)+max(M,N)*NB, where NB is an upper bound for the optimal blocksizes for ZGEQRF, CGERQF, ZUNMQR and CUNMRQ.

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.

Parameters
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value.