MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages

Functions

magma_int_t magma_cpotrf_batched (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex **dA_array, magma_int_t ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 CPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA. More...
 
magma_int_t magma_cpotrf_vbatched (magma_uplo_t uplo, magma_int_t *n, magmaFloatComplex **dA_array, magma_int_t *ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 CPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA. More...
 
magma_int_t magma_dpotrf_batched (magma_uplo_t uplo, magma_int_t n, double **dA_array, magma_int_t ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 DPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA. More...
 
magma_int_t magma_dpotrf_vbatched (magma_uplo_t uplo, magma_int_t *n, double **dA_array, magma_int_t *ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 DPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA. More...
 
magma_int_t magma_spotrf_batched (magma_uplo_t uplo, magma_int_t n, float **dA_array, magma_int_t ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 SPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA. More...
 
magma_int_t magma_spotrf_vbatched (magma_uplo_t uplo, magma_int_t *n, float **dA_array, magma_int_t *ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 SPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA. More...
 
magma_int_t magma_zpotrf_batched (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex **dA_array, magma_int_t ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 ZPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA. More...
 
magma_int_t magma_zpotrf_vbatched (magma_uplo_t uplo, magma_int_t *n, magmaDoubleComplex **dA_array, magma_int_t *ldda, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue)
 ZPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA. More...
 

Detailed Description

Function Documentation

magma_int_t magma_cpotrf_batched ( magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex **  dA_array,
magma_int_t  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

CPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the fixed size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER The order of the matrix dA. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array on the GPU, dimension (LDDA,N) On entry, each pointer is a Hermitian matrix dA. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of dA contains the upper triangular part of the matrix dA, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of dA contains the lower triangular part of the matrix dA, and the strictly upper triangular part of dA is not referenced.
On exit, if corresponding entry in info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER The leading dimension of each array dA. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_cpotrf_vbatched ( magma_uplo_t  uplo,
magma_int_t *  n,
magmaFloatComplex **  dA_array,
magma_int_t *  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

CPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the variable size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER array, dimension(batchCount + 1). Each element N specifies the order of each matrix A. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array A 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 A 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 info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER array, dimension(batchCount + 1). Each element LDDA specifies the leading dimension of each array A. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dpotrf_batched ( magma_uplo_t  uplo,
magma_int_t  n,
double **  dA_array,
magma_int_t  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

DPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the fixed size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER The order of the matrix dA. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array on the GPU, dimension (LDDA,N) On entry, each pointer is a symmetric matrix dA. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of dA contains the upper triangular part of the matrix dA, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of dA contains the lower triangular part of the matrix dA, and the strictly upper triangular part of dA is not referenced.
On exit, if corresponding entry in info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER The leading dimension of each array dA. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dpotrf_vbatched ( magma_uplo_t  uplo,
magma_int_t *  n,
double **  dA_array,
magma_int_t *  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

DPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the variable size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER array, dimension(batchCount + 1). Each element N specifies the order of each matrix A. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array A 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 A 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 info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER array, dimension(batchCount + 1). Each element LDDA specifies the leading dimension of each array A. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_spotrf_batched ( magma_uplo_t  uplo,
magma_int_t  n,
float **  dA_array,
magma_int_t  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

SPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the fixed size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER The order of the matrix dA. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a REAL array on the GPU, dimension (LDDA,N) On entry, each pointer is a symmetric matrix dA. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of dA contains the upper triangular part of the matrix dA, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of dA contains the lower triangular part of the matrix dA, and the strictly upper triangular part of dA is not referenced.
On exit, if corresponding entry in info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER The leading dimension of each array dA. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_spotrf_vbatched ( magma_uplo_t  uplo,
magma_int_t *  n,
float **  dA_array,
magma_int_t *  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

SPOTRF computes the Cholesky factorization of a real symmetric positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the variable size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER array, dimension(batchCount + 1). Each element N specifies the order of each matrix A. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a REAL array A 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 A 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 info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER array, dimension(batchCount + 1). Each element LDDA specifies the leading dimension of each array A. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpotrf_batched ( magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the fixed size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER The order of the matrix dA. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array on the GPU, dimension (LDDA,N) On entry, each pointer is a Hermitian matrix dA. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of dA contains the upper triangular part of the matrix dA, and the strictly lower triangular part of dA is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of dA contains the lower triangular part of the matrix dA, and the strictly upper triangular part of dA is not referenced.
On exit, if corresponding entry in info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER The leading dimension of each array dA. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_zpotrf_vbatched ( magma_uplo_t  uplo,
magma_int_t *  n,
magmaDoubleComplex **  dA_array,
magma_int_t *  ldda,
magma_int_t *  info_array,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZPOTRF computes the Cholesky factorization of a complex Hermitian positive definite matrix dA.

The factorization has the form dA = U**H * U, if UPLO = MagmaUpper, or dA = L * L**H, if UPLO = MagmaLower, where U is an upper triangular matrix and L is lower triangular.

This is the block version of the algorithm, calling Level 3 BLAS. This is the variable size batched version of the operation.

Parameters
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of dA is stored;
  • = MagmaLower: Lower triangle of dA is stored. Only MagmaLower is supported.
[in]nINTEGER array, dimension(batchCount + 1). Each element N specifies the order of each matrix A. N >= 0.
[in,out]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array A 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 A 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 info_array = 0, each pointer is the factor U or L from the Cholesky factorization dA = U**H * U or dA = L * L**H.
[in]lddaINTEGER array, dimension(batchCount + 1). Each element LDDA specifies the leading dimension of each array A. LDDA >= max(1,N). To benefit from coalescent memory accesses LDDA must be divisible by 16.
[out]info_arrayArray of INTEGERs, dimension (batchCount), for corresponding matrices.
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.