![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_cposv_batched (magma_uplo_t uplo, magma_int_t n, magma_int_t nrhs, magmaFloatComplex **dA_array, magma_int_t ldda, magmaFloatComplex **dB_array, magma_int_t lddb, magma_int_t *dinfo_array, magma_int_t batchCount, magma_queue_t queue) |
CPOSV computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices. | |
magma_int_t | magma_dposv_batched (magma_uplo_t uplo, 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 *dinfo_array, magma_int_t batchCount, magma_queue_t queue) |
DPOSV computes the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric positive definite matrix and X and B are N-by-NRHS matrices. | |
magma_int_t | magma_sposv_batched (magma_uplo_t uplo, magma_int_t n, magma_int_t nrhs, float **dA_array, magma_int_t ldda, float **dB_array, magma_int_t lddb, magma_int_t *dinfo_array, magma_int_t batchCount, magma_queue_t queue) |
SPOSV computes the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric positive definite matrix and X and B are N-by-NRHS matrices. | |
magma_int_t | magma_zposv_batched (magma_uplo_t uplo, magma_int_t n, magma_int_t nrhs, magmaDoubleComplex **dA_array, magma_int_t ldda, magmaDoubleComplex **dB_array, magma_int_t lddb, magma_int_t *dinfo_array, magma_int_t batchCount, magma_queue_t queue) |
ZPOSV computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices. | |
magma_int_t magma_cposv_batched | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magma_int_t | nrhs, | ||
magmaFloatComplex ** | dA_array, | ||
magma_int_t | ldda, | ||
magmaFloatComplex ** | dB_array, | ||
magma_int_t | lddb, | ||
magma_int_t * | dinfo_array, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue ) |
CPOSV computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices.
The Cholesky decomposition is used to factor A as A = U**H * U, if UPLO = MagmaUpper, or A = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations A * X = B.
[in] | uplo | magma_uplo_t
|
[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_array | Array of pointers, dimension (batchCount). Each is a COMPLEX array on the GPU, dimension (LDDA,N) On entry, each pointer is a Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if corresponding entry in dinfo_array = 0, each pointer is the factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H. |
[in] | ldda | INTEGER The leading dimension of each array A. LDA >= max(1,N). |
[in,out] | dB_array | Array of pointers, dimension (batchCount). Each is a COMPLEX array on the GPU, dimension (LDB,NRHS) On entry, each pointer is a right hand side matrix B. On exit, each pointer is the corresponding solution matrix X. |
[in] | lddb | INTEGER The leading dimension of each array B. LDB >= max(1,N). |
[out] | dinfo_array | Array of INTEGERs, dimension (batchCount), for corresponding matrices.
|
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_dposv_batched | ( | magma_uplo_t | uplo, |
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 * | dinfo_array, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue ) |
DPOSV computes the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric positive definite matrix and X and B are N-by-NRHS matrices.
The Cholesky decomposition is used to factor A as A = U**H * U, if UPLO = MagmaUpper, or A = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations A * X = B.
[in] | uplo | magma_uplo_t
|
[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_array | Array of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array on the GPU, dimension (LDDA,N) On entry, each pointer is a symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if corresponding entry in dinfo_array = 0, each pointer is the factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H. |
[in] | ldda | INTEGER The leading dimension of each array A. LDA >= max(1,N). |
[in,out] | dB_array | Array of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array on the GPU, dimension (LDB,NRHS) On entry, each pointer is a right hand side matrix B. On exit, each pointer is the corresponding solution matrix X. |
[in] | lddb | INTEGER The leading dimension of each array B. LDB >= max(1,N). |
[out] | dinfo_array | Array of INTEGERs, dimension (batchCount), for corresponding matrices.
|
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_sposv_batched | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magma_int_t | nrhs, | ||
float ** | dA_array, | ||
magma_int_t | ldda, | ||
float ** | dB_array, | ||
magma_int_t | lddb, | ||
magma_int_t * | dinfo_array, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue ) |
SPOSV computes the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric positive definite matrix and X and B are N-by-NRHS matrices.
The Cholesky decomposition is used to factor A as A = U**H * U, if UPLO = MagmaUpper, or A = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations A * X = B.
[in] | uplo | magma_uplo_t
|
[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_array | Array of pointers, dimension (batchCount). Each is a REAL array on the GPU, dimension (LDDA,N) On entry, each pointer is a symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if corresponding entry in dinfo_array = 0, each pointer is the factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H. |
[in] | ldda | INTEGER The leading dimension of each array A. LDA >= max(1,N). |
[in,out] | dB_array | Array of pointers, dimension (batchCount). Each is a REAL array on the GPU, dimension (LDB,NRHS) On entry, each pointer is a right hand side matrix B. On exit, each pointer is the corresponding solution matrix X. |
[in] | lddb | INTEGER The leading dimension of each array B. LDB >= max(1,N). |
[out] | dinfo_array | Array of INTEGERs, dimension (batchCount), for corresponding matrices.
|
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_zposv_batched | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magma_int_t | nrhs, | ||
magmaDoubleComplex ** | dA_array, | ||
magma_int_t | ldda, | ||
magmaDoubleComplex ** | dB_array, | ||
magma_int_t | lddb, | ||
magma_int_t * | dinfo_array, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue ) |
ZPOSV computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices.
The Cholesky decomposition is used to factor A as A = U**H * U, if UPLO = MagmaUpper, or A = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations A * X = B.
[in] | uplo | magma_uplo_t
|
[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_array | Array of pointers, dimension (batchCount). Each is a COMPLEX_16 array on the GPU, dimension (LDDA,N) On entry, each pointer is a Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if corresponding entry in dinfo_array = 0, each pointer is the factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H. |
[in] | ldda | INTEGER The leading dimension of each array A. LDA >= max(1,N). |
[in,out] | dB_array | Array of pointers, dimension (batchCount). Each is a COMPLEX_16 array on the GPU, dimension (LDB,NRHS) On entry, each pointer is a right hand side matrix B. On exit, each pointer is the corresponding solution matrix X. |
[in] | lddb | INTEGER The leading dimension of each array B. LDB >= max(1,N). |
[out] | dinfo_array | Array of INTEGERs, dimension (batchCount), for corresponding matrices.
|
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |