MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
gesv: Solves Ax = b using LU factorization (driver)

Functions

magma_int_t magma_cgesv (magma_int_t n, magma_int_t nrhs, magmaFloatComplex *A, magma_int_t lda, magma_int_t *ipiv, magmaFloatComplex *B, magma_int_t ldb, magma_int_t *info)
 CGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_cgesv_gpu (magma_int_t n, magma_int_t nrhs, magmaFloatComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaFloatComplex_ptr dB, magma_int_t lddb, magma_int_t *info)
 CGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_dgesv (magma_int_t n, magma_int_t nrhs, double *A, magma_int_t lda, magma_int_t *ipiv, double *B, magma_int_t ldb, magma_int_t *info)
 DGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_dgesv_gpu (magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaDouble_ptr dB, magma_int_t lddb, magma_int_t *info)
 DGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_dsgesv_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaInt_ptr dipiv, magmaDouble_ptr dB, magma_int_t lddb, magmaDouble_ptr dX, magma_int_t lddx, magmaDouble_ptr dworkd, magmaFloat_ptr dworks, magma_int_t *iter, magma_int_t *info)
 DSGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_dxgesv_gmres_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaInt_ptr dipiv, magmaDouble_ptr dB, magma_int_t lddb, magmaDouble_ptr dX, magma_int_t lddx, magmaDouble_ptr dworkd, magmaFloat_ptr dworks, magma_refinement_t facto_type, magma_refinement_t solver_type, magma_int_t *iter, magma_int_t *info, real_Double_t *facto_time)
 DSGESV or DHGESV expert interface. More...
 
magma_int_t magma_dsgesv_iteref_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaInt_ptr dipiv, magmaDouble_ptr dB, magma_int_t lddb, magmaDouble_ptr dX, magma_int_t lddx, magmaDouble_ptr dworkd, magmaFloat_ptr dworks, magma_int_t *iter, magma_int_t *info)
 DSGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_dhgesv_iteref_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaInt_ptr dipiv, magmaDouble_ptr dB, magma_int_t lddb, magmaDouble_ptr dX, magma_int_t lddx, magmaDouble_ptr dworkd, magmaFloat_ptr dworks, magma_int_t *iter, magma_int_t *info)
 DHGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_sgesv (magma_int_t n, magma_int_t nrhs, float *A, magma_int_t lda, magma_int_t *ipiv, float *B, magma_int_t ldb, magma_int_t *info)
 SGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_sgesv_gpu (magma_int_t n, magma_int_t nrhs, magmaFloat_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaFloat_ptr dB, magma_int_t lddb, magma_int_t *info)
 SGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_zcgesv_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaInt_ptr dipiv, magmaDoubleComplex_ptr dB, magma_int_t lddb, magmaDoubleComplex_ptr dX, magma_int_t lddx, magmaDoubleComplex_ptr dworkd, magmaFloatComplex_ptr dworks, magma_int_t *iter, magma_int_t *info)
 ZCGESV computes the solution to a complex system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_zgesv (magma_int_t n, magma_int_t nrhs, magmaDoubleComplex *A, magma_int_t lda, magma_int_t *ipiv, magmaDoubleComplex *B, magma_int_t ldb, magma_int_t *info)
 ZGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 
magma_int_t magma_zgesv_gpu (magma_int_t n, magma_int_t nrhs, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_int_t *info)
 ZGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices. More...
 

Detailed Description

Function Documentation

magma_int_t magma_cgesv ( magma_int_t  n,
magma_int_t  nrhs,
magmaFloatComplex *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
magmaFloatComplex *  B,
magma_int_t  ldb,
magma_int_t *  info 
)

CGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]ACOMPLEX array, dimension (LDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]BCOMPLEX array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_cgesv_gpu ( magma_int_t  n,
magma_int_t  nrhs,
magmaFloatComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaFloatComplex_ptr  dB,
magma_int_t  lddb,
magma_int_t *  info 
)

CGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dACOMPLEX array on the GPU, dimension (LDDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]dBCOMPLEX array on the GPU, dimension (LDDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_dgesv ( magma_int_t  n,
magma_int_t  nrhs,
double *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
double *  B,
magma_int_t  ldb,
magma_int_t *  info 
)

DGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]ADOUBLE PRECISION array, dimension (LDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]BDOUBLE PRECISION array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_dgesv_gpu ( magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magma_int_t *  info 
)

DGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dADOUBLE PRECISION array on the GPU, dimension (LDDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]dBDOUBLE PRECISION array on the GPU, dimension (LDDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_dsgesv_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaInt_ptr  dipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magmaDouble_ptr  dX,
magma_int_t  lddx,
magmaDouble_ptr  dworkd,
magmaFloat_ptr  dworks,
magma_int_t *  iter,
magma_int_t *  info 
)

DSGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.

DSGESV first attempts to factorize the matrix in real SINGLE PRECISION and use this factorization within an iterative refinement procedure to produce a solution with real DOUBLE PRECISION norm-wise backward error quality (see below). If the approach fails the method switches to a real DOUBLE PRECISION factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio real SINGLE PRECISION performance over real DOUBLE PRECISION performance is too small. A reasonable strategy should take the number of right-hand sides and the size of the matrix into account. This might be done with a call to ILAENV in the future. Up to now, we always try iterative refinement.

The iterative refinement process is stopped if ITER > ITERMAX or for all the RHS we have: RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX where o ITER is the number of the current iteration in the iterative refinement process o RNRM is the infinity-norm of the residual o XNRM is the infinity-norm of the solution o ANRM is the infinity-operator-norm of the matrix A o EPS is the machine epsilon returned by DLAMCH('Epsilon') The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00 respectively.

Parameters
[in]transmagma_trans_t Specifies the form of the system of equations:
  • = MagmaNoTrans: A * X = B (No transpose)
  • = MagmaTrans: A**T * X = B (Transpose)
  • = MagmaConjTrans: A**H * X = B (Conjugate transpose)
[in]nINTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dADOUBLE PRECISION array on the GPU, dimension (ldda,N) On entry, the N-by-N coefficient matrix A. On exit, if iterative refinement has been successfully used (info.EQ.0 and ITER.GE.0, see description below), A is unchanged. If double precision factorization has been used (info.EQ.0 and ITER.LT.0, see description below), then the array dA contains the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array dA. ldda >= max(1,N).
[out]ipivINTEGER array, dimension (N) The pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i). Corresponds either to the single precision factorization (if info.EQ.0 and ITER.GE.0) or the double precision factorization (if info.EQ.0 and ITER.LT.0).
[out]dipivINTEGER array on the GPU, dimension (N) The pivot indices; for 1 <= i <= N, after permuting, row i of the matrix was moved to row dIPIV(i). Note this is different than IPIV, where interchanges are applied one-after-another.
[in]dBDOUBLE PRECISION array on the GPU, dimension (lddb,NRHS) The N-by-NRHS right hand side matrix B.
[in]lddbINTEGER The leading dimension of the array dB. lddb >= max(1,N).
[out]dXDOUBLE PRECISION array on the GPU, dimension (lddx,NRHS) If info = 0, the N-by-NRHS solution matrix X.
[in]lddxINTEGER The leading dimension of the array dX. lddx >= max(1,N).
dworkd(workspace) DOUBLE PRECISION array on the GPU, dimension (N*NRHS) This array is used to hold the residual vectors.
dworks(workspace) SINGLE PRECISION array on the GPU, dimension (N*(N+NRHS)) This array is used to store the real single precision matrix and the right-hand sides or solutions in single precision.
[out]iterINTEGER
  • < 0: iterative refinement has failed, double precision factorization has been performed
    • -1 : the routine fell back to full precision for implementation- or machine-specific reasons
    • -2 : narrowing the precision induced an overflow, the routine fell back to full precision
    • -3 : failure of SGETRF
    • -31: stop the iterative refinement after the 30th iteration
  • > 0: iterative refinement has been successfully used. Returns the number of iterations
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if info = -i, the i-th argument had an illegal value
  • > 0: if info = i, U(i,i) computed in DOUBLE PRECISION is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed.
magma_int_t magma_dxgesv_gmres_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaInt_ptr  dipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magmaDouble_ptr  dX,
magma_int_t  lddx,
magmaDouble_ptr  dworkd,
magmaFloat_ptr  dworks,
magma_refinement_t  facto_type,
magma_refinement_t  solver_type,
magma_int_t *  iter,
magma_int_t *  info,
real_Double_t *  facto_time 
)

DSGESV or DHGESV expert interface.

It computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. the accomodate the Single Precision DSGESV and the Half precision dhgesv API. precision and iterative refinement solver are specified by facto_type, solver_type. For other API parameter please refer to the corresponding dsgesv or dhgesv.

