![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_cgetri_expert_gpu_work (magma_int_t n, magmaFloatComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info, magma_mode_t mode, void *host_work, magma_int_t *lwork_host, void *device_work, magma_int_t *lwork_device, magma_queue_t queues[2]) |
CGETRI computes the inverse of a matrix using the LU factorization computed by CGETRF. | |
magma_int_t | magma_cgetri_gpu (magma_int_t n, magmaFloatComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaFloatComplex_ptr dwork, magma_int_t lwork, magma_int_t *info) |
CGETRI computes the inverse of a matrix using the LU factorization computed by CGETRF. | |
magma_int_t | magma_dgetri_expert_gpu_work (magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info, magma_mode_t mode, void *host_work, magma_int_t *lwork_host, void *device_work, magma_int_t *lwork_device, magma_queue_t queues[2]) |
DGETRI computes the inverse of a matrix using the LU factorization computed by DGETRF. | |
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. | |
magma_int_t | magma_sgetri_expert_gpu_work (magma_int_t n, magmaFloat_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info, magma_mode_t mode, void *host_work, magma_int_t *lwork_host, void *device_work, magma_int_t *lwork_device, magma_queue_t queues[2]) |
SGETRI computes the inverse of a matrix using the LU factorization computed by SGETRF. | |
magma_int_t | magma_sgetri_gpu (magma_int_t n, magmaFloat_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaFloat_ptr dwork, magma_int_t lwork, magma_int_t *info) |
SGETRI computes the inverse of a matrix using the LU factorization computed by SGETRF. | |
magma_int_t | magma_zgetri_expert_gpu_work (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info, magma_mode_t mode, void *host_work, magma_int_t *lwork_host, void *device_work, magma_int_t *lwork_device, magma_queue_t queues[2]) |
ZGETRI computes the inverse of a matrix using the LU factorization computed by ZGETRF. | |
magma_int_t | magma_zgetri_gpu (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *ipiv, magmaDoubleComplex_ptr dwork, magma_int_t lwork, magma_int_t *info) |
ZGETRI computes the inverse of a matrix using the LU factorization computed by ZGETRF. | |
magma_int_t magma_cgetri_expert_gpu_work | ( | magma_int_t | n, |
magmaFloatComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info, | ||
magma_mode_t | mode, | ||
void * | host_work, | ||
magma_int_t * | lwork_host, | ||
void * | device_work, | ||
magma_int_t * | lwork_device, | ||
magma_queue_t | queues[2] ) |
CGETRI computes the inverse of a matrix using the LU factorization computed by CGETRF.
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 CGESV, or CGETRF and CGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by CGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from CGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | info | INTEGER
|
[in] | mode | magma_mode_t
|
[in,out] | host_work | Workspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned() |
[in,out] | lwork_host | INTEGER pointer The size of the workspace (host_work) in bytes
|
[in,out] | device_work | Workspace, allocated on device (GPU) memory. |
[in,out] | lwork_device | INTEGER pointer The size of the workspace (device_work) in bytes
|
[in] | queue | magma_queue_t
|
magma_int_t magma_cgetri_gpu | ( | magma_int_t | n, |
magmaFloatComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magmaFloatComplex_ptr | dwork, | ||
magma_int_t | lwork, | ||
magma_int_t * | info ) |
CGETRI computes the inverse of a matrix using the LU factorization computed by CGETRF.
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 CGESV, or CGETRF and CGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by CGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from CGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | dwork | (workspace) COMPLEX array on the GPU, dimension (MAX(1,LWORK)) |
[in] | lwork | INTEGER The dimension of the array DWORK. LWORK >= N*NB, where NB is the optimal blocksize returned by magma_get_cgetri_nb(n). Unlike LAPACK, this version does not currently support a workspace query, because the workspace is on the GPU. |
[out] | info | INTEGER
|
magma_int_t magma_dgetri_expert_gpu_work | ( | magma_int_t | n, |
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info, | ||
magma_mode_t | mode, | ||
void * | host_work, | ||
magma_int_t * | lwork_host, | ||
void * | device_work, | ||
magma_int_t * | lwork_device, | ||
magma_queue_t | queues[2] ) |
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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE 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] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from DGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | info | INTEGER
|
[in] | mode | magma_mode_t
|
[in,out] | host_work | Workspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned() |
[in,out] | lwork_host | INTEGER pointer The size of the workspace (host_work) in bytes
|
[in,out] | device_work | Workspace, allocated on device (GPU) memory. |
[in,out] | lwork_device | INTEGER pointer The size of the workspace (device_work) in bytes
|
[in] | queue | magma_queue_t
|
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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE 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] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER 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] | lwork | INTEGER 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] | info | INTEGER
|
magma_int_t magma_sgetri_expert_gpu_work | ( | magma_int_t | n, |
magmaFloat_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info, | ||
magma_mode_t | mode, | ||
void * | host_work, | ||
magma_int_t * | lwork_host, | ||
void * | device_work, | ||
magma_int_t * | lwork_device, | ||
magma_queue_t | queues[2] ) |
SGETRI computes the inverse of a matrix using the LU factorization computed by SGETRF.
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 SGESV, or SGETRF and SGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | REAL array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by SGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from SGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | info | INTEGER
|
[in] | mode | magma_mode_t
|
[in,out] | host_work | Workspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned() |
[in,out] | lwork_host | INTEGER pointer The size of the workspace (host_work) in bytes
|
[in,out] | device_work | Workspace, allocated on device (GPU) memory. |
[in,out] | lwork_device | INTEGER pointer The size of the workspace (device_work) in bytes
|
[in] | queue | magma_queue_t
|
magma_int_t magma_sgetri_gpu | ( | magma_int_t | n, |
magmaFloat_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magmaFloat_ptr | dwork, | ||
magma_int_t | lwork, | ||
magma_int_t * | info ) |
SGETRI computes the inverse of a matrix using the LU factorization computed by SGETRF.
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 SGESV, or SGETRF and SGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | REAL array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by SGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from SGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | dwork | (workspace) REAL array on the GPU, dimension (MAX(1,LWORK)) |
[in] | lwork | INTEGER The dimension of the array DWORK. LWORK >= N*NB, where NB is the optimal blocksize returned by magma_get_sgetri_nb(n). Unlike LAPACK, this version does not currently support a workspace query, because the workspace is on the GPU. |
[out] | info | INTEGER
|
magma_int_t magma_zgetri_expert_gpu_work | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info, | ||
magma_mode_t | mode, | ||
void * | host_work, | ||
magma_int_t * | lwork_host, | ||
void * | device_work, | ||
magma_int_t * | lwork_device, | ||
magma_queue_t | queues[2] ) |
ZGETRI computes the inverse of a matrix using the LU factorization computed by ZGETRF.
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 ZGESV, or ZGETRF and ZGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX_16 array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by ZGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | info | INTEGER
|
[in] | mode | magma_mode_t
|
[in,out] | host_work | Workspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned() |
[in,out] | lwork_host | INTEGER pointer The size of the workspace (host_work) in bytes
|
[in,out] | device_work | Workspace, allocated on device (GPU) memory. |
[in,out] | lwork_device | INTEGER pointer The size of the workspace (device_work) in bytes
|
[in] | queue | magma_queue_t
|
magma_int_t magma_zgetri_gpu | ( | magma_int_t | n, |
magmaDoubleComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magmaDoubleComplex_ptr | dwork, | ||
magma_int_t | lwork, | ||
magma_int_t * | info ) |
ZGETRI computes the inverse of a matrix using the LU factorization computed by ZGETRF.
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 ZGESV, or ZGETRF and ZGETRS, 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.
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX_16 array on the GPU, dimension (LDDA,N) On entry, the factors L and U from the factorization A = P*L*U as computed by ZGETRF_GPU. On exit, if INFO = 0, the inverse of the original matrix A. |
[in] | ldda | INTEGER The leading dimension of the array A. LDDA >= max(1,N). |
[in] | ipiv | INTEGER array, dimension (N) The pivot indices from ZGETRF; for 1 <= i <= N, row i of the matrix was interchanged with row IPIV(i). |
[out] | dwork | (workspace) COMPLEX_16 array on the GPU, dimension (MAX(1,LWORK)) |
[in] | lwork | INTEGER The dimension of the array DWORK. LWORK >= N*NB, where NB is the optimal blocksize returned by magma_get_zgetri_nb(n). Unlike LAPACK, this version does not currently support a workspace query, because the workspace is on the GPU. |
[out] | info | INTEGER
|