LU Factorization and Solve#
GESV (GEneral linear SolVe) function solves the system
by factoring A
with LU factorization and overwriting B
with the solution X
. The function GESV is equivalent to executing a GETRS immediately followed by GETRF.
cuSolverDx provides separate function operators depending on whether pivoting is to be performed:
cusolverdx::function::gesv_no_pivot
: Linear system solve after LU factorization with no pivotingcusolverdx::function::gesv_partial_pivot
: Linear system solve after LU factorization with partial pivoting
Note
If a nonsingular matrix A
is diagonal dominant, then it is safe to compute \(A = L * U\) without pivoting. If a matrix is not diagonal dominant, then pivoting is usually required to ensure numerical stability.
cuSolverDx gesv_no_pivot
device functions are (see Execution Methods):
__device__ void execute(data_type* A, data_type* B, status_type* info);
// with runtime lda and ldb
__device__ void execute(data_type* A, const unsigned int lda,
data_type* B, status_type* info);
__device__ void execute(data_type* A,
data_type* B, const unsigned int ldb, status_type* info);
__device__ void execute(data_type* A, const unsigned int lda,
data_type* B, const unsigned int ldb, status_type* info);
cuSolverDx gesv_partial_pivot
device functions are:
__device__ void execute(data_type* A, int* ipiv, data_type* B, status_type* info);
// with runtime lda and ldb
__device__ void execute(data_type* A, const unsigned int lda,
int* ipiv,
data_type* B, status_type* info);
__device__ void execute(data_type* A,
int* ipiv,
data_type* B, const unsigned int ldb, status_type* info);
__device__ void execute(data_type* A, const unsigned int lda, int* ipiv,
data_type* B, const unsigned int ldb, status_type* info);
After the function, input A
is replaced by the lower triangular LU factor L
, and the upper triangular LU factor U
.
The output status parameter, info
, is an array of batch size. If LU factorization is successful, every element of info
is zero.
If LU factorization failed for any batches, i.e. A
(U
) is singular, info[batch_id] = i
indicates U(i,i) = 0
.
For getrf_partial_pivot
, the output ipiv
is an array of size min(M, N)
for each batch, and ipiv[batch_id, i]
indicates the row i
interchanges with row ipiv[batch_id, i]
on the batch_id
-th batch of A
.
B
is a batched N x K
righ-hand side matrix. After the function, matrix X
overwrites matrix B
with the same leading dimension ldb
. The leading dimension of B
is ldb >= K
if B
is column-major, or ldb >= N
if B
is row-major.
The functions support:
A
andB
either being the same or different column- or row-major layouts, see Arrangement operator, and\(op(A)\) either being
non_transposed
,transposed
for real data type, orconj_transposed
for complex data type, see TransposeMode operator.