Linear solve¶
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
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 |
|
|
|
|
|
|
|
|
C++ API¶
-
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
(AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B)¶