17#ifdef MFEM_USE_STRUMPACK
25#include "StrumpackSparseSolverMPIDist.hpp"
26#include "StrumpackSparseSolverMixedPrecisionMPIDist.hpp"
44 bool sym_sparse =
false);
57 MFEM_ABORT(
"STRUMPACKRowLocMatrix::Mult: Matrix vector products are not "
62 MPI_Comm
GetComm()
const {
return A_->comm(); }
65 strumpack::CSRMatrixMPI<double, HYPRE_BigInt> *
GetA()
const {
return A_; }
68 strumpack::CSRMatrixMPI<double, HYPRE_BigInt> *A_;
77template <
typename STRUMPACKSolverType>
251#if STRUMPACK_VERSION_MAJOR >= 5
268 void FactorInternal()
const;
284 SparseSolverMPIDist<double, HYPRE_BigInt>>
314#if STRUMPACK_VERSION_MAJOR >= 7
317 SparseSolverMixedPrecisionMPIDist<float, double, HYPRE_BigInt>>
340 int argc,
char *argv[]);
STRUMPACKMixedPrecisionSolver(MPI_Comm comm)
Constructor with MPI_Comm parameter.
~STRUMPACKMixedPrecisionSolver()
Destructor.
STRUMPACKRowLocMatrix(MPI_Comm comm, int num_loc_rows, HYPRE_BigInt first_loc_row, HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols, int *I, HYPRE_BigInt *J, double *data, bool sym_sparse=false)
Creates a general parallel matrix from a local CSR matrix on each processor.
void Mult(const Vector &x, Vector &y) const
Matrix vector products are not supported for this type of matrix.
strumpack::CSRMatrixMPI< double, HYPRE_BigInt > * GetA() const
Get access to the internal CSR matrix.
MPI_Comm GetComm() const
Get the MPI Comm being used by the parallel matrix.
void SetOperator(const Operator &op)
Set the operator/matrix.
void SetAbsTol(double atol)
Set the absolute tolerance for iterative solvers.
void ArrayMult(const Array< const Vector * > &X, Array< Vector * > &Y) const
Factor and solve the linear systems across the array of vectors.
void SetMaxIter(int max_it)
Set the maximum number of iterations for iterative solvers.
STRUMPACKSolverType * solver_
void SetMatching(strumpack::MatchingJob job)
Configure static pivoting for stability.
void SetCompression(strumpack::CompressionType type)
Select compression for sparse data types.
STRUMPACKSolverBase(MPI_Comm comm, int argc, char *argv[])
Constructor with MPI_Comm parameter and command line arguments.
void SetPrintFactorStatistics(bool print_stat)
Set up verbose printing during the factor step.
void SetKrylovSolver(strumpack::KrylovSolver method)
Set the Krylov solver method to use.
void Mult(const Vector &x, Vector &y) const
Factor and solve the linear system .
void SetCompressionButterflyLevels(int levels)
Set the number of butterfly levels for the HODLR compression option.
void SetCompressionAbsTol(double atol)
Set the absolute tolerance for low rank compression methods.
void EnableGPU()
Enable GPU off-loading available if STRUMPACK was compiled with CUDA.
virtual ~STRUMPACKSolverBase()
Default destructor.
void DisableGPU()
Disable GPU off-loading available if STRUMPACK was compiled with CUDA.
const STRUMPACKRowLocMatrix * APtr_
void SetFromCommandLine()
Set options that were captured from the command line.
void SetRelTol(double rtol)
Set the relative tolerance for interative solvers.
void SetCompressionLossyPrecision(int precision)
Set the precision for the lossy compression option.
void SetReorderingReuse(bool reuse)
Set the flag controlling reuse of the symbolic factorization for multiple operators.
void SetPrintSolveStatistics(bool print_stat)
Set up verbose printing during the solve step.
void SetCompressionRelTol(double rtol)
Set the relative tolerance for low rank compression methods.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Set matrix reordering strategy.
STRUMPACKSolver(MPI_Comm comm)
Constructor with the MPI Comm parameter.
~STRUMPACKSolver()
Destructor.
MFEM_DEPRECATED STRUMPACKSolver(int argc, char *argv[], MPI_Comm comm)