Linear solvers

HPD solve

Solves either \(AX=B\) or \(A^T X=B\) for \(X\) given Hermitian positive-definite (HPD) \(A\) and right-hand side matrix \(B\). The solution is computed by first finding the Cholesky factorization of \(A\) and then performing two successive triangular solves against \(B\).

void HPDSolve(UpperOrLower uplo, Orientation orientation, Matrix<F> &A, Matrix<F> &B)
void HPDSolve(UpperOrLower uplo, Orientation orientation, DistMatrix<F> &A, DistMatrix<F> &B)

Overwrite B with the solution to \(AX=B\) or \(A^T X=B\), where A is Hermitian positive-definite and only the triangle of A specified by uplo is accessed.

Gaussian elimination

Solves \(AX=B\) for \(X\) given a general square nonsingular matrix \(A\) and right-hand side matrix \(B\). The solution is computed through (partially pivoted) Gaussian elimination.

void GaussianElimination(Matrix<F> &A, Matrix<F> &B)
void GaussianElimination(DistMatrix<F> &A, DistMatrix<F> &B)

Upon completion, \(A\) will have been overwritten with Gaussian elimination and \(B\) will be overwritten with \(X\).

Least-squares

Solves \(AX=B\) or \(A^H X = B\) for \(X\) in a least-squares sense given a general full-rank matrix \(A \in \mathbb{F}^{m \times n}\). If \(m \ge n\), then the first step is to form the QR factorization of \(A\), otherwise the LQ factorization is computed.

  • If solving \(AX=B\), then either \(X=R^{-1} Q^H B\) or \(X=Q^H L^{-1} B\).

  • If solving \(A^H X=B\), then either \(X=Q R^{-H} B\) or \(X=L^{-H} Q B\).

void LeastSquares(Orientation orientation, Matrix<F> &A, const Matrix<F> &B, Matrix<F> &X)
void LeastSquares(Orientation orientation, DistMatrix<F> &A, const DistMatrix<F> &B, DistMatrix<F> &X)

If orientation is set to NORMAL, then solve \(AX=B\), otherwise orientation must be equal to ADJOINT and \(A^H X=B\) will be solved. Upon completion, \(A\) is overwritten with its QR or LQ factorization, and \(X\) is overwritten with the solution.

Solve after Cholesky

Uses an existing in-place Cholesky factorization to solve against one or more right-hand sides.

void cholesky::SolveAfter(UpperOrLower uplo, Orientation orientation, const Matrix<F> &A, Matrix<F> &B)
void cholesky::SolveAfter(UpperOrLower uplo, Orientation orientation, const DistMatrix<F> &A, DistMatrix<F> &B)

Update \(B := A^{-1} B\), \(B := A^{-T} B\), or \(B := A^{-H} B\), where one triangle of \(A\) has been overwritten with its Cholesky factor.

Solve after LU

Uses an existing in-place LU factorization (with or without partial pivoting) to solve against one or more right-hand sides.

void lu::SolveAfter(Orientation orientation, const Matrix<F> &A, Matrix<F> &B)
void lu::SolveAfter(Orientation orientation, const DistMatrix<F> &A, DistMatrix<F> &B)

Update \(B := A^{-1} B\), \(B := A^{-T} B\), or \(B := A^{-H} B\), where \(A\) has been overwritten with its LU factors (without partial pivoting).

void lu::SolveAfter(Orientation orientation, const Matrix<F> &A, const Matrix<int> &p, Matrix<F> &B)
void lu::SolveAfter(Orientation orientation, const DistMatrix<F> &A, const DistMatrix<int, VC, STAR> &p, DistMatrix<F> &B)

Update \(B := A^{-1} B\), \(B := A^{-T} B\), or \(B := A^{-H} B\), where \(A\) has been overwritten with its LU factors with partial pivoting, which satisfy \(P A = L U\), where the permutation matrix \(P\) is represented by the pivot vector p.

void lu::SolveAfter(Orientation orientation, const Matrix<F> &A, const Matrix<int> &p, const Matrix<int> &q, Matrix<F> &B)
void lu::SolveAfter(Orientation orientation, const DistMatrix<F> &A, const DistMatrix<int, VC, STAR> &p, const DistMatrix<int, VC, STAR> &q, DistMatrix<F> &B)

Update \(B := A^{-1} B\), \(B := A^{-T} B\), or \(B := A^{-H} B\), where \(A\) has been overwritten with its LU factors with full pivoting, which satisfy \(P A Q = L U\), where the permutation matrices \(P\) and \(Q\) are represented by the pivot vector p and q, respectively.