Symmetric solve

Implementation

Solve \(AX=B\), \(A^T X = B\), or \(A^H X = B\) for \(X\) given a Hermitian matrix \(A\) and a right-hand side matrix \(B\). When \(A\) is dense, the default algorithm is Bunch-Kaufman, whereas, when \(A\) is sparse, the default approach embeds the problem in the same manner as LinearSolve(), though it is possible to override this behaviour and attempt a sparse \(LDL\) factorization which only dynamically pivots within supernodes.

Note

Only the lower-triangular storage case (uplo=LOWER) is supported by the following routines.

Python API

SymmetricSolve(A, B[, tryLDL=False, conjugate=False, uplo=LOWER, orient=NORMAL])
Parameters
  • A – Dense or sparse symmetric matrix to solve against

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

  • tryLDL – (optional) if a sparse \(LDL^T\) or \(LDL^H\) factorization without pivoting should be attempted instead of embedding in a QSD system

  • conjugate – (optional) whether the matrix is equal to its transpose or conjugate-transpose

  • uplo – (optional) which triangle the data is explicitly stored in

  • orient – (optional) whether to use \(A\), \(A^T\), or \(A^H\)

Type of \(A\)

Type of B

Matrix

Matrix

AbstractDistMatrix

AbstractDistMatrix

SparseMatrix

Matrix

DistSparseMatrix

DistMultiVec

C++ API

void SymmetricSolve(UpperOrLower uplo, Orientation orientation, const Matrix<F> &A, Matrix<F> &B, bool conjugate = false, const LDLPivotCtrl<Base<F>> &ctrl = LDLPivotCtrl<Base<F>>())
void SymmetricSolve(UpperOrLower uplo, Orientation orientation, const AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B, bool conjugate = false, const LDLPivotCtrl<Base<F>> &ctrl = LDLPivotCtrl<Base<F>>())
void SymmetricSolve(const SparseMatrix<F> &A, Matrix<F> &B, bool conjugate = false, bool tryLDL = false, const BisectCtrl &ctrl = BisectCtrl())
void SymmetricSolve(const DistSparseMatrix<F> &A, DistMultiVec<F> &B, bool conjugate = false, bool tryLDL = false, const BisectCtrl &ctrl = BisectCtrl())

Dense versions which factor in-place

void symm_solve::Overwrite(UpperOrLower uplo, Orientation orientation, Matrix<F> &A, Matrix<F> &B, bool conjugate = false, const LDLPivotCtrl<Base<F>> &ctrl = LDLPivotCtrl<Base<F>>())
void symm_solve::Overwrite(UpperOrLower uplo, Orientation orientation, AbstractDistMatrix<F> &A, AbstractDistMatrix<F> &B, bool conjugate = false, const LDLPivotCtrl<Base<F>> &ctrl = LDLPivotCtrl<Base<F>>())

C API

Single-precision

ElError ElSymmetricSolve_s(ElUpperOrLower uplo, ElOrientation orientation, ElConstMatrix_s A, ElMatrix_s B)
ElError ElSymmetricSolveDist_s(ElUpperOrLower uplo, ElOrientation orientation, ElConstDistMatrix_s A, ElDistMatrix_s B)
ElError ElSymmetricSolveSparse_s(ElConstSparseMatrix_s A, ElMatrix_s B, bool tryLDL)
ElError ElSymmetricSolveDistSparse_s(ElConstDistSparseMatrix_s A, ElDistMultiVec_s B, bool tryLDL)

Double-precision

ElError ElSymmetricSolve_d(ElUpperOrLower uplo, ElOrientation orientation, ElConstMatrix_d A, ElMatrix_d B)
ElError ElSymmetricSolveDist_d(ElUpperOrLower uplo, ElOrientation orientation, ElConstDistMatrix_d A, ElDistMatrix_d B)
ElError ElSymmetricSolveSparse_d(ElConstSparseMatrix_d A, ElMatrix_d B, bool tryLDL)
ElError ElSymmetricSolveDistSparse_d(ElConstDistSparseMatrix_d A, ElDistMultiVec_d B, bool tryLDL)

Single-precision complex

ElError ElSymmetricSolve_c(ElUpperOrLower uplo, ElOrientation orientation, ElConstMatrix_c A, ElMatrix_c B)
ElError ElSymmetricSolveDist_c(ElUpperOrLower uplo, ElOrientation orientation, ElConstDistMatrix_c A, ElDistMatrix_c B)
ElError ElSymmetricSolveSparse_c(ElConstSparseMatrix_c A, ElMatrix_c B, bool tryLDL)
ElError ElSymmetricSolveDistSparse_c(ElConstDistSparseMatrix_c A, ElDistMultiVec_c B, bool tryLDL)

Double-precision complex

ElError ElSymmetricSolve_z(ElUpperOrLower uplo, ElOrientation orientation, ElConstMatrix_z A, ElMatrix_z B)
ElError ElSymmetricSolveDist_z(ElUpperOrLower uplo, ElOrientation orientation, ElConstDistMatrix_z A, ElDistMatrix_z B)
ElError ElSymmetricSolveSparse_z(ElConstSparseMatrix_z A, ElMatrix_z B, bool tryLDL)
ElError ElSymmetricSolveDistSparse_z(ElConstDistSparseMatrix_z A, ElDistMultiVec_z B, bool tryLDL)