MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
Loading...
Searching...
No Matches

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.
 

Detailed Description

Function Documentation

◆ magma_cgetri_expert_gpu_work()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX 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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from CGETRF; for 1 <= i <= 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
  • > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its cannot be computed.
[in]modemagma_mode_t
  • specifies execution mode (hybrid vs. native)
  • currently ignored, reserved for future use
[in,out]host_workWorkspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned()
[in,out]lwork_hostINTEGER pointer The size of the workspace (host_work) in bytes
  • lwork_host[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_host. The workspace itself is not referenced, and no factorization is performed.
  • lwork_host[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_host.
Parameters
[in,out]device_workWorkspace, allocated on device (GPU) memory.
[in,out]lwork_deviceINTEGER pointer The size of the workspace (device_work) in bytes
  • lwork_device[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_device. The workspace itself is not referenced, and no factorization is performed.
  • lwork_device[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_device.
[in]queuemagma_queue_t
  • created/destroyed by the user outside the routine
  • Used for kernel execution on the GPU

◆ magma_cgetri_gpu()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX 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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER 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]lworkINTEGER 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]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_dgetri_expert_gpu_work()

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.

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]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.
[in]modemagma_mode_t
  • specifies execution mode (hybrid vs. native)
  • currently ignored, reserved for future use
[in,out]host_workWorkspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned()
[in,out]lwork_hostINTEGER pointer The size of the workspace (host_work) in bytes
  • lwork_host[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_host. The workspace itself is not referenced, and no factorization is performed.
  • lwork_host[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_host.
Parameters
[in,out]device_workWorkspace, allocated on device (GPU) memory.
[in,out]lwork_deviceINTEGER pointer The size of the workspace (device_work) in bytes
  • lwork_device[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_device. The workspace itself is not referenced, and no factorization is performed.
  • lwork_device[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_device.
[in]queuemagma_queue_t
  • created/destroyed by the user outside the routine
  • Used for kernel execution on the GPU

◆ magma_dgetri_gpu()

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_sgetri_expert_gpu_work()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dAREAL 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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from SGETRF; for 1 <= i <= 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
  • > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its cannot be computed.
[in]modemagma_mode_t
  • specifies execution mode (hybrid vs. native)
  • currently ignored, reserved for future use
[in,out]host_workWorkspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned()
[in,out]lwork_hostINTEGER pointer The size of the workspace (host_work) in bytes
  • lwork_host[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_host. The workspace itself is not referenced, and no factorization is performed.
  • lwork_host[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_host.
Parameters
[in,out]device_workWorkspace, allocated on device (GPU) memory.
[in,out]lwork_deviceINTEGER pointer The size of the workspace (device_work) in bytes
  • lwork_device[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_device. The workspace itself is not referenced, and no factorization is performed.
  • lwork_device[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_device.
[in]queuemagma_queue_t
  • created/destroyed by the user outside the routine
  • Used for kernel execution on the GPU

◆ magma_sgetri_gpu()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dAREAL 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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER 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]lworkINTEGER 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]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_zgetri_expert_gpu_work()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX_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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER array, dimension (N) The pivot indices from ZGETRF; for 1 <= i <= 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
  • > 0: if INFO = i, U(i,i) is exactly zero; the matrix is singular and its cannot be computed.
[in]modemagma_mode_t
  • specifies execution mode (hybrid vs. native)
  • currently ignored, reserved for future use
[in,out]host_workWorkspace, allocated on host (CPU) memory. For faster CPU-GPU communication, user can allocate it as pinned memory using magma_malloc_pinned()
[in,out]lwork_hostINTEGER pointer The size of the workspace (host_work) in bytes
  • lwork_host[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_host. The workspace itself is not referenced, and no factorization is performed.
  • lwork_host[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_host.
Parameters
[in,out]device_workWorkspace, allocated on device (GPU) memory.
[in,out]lwork_deviceINTEGER pointer The size of the workspace (device_work) in bytes
  • lwork_device[0] < 0: a workspace query is assumed, the routine calculates the required amount of workspace and returns it in lwork_device. The workspace itself is not referenced, and no factorization is performed.
  • lwork_device[0] >= 0: the routine assumes that the user has provided a workspace with the size in lwork_device.
[in]queuemagma_queue_t
  • created/destroyed by the user outside the routine
  • Used for kernel execution on the GPU

◆ magma_zgetri_gpu()

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.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX_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]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,N).
[in]ipivINTEGER 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]lworkINTEGER 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]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.