Linear solve

Implementation

Solves \(AX=B\) for \(X\) given a general square nonsingular matrix \(A\) and right-hand side matrix \(B\). In all cases, on exit, the following routines overwrite \(B\) with \(\text{inv}(A) B\).

Dense algorithm

For dense matrices, the solution is computed through Gaussian elimination with partial pivoting.

Sparse-direct algorithm

For sparse matrices, the original problem is embedded into the augmented system

\[\begin{split}\begin{pmatrix} \alpha I & A \\ A^H & 0 \end{pmatrix} \begin{pmatrix} R/\alpha \\ X \end{pmatrix} = \begin{pmatrix} B \\ 0 \end{pmatrix},\end{split}\]

where \(\alpha \approx \sigma_{\text{min}}(A)\) is known to be near-optimal with respect to minimizing the condition number of the augmented system. A priori regularization is added to the linear system, a sparse-direct Cholesky-like factorization is performed, and the factorization is used as a preconditioner for Flexible GMRES.

It is worth noting that this is a special case of LeastSquares(), and that, unlike cases where \(A\) is not invertible, it is possible for LinearSolve() to choose \(\alpha=0\).

Python API

LinearSolve(A, B[, ctrl=None])
Parameters
  • A – Dense or sparse matrix to solve against

  • B – Right-hand sides (which will be overwritten with the solution).

  • ctrl – (optional) sparse-direct Least Squares control structure

Type of \(A\)

Type of B

Matrix

Matrix

AbstractDistMatrix

AbstractDistMatrix

SparseMatrix

Matrix

DistSparseMatrix

DistMultiVec

C++ API

void LinearSolve(const Matrix<F> &A, Matrix<F> &B)
void LinearSolve(const AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B)
void LinearSolve(const SparseMatrix<F> &A, Matrix<F> &B, const LeastSquaresCtrl<Base<F>> &ctrl = LeastSquaresCtrl<Base<F>>())
void LinearSolve(const DistSparseMatrix<F> &A, DistMultiVec<F> &B, const LeastSquaresCtrl<Base<F>> &ctrl = LeastSquaresCtrl<Base<F>>())

Dense in-place alternatives

The following routines factor \(A\) in place.

lin_solve::Overwrite(Matrix<F> &A, Matrix<F> &B)
lin_solve::Overwrite(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B)

C API

Standard interface

Single-precision

ElError ElLinearSolve_s(ElConstMatrix_s A, ElMatrix_s B)
ElError ElLinearSolveDist_s(ElConstDistMatrix_s A, ElDistMatrix_s B)
ElError ElLinearSolveSparse_s(ElConstSparseMatrix_s A, ElMatrix_s B)
ElError ElLinearSolveDistSparse_s(ElConstDistSparseMatrix_s A, ElDistMultiVec_s B)

Double-precision

ElError ElLinearSolve_d(ElConstMatrix_d A, ElMatrix_d B)
ElError ElLinearSolveDist_d(ElConstDistMatrix_d A, ElDistMatrix_d B)
ElError ElLinearSolveSparse_d(ElConstSparseMatrix_d A, ElMatrix_d B)
ElError ElLinearSolveDistSparse_d(ElConstDistSparseMatrix_d A, ElDistMultiVec_d B)

Single-precision complex

ElError ElLinearSolve_c(ElConstMatrix_c A, ElMatrix_c B)
ElError ElLinearSolveDist_c(ElConstDistMatrix_c A, ElDistMatrix_c B)
ElError ElLinearSolveSparse_c(ElConstSparseMatrix_c A, ElMatrix_c B)
ElError ElLinearSolveDistSparse_c(ElConstDistSparseMatrix_c A, ElDistMultiVec_c B)

Double-precision complex

ElError ElLinearSolve_z(ElConstMatrix_z A, ElMatrix_z B)
ElError ElLinearSolveDist_z(ElConstDistMatrix_z A, ElDistMatrix_z B)
ElError ElLinearSolveSparse_z(ElConstSparseMatrix_z A, ElMatrix_z B)
ElError ElLinearSolveDistSparse_z(ElConstDistSparseMatrix_z A, ElDistMultiVec_z B)

Expert versions

Single-precision

ElError ElLinearSolveXSparse_s(ElConstSparseMatrix_s A, ElMatrix_s B, ElLeastSquaresCtrl_s ctrl)
ElError ElLinearSolveXDistSparse_s(ElConstDistSparseMatrix_s A, ElDistMultiVec_s B, ElLeastSquaresCtrl_s ctrl)

Double-precision

ElError ElLinearSolveXSparse_d(ElConstSparseMatrix_d A, ElMatrix_d B, ElLeastSquaresCtrl_d ctrl)
ElError ElLinearSolveXDistSparse_d(ElConstDistSparseMatrix_d A, ElDistMultiVec_d B, ElLeastSquaresCtrl_d ctrl)

Single-precision complex

ElError ElLinearSolveXSparse_c(ElConstSparseMatrix_c A, ElMatrix_c B, ElLeastSquaresCtrl_s ctrl)
ElError ElLinearSolveXDistSparse_c(ElConstDistSparseMatrix_c A, ElDistMultiVec_c B, ElLeastSquaresCtrl_s ctrl)

Double-precision complex

ElError ElLinearSolveXSparse_z(ElConstSparseMatrix_z A, ElMatrix_z B, ElLeastSquaresCtrl_d ctrl)
ElError ElLinearSolveXDistSparse_z(ElConstDistSparseMatrix_z A, ElDistMultiVec_z B, ElLeastSquaresCtrl_d ctrl)