![]() |
MAGMA
1.5.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
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) |
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, double *dA, magma_int_t ldda, magma_int_t *ipiv, double *dB, magma_int_t lddb, magma_int_t *info) |
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, double *dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *dipiv, double *dB, magma_int_t lddb, double *dX, magma_int_t lddx, double *dworkd, float *dworks, magma_int_t *iter, magma_int_t *info) |
DSGESV computes the solution to a real system of linear equations A * X = B or A' * X = B where A is an 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 | ||
) |
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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in] | nrhs | INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in,out] | A | DOUBLE_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] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER 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] | B | DOUBLE_PRECISION array, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | info | INTEGER
|
magma_int_t magma_dgesv_gpu | ( | magma_int_t | n, |
magma_int_t | nrhs, | ||
double * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
double * | dB, | ||
magma_int_t | lddb, | ||
magma_int_t * | info | ||
) |
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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in] | nrhs | INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in,out] | dA | DOUBLE_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] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER 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] | dB | DOUBLE_PRECISION array on the GPU, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X. |
[in] | lddb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | info | INTEGER
|
magma_int_t magma_dsgesv_gpu | ( | magma_trans_t | trans, |
magma_int_t | n, | ||
magma_int_t | nrhs, | ||
double * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | dipiv, | ||
double * | dB, | ||
magma_int_t | lddb, | ||
double * | dX, | ||
magma_int_t | lddx, | ||
double * | dworkd, | ||
float * | dworks, | ||
magma_int_t * | iter, | ||
magma_int_t * | info | ||
) |
DSGESV computes the solution to a real system of linear equations A * X = B or A' * 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.
[in] | trans | magma_trans_t Specifies the form of the system of equations:
|
[in] | n | INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | nrhs | INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in,out] | dA | DOUBLE 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] | ldda | INTEGER The leading dimension of the array dA. ldda >= max(1,N). |
[out] | ipiv | INTEGER 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] | dipiv | INTEGER 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] | dB | DOUBLE PRECISION array on the GPU, dimension (lddb,NRHS) The N-by-NRHS right hand side matrix B. |
[in] | lddb | INTEGER The leading dimension of the array dB. lddb >= max(1,N). |
[out] | dX | DOUBLE PRECISION array on the GPU, dimension (lddx,NRHS) If info = 0, the N-by-NRHS solution matrix X. |
[in] | lddx | INTEGER 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] | iter | INTEGER
|
[out] | info | INTEGER
|