MAGMA  1.6.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
double precision

Functions

magma_int_t magma_dgerbt_batched (magma_bool_t gen, magma_int_t n, magma_int_t nrhs, double **dA_array, magma_int_t ldda, double **dB_array, magma_int_t lddb, double *U, double *V, magma_int_t *info, magma_int_t batchCount, magma_queue_t queue)
 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_dgerbt_gpu (magma_bool_t gen, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dB, magma_int_t lddb, double *U, double *V, 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_rbt_batched (magma_int_t n, magma_int_t nrhs, double **dA_array, magma_int_t ldda, double **dB_array, magma_int_t lddb, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU. More...
 
magma_int_t magma_dgetrf (magma_int_t m, magma_int_t n, double *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf2_mgpu (magma_int_t ngpu, magma_int_t m, magma_int_t n, magma_int_t nb, magma_int_t offset, magmaDouble_ptr d_lAT[], magma_int_t lddat, magma_int_t *ipiv, magmaDouble_ptr d_lAP[], double *W, magma_int_t ldw, magma_queue_t queues[][2], magma_int_t *info)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_batched (magma_int_t m, magma_int_t n, double **dA_array, magma_int_t ldda, magma_int_t **ipiv_array, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_m (magma_int_t ngpu, magma_int_t m, magma_int_t n, double *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info)
 DGETRF_m computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_mgpu (magma_int_t ngpu, magma_int_t m, magma_int_t n, magmaDouble_ptr d_lA[], magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_nopiv (magma_int_t m, magma_int_t n, double *A, magma_int_t lda, magma_int_t *info)
 DGETRF_NOPIV computes an LU factorization of a general M-by-N matrix A without pivoting. More...
 
magma_int_t magma_dgetrf_nopiv_batched (magma_int_t m, magma_int_t n, double **dA_array, magma_int_t lda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges. More...
 
magma_int_t magma_dgetrf_nopiv_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *info)
 DGETRF_NOPIV_GPU computes an LU factorization of a general M-by-N matrix A without any pivoting. More...
 
magma_int_t magma_dgetri_gpu (magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaDouble_ptr dwork, magma_int_t lwork, magma_int_t *info)
 DGETRI computes the inverse of a matrix using the LU factorization computed by DGETRF. More...
 
magma_int_t magma_dgetri_outofplace_batched (magma_int_t n, double **dA_array, magma_int_t ldda, magma_int_t **dipiv_array, double **dinvA_array, magma_int_t lddia, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 DGETRI computes the inverse of a matrix using the LU factorization computed by DGETRF. More...
 
magma_int_t magma_dgetrs_batched (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, double **dA_array, magma_int_t ldda, magma_int_t **dipiv_array, double **dB_array, magma_int_t lddb, magma_int_t batchCount, magma_queue_t queue)
 Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU. More...
 
magma_int_t magma_dgetrs_gpu (magma_trans_t trans, 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)
 Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU. More...
 
magma_int_t magma_dgetrs_nopiv_batched (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, double **dA_array, magma_int_t ldda, double **dB_array, magma_int_t lddb, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU. More...
 
magma_int_t magma_dgetrs_nopiv_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dB, magma_int_t lddb, magma_int_t *info)
 Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_NOPIV_GPU. More...
 
magma_int_t magma_dsgetrs_gpu (magma_trans_t trans, magma_int_t n, magma_int_t nrhs, magmaFloat_ptr dA, magma_int_t ldda, magmaInt_ptr dipiv, magmaDouble_ptr dB, magma_int_t lddb, magmaDouble_ptr dX, magma_int_t lddx, magmaFloat_ptr dSX, magma_int_t *info)
 DSGETRS solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by MAGMA_SGETRF_GPU. More...
 

Detailed Description

Function Documentation

magma_int_t magma_dgerbt_batched ( magma_bool_t  gen,
magma_int_t  n,
magma_int_t  nrhs,
double **  dA_array,
magma_int_t  ldda,
double **  dB_array,
magma_int_t  lddb,
double *  U,
double *  V,
magma_int_t *  info,
magma_int_t  batchCount,
magma_queue_t  queue 
)

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.

Random Butterfly Tranformation is applied on A and B, then the LU decomposition with no pivoting is used to factor A as A = L * U, where 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]genmagma_bool_t
  • = MagmaTrue: new matrices are generated for U and V
  • = MagmaFalse: matrices U and V given as parameter are used
[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, dimension (LDA,n). On entry, the M-by-n matrix to be factored. On exit, the factors L and U from the factorization A = L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,n).
[in,out]dBDOUBLE_PRECISION array, dimension (LDB,nrhs) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDB >= max(1,n).
[in,out]UDOUBLE_PRECISION array, dimension (2,n) Random butterfly matrix, if gen = MagmaTrue U is generated and returned as output; else we use U given as input. CPU memory
[in,out]VDOUBLE_PRECISION array, dimension (2,n) Random butterfly matrix, if gen = MagmaTrue V is generated and returned as output; else we use U given as input. CPU memory
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
magma_int_t magma_dgerbt_gpu ( magma_bool_t  gen,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dB,
magma_int_t  lddb,
double *  U,
double *  V,
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.

Random Butterfly Tranformation is applied on A and B, then the LU decomposition with no pivoting is used to factor A as A = L * U, where 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]genmagma_bool_t
  • = MagmaTrue: new matrices are generated for U and V
  • = MagmaFalse: matrices U and V given as parameter are used
[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, dimension (LDA,n). On entry, the M-by-n matrix to be factored. On exit, the factors L and U from the factorization A = L*U; the unit diagonal elements of L are not stored.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,n).
[in,out]dBDOUBLE_PRECISION array, dimension (LDB,nrhs) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER The leading dimension of the array B. LDB >= max(1,n).
[in,out]UDOUBLE_PRECISION array, dimension (2,n) Random butterfly matrix, if gen = MagmaTrue U is generated and returned as output; else we use U given as input. CPU memory
[in,out]VDOUBLE_PRECISION array, dimension (2,n) Random butterfly matrix, if gen = MagmaTrue V is generated and returned as output; else we use U given as input. CPU memory
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
magma_int_t magma_dgesv_rbt_batched ( magma_int_t  n,
magma_int_t  nrhs,
double **  dA_array,
magma_int_t  ldda,
double **  dB_array,
magma_int_t  lddb,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU.

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 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]dADOUBLE_PRECISION array on the GPU, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[in,out]dBDOUBLE_PRECISION array on the GPU, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER 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_dgetrf ( magma_int_t  m,
magma_int_t  n,
double *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
magma_int_t *  info 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

This version does not require work space on the GPU passed as input. GPU memory is allocated in the routine.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

If the current stream is NULL, this version replaces it with a new stream to overlap computation with communication.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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.
Higher performance is achieved if A is in pinned memory, e.g. allocated using magma_malloc_pinned.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf2_mgpu ( magma_int_t  ngpu,
magma_int_t  m,
magma_int_t  n,
magma_int_t  nb,
magma_int_t  offset,
magmaDouble_ptr  d_lAT[],
magma_int_t  lddat,
magma_int_t *  ipiv,
magmaDouble_ptr  d_lAP[],
double *  W,
magma_int_t  ldw,
magma_queue_t  queues[][2],
magma_int_t *  info 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm. Use two buffer to send panels.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 0.
[in,out]d_lATDOUBLE_PRECISION array of pointers on the GPU, dimension (ngpu). On entry, the M-by-N matrix A distributed over GPUs (d_lAT[d] points to the local matrix on d-th GPU). It uses a 1D block column cyclic format (with the block size nb), and each local matrix is stored by row. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in]lddatINTEGER The leading dimension of the array d_lAT[d]. LDDA >= max(1,M).
[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).
(workspace)on device d_lAP DOUBLE_PRECISION array of pointers on the GPU, dimension (ngpu). d_lAP[d] is the workspace on d-th GPU. Each local workspace must be of size (3+ngpu)*nb*maxm, where maxm is m rounded up to a multiple of 32 and nb is the block size.
(workspace)W DOUBLE_PRECISION array, dimension (ngpu*nb*maxm). It is used to store panel on CPU.
[in]ldwINTEGER The leading dimension of the workspace w.
[in]queuesmagma_queue_t queues[d] points to the streams for the d-th GPU to execute in. Each GPU require two streams.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_batched ( magma_int_t  m,
magma_int_t  n,
double **  dA_array,
magma_int_t  ldda,
magma_int_t **  ipiv_array,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

If the current stream is NULL, this version replaces it with a new stream to overlap computation with communication.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_gpu ( magma_int_t  m,
magma_int_t  n,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magma_int_t *  info 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

If the current stream is NULL, this version replaces it with a new stream to overlap computation with communication.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_m ( magma_int_t  ngpu,
magma_int_t  m,
magma_int_t  n,
double *  A,
magma_int_t  lda,
magma_int_t *  ipiv,
magma_int_t *  info 
)

DGETRF_m computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

This version does not require work space on the GPU passed as input. GPU memory is allocated in the routine. The matrix may exceed the GPU memory.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

Note: The factorization of big panel is done calling multiple-gpu-interface. Pivots are applied on GPU within the big panel.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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.
Higher performance is achieved if A is in pinned memory, e.g. allocated using magma_malloc_pinned.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_mgpu ( magma_int_t  ngpu,
magma_int_t  m,
magma_int_t  n,
magmaDouble_ptr  d_lA[],
magma_int_t  ldda,
magma_int_t *  ipiv,
magma_int_t *  info 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 0.
[in,out]d_lADOUBLE_PRECISION array of pointers on the GPU, dimension (ngpu). On entry, the M-by-N matrix A distributed over GPUs (d_lA[d] points to the local matrix on d-th GPU). It uses 1D block column cyclic format with the block size of nb, and each local matrix is stored by column. 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 d_lA. LDDA >= max(1,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_nopiv ( magma_int_t  m,
magma_int_t  n,
double *  A,
magma_int_t  lda,
magma_int_t *  info 
)

DGETRF_NOPIV computes an LU factorization of a general M-by-N matrix A without pivoting.

The factorization has the form A = L * U where L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_nopiv_batched ( magma_int_t  m,
magma_int_t  n,
double **  dA_array,
magma_int_t  lda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

DGETRF computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

If the current stream is NULL, this version replaces it with a new stream to overlap computation with communication.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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,M).
[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).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetrf_nopiv_gpu ( magma_int_t  m,
magma_int_t  n,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  info 
)

DGETRF_NOPIV_GPU computes an LU factorization of a general M-by-N matrix A without any pivoting.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking Level 3 BLAS version of the algorithm.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 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,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value or another error occured, such as memory allocation failed.
  • > 0: if INFO = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.
magma_int_t magma_dgetri_gpu ( magma_int_t  n,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magma_int_t *  ipiv,
magmaDouble_ptr  dwork,
magma_int_t  lwork,
magma_int_t *  info 
)

DGETRI computes the inverse of a matrix using the LU factorization computed by DGETRF.

This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).

Note that it is generally both faster and more accurate to use DGESV, or DGETRF and DGETRS, to solve the system AX = B, rather than inverting the matrix and multiplying to form X = inv(A)*B. Only in special instances should an explicit inverse be computed with this routine.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dADOUBLE_PRECISION array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[out]dwork(workspace) DOUBLE_PRECISION array on the GPU, dimension (MAX(1,LWORK))
[in]lworkINTEGER The dimension of the array DWORK. LWORK >= N*NB, where NB is the optimal blocksize returned by magma_get_dgetri_nb(n).
Unlike LAPACK, this version does not currently support a workspace query, because the workspace is on the GPU.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its cannot be computed.
magma_int_t magma_dgetri_outofplace_batched ( magma_int_t  n,
double **  dA_array,
magma_int_t  ldda,
magma_int_t **  dipiv_array,
double **  dinvA_array,
magma_int_t  lddia,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

DGETRI computes the inverse of a matrix using the LU factorization computed by DGETRF.

This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).

Note that it is generally both faster and more accurate to use DGESV, or DGETRF and DGETRS, to solve the system AX = B, rather than inverting the matrix and multiplying to form X = inv(A)*B. Only in special instances should an explicit inverse be computed with this routine.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dADOUBLE_PRECISION array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[out]dwork(workspace) DOUBLE_PRECISION array on the GPU, dimension (MAX(1,LWORK))
[in]lworkINTEGER The dimension of the array DWORK. LWORK >= N*NB, where NB is the optimal blocksize returned by magma_get_dgetri_nb(n).
Unlike LAPACK, this version does not currently support a workspace query, because the workspace is on the GPU.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its cannot be computed.
magma_int_t magma_dgetrs_batched ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
double **  dA_array,
magma_int_t  ldda,
magma_int_t **  dipiv_array,
double **  dB_array,
magma_int_t  lddb,
magma_int_t  batchCount,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU.

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 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]dADOUBLE_PRECISION array on the GPU, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[in,out]dBDOUBLE_PRECISION array on the GPU, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER 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_dgetrs_gpu ( magma_trans_t  trans,
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 
)

Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU.

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 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]dADOUBLE_PRECISION array on the GPU, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[in,out]dBDOUBLE_PRECISION array on the GPU, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER 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_dgetrs_nopiv_batched ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
double **  dA_array,
magma_int_t  ldda,
double **  dB_array,
magma_int_t  lddb,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_GPU.

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 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]dADOUBLE_PRECISION array on the GPU, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i).
[in,out]dBDOUBLE_PRECISION array on the GPU, dimension (LDB,NRHS) On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]lddbINTEGER 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_dgetrs_nopiv_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magma_int_t *  info 
)

Solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by DGETRF_NOPIV_GPU.

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 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]dADOUBLE_PRECISION array on the GPU, dimension (LDA,N) The factors L and U from the factorization A = P*L*U as computed by DGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,N).

param[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.

Parameters
[in]lddbINTEGER 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_dsgetrs_gpu ( magma_trans_t  trans,
magma_int_t  n,
magma_int_t  nrhs,
magmaFloat_ptr  dA,
magma_int_t  ldda,
magmaInt_ptr  dipiv,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magmaDouble_ptr  dX,
magma_int_t  lddx,
magmaFloat_ptr  dSX,
magma_int_t *  info 
)

DSGETRS solves a system of linear equations A * X = B, A**T * X = B, or A**H * X = B with a general N-by-N matrix A using the LU factorization computed by MAGMA_SGETRF_GPU.

B and X are in DOUBLE PRECISION, and A is in SINGLE PRECISION. This routine is used in the mixed precision iterative solver magma_dsgesv.

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 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]dASINGLE PRECISION array on the GPU, dimension (LDDA,N) The factors L and U from the factorization A = P*L*U as computed by CGETRF_GPU.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,N).
[in]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 from DGETRF, where interchanges are applied one-after-another.
[in]dBDOUBLE PRECISION array on the GPU, dimension (LDDB,NRHS) On entry, the right hand side matrix B.
[in]lddbINTEGER The leading dimension of the arrays X and B. LDDB >= max(1,N).
[out]dXDOUBLE PRECISION array on the GPU, dimension (LDDX, NRHS) On exit, the solution matrix dX.
[in]lddxINTEGER The leading dimension of the array dX, LDDX >= max(1,N).
dSX(workspace) SINGLE PRECISION array on the GPU used as workspace, dimension (N, NRHS)
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value