Parameters
[in]facto_typemagma_refinement_t Specify the mixed precision factorization algorithm. Magma_PREC_SS for FP32 Magma_PREC_SHT for FP16 Tensor Cores More details will be released soon.
[in]solver_typemagma_refinement_t Specify the iterative refinement technique to be used. classical IR or GMRES etc. More details will be released soon.

More details can be found in Azzam Haidar, Stanimire Tomov, Jack Dongarra, and Nicholas J. Higham. 2018. Harnessing GPU tensor cores for fast FP16 arithmetic to speed up mixed-precision iterative refinement solvers. In Proceedings of the International Conference for High Performance Computing, Networking, Storage, and Analysis (SC '18). IEEE Press, Piscataway, NJ, USA, Article 47, 11 pages.

magma_int_t magma_dsgesv_iteref_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaInt_ptr  dipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magmaDouble_ptr  dX,
magma_int_t  lddx,
magmaDouble_ptr  dworkd,
magmaFloat_ptr  dworks,
magma_int_t *  iter,
magma_int_t *  info 
)

DSGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.

DSGESV first attempts to factorize the matrix in real SINGLE PRECISION and use this factorization within an iterative refinement procedure to produce a solution with real DOUBLE PRECISION norm-wise backward error quality (see below). If the approach fails the method switches to a real DOUBLE PRECISION factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio real SINGLE PRECISION performance over real DOUBLE PRECISION performance is too small. A reasonable strategy should take the number of right-hand sides and the size of the matrix into account. This might be done with a call to ILAENV in the future. Up to now, we always try iterative refinement.

The iterative refinement process is stopped if ITER > solver_outer_itermax or for all the RHS we have: RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX where o ITER is the number of the current iteration in the iterative refinement process o RNRM is the infinity-norm of the residual o XNRM is the infinity-norm of the solution o ANRM is the infinity-operator-norm of the matrix A o EPS is the machine epsilon returned by DLAMCH('Epsilon') The value solver_outer_itermax and BWDMAX are fixed to 30 and 1.0D+00 respectively.

Parameters
[in]transmagma_trans_t Specifies the form of the system of equations:
  • = MagmaNoTrans: A * X = B (No transpose)
  • = MagmaTrans: A**T * X = B (Transpose)
  • = MagmaConjTrans: A**H * X = B (Conjugate transpose)
[in]nINTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dADOUBLE PRECISION array on the GPU, dimension (ldda,N) On entry, the N-by-N coefficient matrix A. On exit, if iterative refinement has been successfully used (info.EQ.0 and ITER.GE.0, see description below), A is unchanged. If double precision factorization has been used (info.EQ.0 and ITER.LT.0, see description below), then the array dA contains the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array dA. ldda >= max(1,N).
[out]ipivINTEGER array, dimension (N) The pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i). Corresponds either to the single precision factorization (if info.EQ.0 and ITER.GE.0) or the double precision factorization (if info.EQ.0 and ITER.LT.0).
[out]dipivINTEGER array on the GPU, dimension (N) The pivot indices; for 1 <= i <= N, after permuting, row i of the matrix was moved to row dIPIV(i). Note this is different than IPIV, where interchanges are applied one-after-another.
[in]dBDOUBLE PRECISION array on the GPU, dimension (lddb,NRHS) The N-by-NRHS right hand side matrix B.
[in]lddbINTEGER The leading dimension of the array dB. lddb >= max(1,N).
[out]dXDOUBLE PRECISION array on the GPU, dimension (lddx,NRHS) If info = 0, the N-by-NRHS solution matrix X.
[in]lddxINTEGER The leading dimension of the array dX. lddx >= max(1,N).
dworkd(workspace) DOUBLE PRECISION array on the GPU, dimension (N*NRHS) This array is used to hold the residual vectors.
dworks(workspace) SINGLE PRECISION array on the GPU, dimension (N*(N+NRHS)) This array is used to store the real single precision matrix and the right-hand sides or solutions in single precision.
[out]iterINTEGER
  • < 0: iterative refinement has failed, double precision factorization has been performed
    • -1 : the routine fell back to full precision for implementation- or machine-specific reasons
    • -2 : narrowing the precision induced an overflow, the routine fell back to full precision
    • -3 : failure of SGETRF
    • -31: stop the iterative refinement after the 30th iteration
  • > 0: iterative refinement has been successfully used. Returns the number of iterations
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if info = -i, the i-th argument had an illegal value
  • > 0: if info = i, U(i,i) computed in DOUBLE PRECISION is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed.
