SLATE 2024.05.31
Software for Linear Algebra Targeting Exascale
|
Solve \(AX = B\). More...
Functions | |
template<typename scalar_t > | |
int64_t | slate::pbsv (HermitianBandMatrix< scalar_t > &A, Matrix< scalar_t > &B, Options const &opts) |
Distributed parallel Cholesky factorization and solve. | |
Solve \(AX = B\).
int64_t slate::pbsv | ( | HermitianBandMatrix< scalar_t > & | A, |
Matrix< scalar_t > & | B, | ||
Options const & | opts | ||
) |
Distributed parallel Cholesky factorization and solve.
Computes the solution to a system of linear equations
\[ A X = B, \]
where \(A\) is an n-by-n Hermitian positive definite band matrix and \(X\) and \(B\) are n-by-nrhs matrices. The Cholesky decomposition is used to factor \(A\) as
\[ A = L L^H, \]
if \(A\) is stored lower, where \(L\) is a lower triangular band matrix, or
\[ A = U^H U, \]
if \(A\) is stored upper, where \(U\) is an upper triangular band matrix. The factored form of \(A\) is then used to solve the system of equations \(A X = B\).
scalar_t | One of float, double, std::complex<float>, std::complex<double>. |
[in,out] | A | On entry, the n-by-n Hermitian positive definite band matrix \(A\). On exit, if return value = 0, overwritten by the factor \(U\) or \(L\) from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\). If scalar_t is real, \(A\) can be a SymmetricMatrix object. |
[in,out] | B | On entry, the n-by-nrhs right hand side matrix \(B\). On exit, if return value = 0, the n-by-nrhs solution matrix \(X\). |
[in] | opts | Additional options, as map of name = value pairs. Possible options:
|