magma_int_t magma_dhgesv_iteref_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaInt_ptr  dipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magmaDouble_ptr  dX,
magma_int_t  lddx,
magmaDouble_ptr  dworkd,
magmaFloat_ptr  dworks,
magma_int_t *  iter,
magma_int_t *  info 
)

DHGESV computes the solution to a real system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.

DHGESV first attempts to factorize the matrix in real HALF PRECISION and use this factorization within an iterative refinement procedure to produce a solution with real DOUBLE PRECISION norm-wise backward error quality (see below). If the approach fails the method switches to a real DOUBLE PRECISION factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio real HALF PRECISION performance over real DOUBLE PRECISION performance is too small. A reasonable strategy should take the number of right-hand sides and the size of the matrix into account. This might be done with a call to ILAENV in the future. Up to now, we always try iterative refinement.

The iterative refinement process is stopped if ITER > solver_outer_itermax or for all the RHS we have: RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX where o ITER is the number of the current iteration in the iterative refinement process o RNRM is the infinity-norm of the residual o XNRM is the infinity-norm of the solution o ANRM is the infinity-operator-norm of the matrix A o EPS is the machine epsilon returned by DLAMCH('Epsilon') The value solver_outer_itermax and BWDMAX are fixed to 30 and 1.0D+00 respectively.

Parameters
[in]transmagma_trans_t Specifies the form of the system of equations:
  • = MagmaNoTrans: A * X = B (No transpose)
  • = MagmaTrans: A**T * X = B (Transpose)
  • = MagmaConjTrans: A**H * X = B (Conjugate transpose)
[in]nINTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dADOUBLE PRECISION array on the GPU, dimension (ldda,N) On entry, the N-by-N coefficient matrix A. On exit, if iterative refinement has been successfully used (info.EQ.0 and ITER.GE.0, see description below), A is unchanged. If double precision factorization has been used (info.EQ.0 and ITER.LT.0, see description below), then the array dA contains the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array dA. ldda >= max(1,N).
[out]ipivINTEGER array, dimension (N) The pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i). Corresponds either to the single precision factorization (if info.EQ.0 and ITER.GE.0) or the double precision factorization (if info.EQ.0 and ITER.LT.0).
[out]dipivINTEGER array on the GPU, dimension (N) The pivot indices; for 1 <= i <= N, after permuting, row i of the matrix was moved to row dIPIV(i). Note this is different than IPIV, where interchanges are applied one-after-another.
[in]dBDOUBLE PRECISION array on the GPU, dimension (lddb,NRHS) The N-by-NRHS right hand side matrix B.
[in]lddbINTEGER The leading dimension of the array dB. lddb >= max(1,N).
[out]dXDOUBLE PRECISION array on the GPU, dimension (lddx,NRHS) If info = 0, the N-by-NRHS solution matrix X.
[in]lddxINTEGER The leading dimension of the array dX. lddx >= max(1,N).
dworkd(workspace) DOUBLE PRECISION array on the GPU, dimension (N*NRHS) This array is used to hold the residual vectors.
dworks(workspace) SINGLE PRECISION array on the GPU, dimension (N*(N+NRHS)) This array is used to store the real single precision matrix and the right-hand sides or solutions in single precision.
[out]iterINTEGER
  • < 0: iterative refinement has failed, double precision factorization has been performed
    • -1 : the routine fell back to full precision for implementation- or machine-specific reasons
    • -2 : narrowing the precision induced an overflow, the routine fell back to full precision
    • -3 : failure of SGETRF
    • -31: stop the iterative refinement after the 30th iteration
  • > 0: iterative refinement has been successfully used. Returns the number of iterations
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if info = -i, the i-th argument had an illegal value
  • > 0: if info = i, U(i,i) computed in DOUBLE PRECISION is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed.
magma_int_t magma_sgesv ( magma_int_t  n,
magma_int_t  nrhs,
float *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
float *  B,
magma_int_t  ldb,
magma_int_t *  info 
)

SGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]AREAL array, dimension (LDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]BREAL array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_sgesv_gpu ( magma_int_t  n,
magma_int_t  nrhs,
magmaFloat_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaFloat_ptr  dB,
magma_int_t  lddb,
magma_int_t *  info 
)

SGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dAREAL array on the GPU, dimension (LDDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]dBREAL array on the GPU, dimension (LDDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_zcgesv_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaInt_ptr  dipiv,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magmaDoubleComplex_ptr  dX,
magma_int_t  lddx,
magmaDoubleComplex_ptr  dworkd,
magmaFloatComplex_ptr  dworks,
magma_int_t *  iter,
magma_int_t *  info 
)

ZCGESV computes the solution to a complex system of linear equations A * X = B, A**T * X = B, or A**H * X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.

ZCGESV first attempts to factorize the matrix in complex SINGLE PRECISION and use this factorization within an iterative refinement procedure to produce a solution with complex DOUBLE PRECISION norm-wise backward error quality (see below). If the approach fails the method switches to a complex DOUBLE PRECISION factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio complex SINGLE PRECISION performance over complex DOUBLE PRECISION performance is too small. A reasonable strategy should take the number of right-hand sides and the size of the matrix into account. This might be done with a call to ILAENV in the future. Up to now, we always try iterative refinement.

The iterative refinement process is stopped if ITER > ITERMAX or for all the RHS we have: RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX where o ITER is the number of the current iteration in the iterative refinement process o RNRM is the infinity-norm of the residual o XNRM is the infinity-norm of the solution o ANRM is the infinity-operator-norm of the matrix A o EPS is the machine epsilon returned by DLAMCH('Epsilon') The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00 respectively.

Parameters
[in]transmagma_trans_t Specifies the form of the system of equations:
  • = MagmaNoTrans: A * X = B (No transpose)
  • = MagmaTrans: A**T * X = B (Transpose)
  • = MagmaConjTrans: A**H * X = B (Conjugate transpose)
[in]nINTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dACOMPLEX_16 array on the GPU, dimension (ldda,N) On entry, the N-by-N coefficient matrix A. On exit, if iterative refinement has been successfully used (info.EQ.0 and ITER.GE.0, see description below), A is unchanged. If double precision factorization has been used (info.EQ.0 and ITER.LT.0, see description below), then the array dA contains the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array dA. ldda >= max(1,N).
[out]ipivINTEGER array, dimension (N) The pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i). Corresponds either to the single precision factorization (if info.EQ.0 and ITER.GE.0) or the double precision factorization (if info.EQ.0 and ITER.LT.0).
[out]dipivINTEGER array on the GPU, dimension (N) The pivot indices; for 1 <= i <= N, after permuting, row i of the matrix was moved to row dIPIV(i). Note this is different than IPIV, where interchanges are applied one-after-another.
[in]dBCOMPLEX_16 array on the GPU, dimension (lddb,NRHS) The N-by-NRHS right hand side matrix B.
[in]lddbINTEGER The leading dimension of the array dB. lddb >= max(1,N).
[out]dXCOMPLEX_16 array on the GPU, dimension (lddx,NRHS) If info = 0, the N-by-NRHS solution matrix X.
[in]lddxINTEGER The leading dimension of the array dX. lddx >= max(1,N).
dworkd(workspace) COMPLEX_16 array on the GPU, dimension (N*NRHS) This array is used to hold the residual vectors.
dworks(workspace) COMPLEX array on the GPU, dimension (N*(N+NRHS)) This array is used to store the complex single precision matrix and the right-hand sides or solutions in single precision.
[out]iterINTEGER
  • < 0: iterative refinement has failed, double precision factorization has been performed
    • -1 : the routine fell back to full precision for implementation- or machine-specific reasons
    • -2 : narrowing the precision induced an overflow, the routine fell back to full precision
    • -3 : failure of SGETRF
    • -31: stop the iterative refinement after the 30th iteration
  • > 0: iterative refinement has been successfully used. Returns the number of iterations
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if info = -i, the i-th argument had an illegal value
  • > 0: if info = i, U(i,i) computed in DOUBLE PRECISION is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed.
magma_int_t magma_zgesv ( magma_int_t  n,
magma_int_t  nrhs,
magmaDoubleComplex *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
magmaDoubleComplex *  B,
magma_int_t  ldb,
magma_int_t *  info 
)

ZGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]ACOMPLEX_16 array, dimension (LDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]BCOMPLEX_16 array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]ldbINTEGER The leading dimension of the array B. LDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
magma_int_t magma_zgesv_gpu ( magma_int_t  n,
magma_int_t  nrhs,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_int_t *  info 
)

ZGESV solves a system of linear equations A * X = B where A is a general N-by-N matrix and X and B are N-by-NRHS matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as A = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]nrhsINTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]dACOMPLEX_16 array on the GPU, dimension (LDDA,N). On entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[out]ipivINTEGER array, dimension (min(M,N)) The pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[in,out]dBCOMPLEX_16 array on the GPU, dimension (LDDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDDB >= max(1,N).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value