MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::STRUMPACKSolverBase< STRUMPACKSolverType > Class Template Reference

#include <strumpack.hpp>

Inheritance diagram for mfem::STRUMPACKSolverBase< STRUMPACKSolverType >:
[legend]
Collaboration diagram for mfem::STRUMPACKSolverBase< STRUMPACKSolverType >:
[legend]

Public Member Functions

virtual ~STRUMPACKSolverBase ()
 Default destructor.
 
void Mult (const Vector &x, Vector &y) const
 Factor and solve the linear system \(y = Op^{-1} x \).
 
void ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Factor and solve the linear systems \( Y_i = Op^{-1} X_i \) across the array of vectors.
 
void SetOperator (const Operator &op)
 Set the operator/matrix.
 
void SetFromCommandLine ()
 Set options that were captured from the command line.
 
void SetPrintFactorStatistics (bool print_stat)
 Set up verbose printing during the factor step.
 
void SetPrintSolveStatistics (bool print_stat)
 Set up verbose printing during the solve step.
 
void SetRelTol (double rtol)
 Set the relative tolerance for interative solvers.
 
void SetAbsTol (double atol)
 Set the absolute tolerance for iterative solvers.
 
void SetMaxIter (int max_it)
 Set the maximum number of iterations for iterative solvers.
 
void SetReorderingReuse (bool reuse)
 Set the flag controlling reuse of the symbolic factorization for multiple operators.
 
void EnableGPU ()
 Enable GPU off-loading available if STRUMPACK was compiled with CUDA.
 
void DisableGPU ()
 Disable GPU off-loading available if STRUMPACK was compiled with CUDA.
 
void SetKrylovSolver (strumpack::KrylovSolver method)
 Set the Krylov solver method to use.
 
void SetReorderingStrategy (strumpack::ReorderingStrategy method)
 Set matrix reordering strategy.
 
void SetMatching (strumpack::MatchingJob job)
 Configure static pivoting for stability.
 
void SetCompression (strumpack::CompressionType type)
 Select compression for sparse data types.
 
void SetCompressionRelTol (double rtol)
 Set the relative tolerance for low rank compression methods.
 
void SetCompressionAbsTol (double atol)
 Set the absolute tolerance for low rank compression methods.
 
void SetCompressionLossyPrecision (int precision)
 Set the precision for the lossy compression option.
 
void SetCompressionButterflyLevels (int levels)
 Set the number of butterfly levels for the HODLR compression option.
 
- Public Member Functions inherited from mfem::Solver
 Solver (int s=0, bool iter_mode=false)
 Initialize a square Solver with size s.
 
 Solver (int h, int w, bool iter_mode=false)
 Initialize a Solver with height h and width w.
 
virtual void SetOperator (const Operator &op)=0
 Set/update the solver for the given operator.
 
- Public Member Functions inherited from mfem::Operator
void InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
 Initializes memory for true vectors of linear system.
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0).
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size).
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows().
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height().
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols().
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width().
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator.
 
virtual void Mult (const Vector &x, Vector &y) const =0
 Operator application: y=A(x).
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error.
 
virtual void AddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator application: y+=A(x) (default) or y+=a*A(x).
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
 
virtual void ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Operator application on a matrix: Y=A(X).
 
virtual void ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X).
 
virtual void ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
 
virtual void ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X).
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error.
 
virtual void AssembleDiagonal (Vector &diag) const
 Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed.
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity.
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach.
 
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
 Form a column-constrained linear system using a matrix-free approach.
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem().
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator.
 
void FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator (including constraints).
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator.
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format.
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format.
 
virtual ~Operator ()
 Virtual destructor.
 
Type GetType () const
 Return the type ID of the Operator class.
 

Protected Member Functions

 STRUMPACKSolverBase (MPI_Comm comm, int argc, char *argv[])
 Constructor with MPI_Comm parameter and command line arguments.
 
 STRUMPACKSolverBase (STRUMPACKRowLocMatrix &A, int argc, char *argv[])
 Constructor with STRUMPACK matrix object and command line arguments.
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator()
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator()
 
OperatorSetupRAP (const Operator *Pi, const Operator *Po)
 Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt".
 

Protected Attributes

const STRUMPACKRowLocMatrixAPtr_
 
STRUMPACKSolverType * solver_
 
bool factor_verbose_
 
bool solve_verbose_
 
bool reorder_reuse_
 
Vector rhs_
 
Vector sol_
 
int nrhs_
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix.
 
int width
 Dimension of the input / number of columns in the matrix.
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  DiagonalPolicy { DIAG_ZERO , DIAG_ONE , DIAG_KEEP }
 Defines operator diagonal policy upon elimination of rows and/or columns. More...
 
enum  Type {
  ANY_TYPE , MFEM_SPARSEMAT , Hypre_ParCSR , PETSC_MATAIJ ,
  PETSC_MATIS , PETSC_MATSHELL , PETSC_MATNEST , PETSC_MATHYPRE ,
  PETSC_MATGENERIC , Complex_Operator , MFEM_ComplexSparseMat , Complex_Hypre_ParCSR ,
  Complex_DenseMat , MFEM_Block_Matrix , MFEM_Block_Operator
}
 Enumeration defining IDs for some classes derived from Operator. More...
 
- Public Attributes inherited from mfem::Solver
bool iterative_mode
 If true, use the second argument of Mult() as an initial guess.
 

Detailed Description

template<typename STRUMPACKSolverType>
class mfem::STRUMPACKSolverBase< STRUMPACKSolverType >

The MFEM STRUMPACK Direct Solver class.

The mfem::STRUMPACKSolver class uses the STRUMPACK library to perform LU factorization of a parallel sparse matrix. The solver is capable of handling double precision types. See http://portal.nersc.gov/project/sparse/strumpack/.

Definition at line 78 of file strumpack.hpp.

Constructor & Destructor Documentation

◆ STRUMPACKSolverBase() [1/2]

template<typename STRUMPACKSolverType >
mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::STRUMPACKSolverBase ( MPI_Comm  comm,
int  argc,
char *  argv[] 
)
protected

Constructor with MPI_Comm parameter and command line arguments.

STRUMPACKSolverBase::SetFromCommandLine must be called for the command line arguments to be used.

Definition at line 120 of file strumpack.cpp.

◆ STRUMPACKSolverBase() [2/2]

template<typename STRUMPACKSolverType >
mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::STRUMPACKSolverBase ( STRUMPACKRowLocMatrix A,
int  argc,
char *  argv[] 
)
protected

Constructor with STRUMPACK matrix object and command line arguments.

STRUMPACKSolverBase::SetFromCommandLine must be called for the command line arguments to be used.

Definition at line 132 of file strumpack.cpp.

◆ ~STRUMPACKSolverBase()

template<typename STRUMPACKSolverType >
mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::~STRUMPACKSolverBase
virtual

Default destructor.

Definition at line 145 of file strumpack.cpp.

Member Function Documentation

◆ ArrayMult()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::ArrayMult ( const Array< const Vector * > &  X,
Array< Vector * > &  Y 
) const
virtual

Factor and solve the linear systems \( Y_i = Op^{-1} X_i \) across the array of vectors.

Reimplemented from mfem::Operator.

Definition at line 373 of file strumpack.cpp.

◆ DisableGPU()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::DisableGPU

Disable GPU off-loading available if STRUMPACK was compiled with CUDA.

Note
Input/Output from MFEM to STRUMPACK is all still through host memory.

Definition at line 209 of file strumpack.cpp.

◆ EnableGPU()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::EnableGPU

Enable GPU off-loading available if STRUMPACK was compiled with CUDA.

Note
Input/Output from MFEM to STRUMPACK is all still through host memory.

Definition at line 202 of file strumpack.cpp.

◆ Mult()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::Mult ( const Vector x,
Vector y 
) const
virtual

Factor and solve the linear system \(y = Op^{-1} x \).

Implements mfem::Operator.

Definition at line 346 of file strumpack.cpp.

◆ SetAbsTol()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetAbsTol ( double  atol)

Set the absolute tolerance for iterative solvers.

Definition at line 181 of file strumpack.cpp.

◆ SetCompression()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetCompression ( strumpack::CompressionType  type)

Select compression for sparse data types.

Enable support for rank-structured data formats, which can be used for compression within the sparse solver.

Supported compression types are:

  • NONE: No compression, purely direct solver (default)
  • HSS: HSS compression of frontal matrices
  • BLR: Block low-rank compression of fronts
  • HODLR: Hierarchically Off-diagonal Low-Rank compression of frontal matrices
  • BLR_HODLR: Block low-rank compression of medium fronts and Hierarchically Off-diagonal Low-Rank compression of large fronts
  • ZFP_BLR_HODLR: ZFP compression for small fronts, Block low-rank compression of medium fronts and Hierarchically Off-diagonal Low-Rank compression of large fronts
  • LOSSLESS: Lossless compression
  • LOSSY: Lossy compression

For versions of STRUMPACK < 5, we support only NONE, HSS, and BLR. BLR_HODLR and ZPR_BLR_HODLR are supported in STRUMPACK >= 6.

Definition at line 236 of file strumpack.cpp.

◆ SetCompressionAbsTol()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetCompressionAbsTol ( double  atol)

Set the absolute tolerance for low rank compression methods.

This currently affects BLR, HSS, and HODLR. Use STRUMPACKSolverBase::SetCompression to set the proper compression type.

Definition at line 275 of file strumpack.cpp.

◆ SetCompressionButterflyLevels()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetCompressionButterflyLevels ( int  levels)

Set the number of butterfly levels for the HODLR compression option.

Use STRUMPACKSolverBase::SetCompression to set the proper compression type.

Definition at line 295 of file strumpack.cpp.

◆ SetCompressionLossyPrecision()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetCompressionLossyPrecision ( int  precision)

Set the precision for the lossy compression option.

Use STRUMPACKSolverBase::SetCompression to set the proper compression type.

Definition at line 288 of file strumpack.cpp.

◆ SetCompressionRelTol()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetCompressionRelTol ( double  rtol)

Set the relative tolerance for low rank compression methods.

This currently affects BLR, HSS, and HODLR. Use STRUMPACKSolverBase::SetCompression to set the proper compression type.

Definition at line 263 of file strumpack.cpp.

◆ SetFromCommandLine()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetFromCommandLine

Set options that were captured from the command line.

These were captured in the constructer STRUMPACKSolverBase. Refer to the STRUMPACK documentation for details.

Definition at line 152 of file strumpack.cpp.

◆ SetKrylovSolver()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetKrylovSolver ( strumpack::KrylovSolver  method)

Set the Krylov solver method to use.

STRUMPACK is an (approximate) direct solver. It can be used as a direct solver or as a preconditioner. To use STRUMPACK as only a preconditioner, set the Krylov solver to DIRECT. STRUMPACK also provides iterative solvers which can use the preconditioner, and these iterative solvers can also be used without preconditioner.

Supported values are:

  • AUTO: Use iterative refinement if no HSS compression is used, otherwise use GMRES
  • DIRECT: No outer iterative solver, just a single application of the multifrontal solver
  • REFINE: Iterative refinement
  • PREC_GMRES: Preconditioned GMRES The preconditioner is the (approx) multifrontal solver
  • GMRES: UN-preconditioned GMRES (for testing mainly)
  • PREC_BICGSTAB: Preconditioned BiCGStab The preconditioner is the (approx) multifrontal solver
  • BICGSTAB: UN-preconditioned BiCGStab. (for testing mainly)

Definition at line 215 of file strumpack.cpp.

◆ SetMatching()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetMatching ( strumpack::MatchingJob  job)

Configure static pivoting for stability.

The static pivoting in STRUMPACK permutes the sparse input matrix in order to get large (nonzero) elements on the diagonal. If the input matrix is already diagonally dominant, this reordering can be disabled.

Supported matching algorithms are:

  • NONE: Don't do anything
  • MAX_CARDINALITY: Maximum cardinality
  • MAX_SMALLEST_DIAGONAL: Maximum smallest diagonal value
  • MAX_SMALLEST_DIAGONAL_2: Same as MAX_SMALLEST_DIAGONAL but different algorithm
  • MAX_DIAGONAL_SUM: Maximum sum of diagonal values
  • MAX_DIAGONAL_PRODUCT_SCALING: Maximum product of diagonal values and row and column scaling (default)
  • COMBBLAS: Use AWPM from CombBLAS (only with version >= 3)

Definition at line 229 of file strumpack.cpp.

◆ SetMaxIter()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetMaxIter ( int  max_it)

Set the maximum number of iterations for iterative solvers.

Definition at line 188 of file strumpack.cpp.

◆ SetOperator()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetOperator ( const Operator op)
virtual

Set the operator/matrix.

Note
A must be a STRUMPACKRowLocMatrix.

Implements mfem::Solver.

Definition at line 303 of file strumpack.cpp.

◆ SetPrintFactorStatistics()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetPrintFactorStatistics ( bool  print_stat)

Set up verbose printing during the factor step.

Definition at line 159 of file strumpack.cpp.

◆ SetPrintSolveStatistics()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetPrintSolveStatistics ( bool  print_stat)

Set up verbose printing during the solve step.

Definition at line 166 of file strumpack.cpp.

◆ SetRelTol()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetRelTol ( double  rtol)

Set the relative tolerance for interative solvers.

Definition at line 174 of file strumpack.cpp.

◆ SetReorderingReuse()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetReorderingReuse ( bool  reuse)

Set the flag controlling reuse of the symbolic factorization for multiple operators.

This method must be called before repeated calls to SetOperator.

Definition at line 195 of file strumpack.cpp.

◆ SetReorderingStrategy()

template<typename STRUMPACKSolverType >
void mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::SetReorderingStrategy ( strumpack::ReorderingStrategy  method)

Set matrix reordering strategy.

Supported reorderings are:

  • NATURAL: Do not reorder the system
  • METIS: Use Metis nested-dissection reordering (default)
  • PARMETIS: Use ParMetis nested-dissection reordering
  • SCOTCH: Use Scotch nested-dissection reordering
  • PTSCOTCH: Use PT-Scotch nested-dissection reordering
  • RCM: Use RCM reordering
  • GEOMETRIC: A simple geometric nested dissection code that only works for regular meshes
  • AMD: Approximate minimum degree
  • MMD: Multiple minimum degree
  • AND: Nested dissection
  • MLF: Minimum local fill
  • SPECTRAL: Spectral nested dissection

Definition at line 222 of file strumpack.cpp.

Member Data Documentation

◆ APtr_

template<typename STRUMPACKSolverType >
const STRUMPACKRowLocMatrix* mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::APtr_
protected

Definition at line 271 of file strumpack.hpp.

◆ factor_verbose_

template<typename STRUMPACKSolverType >
bool mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::factor_verbose_
protected

Definition at line 274 of file strumpack.hpp.

◆ nrhs_

template<typename STRUMPACKSolverType >
int mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::nrhs_
mutableprotected

Definition at line 279 of file strumpack.hpp.

◆ reorder_reuse_

template<typename STRUMPACKSolverType >
bool mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::reorder_reuse_
protected

Definition at line 276 of file strumpack.hpp.

◆ rhs_

template<typename STRUMPACKSolverType >
Vector mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::rhs_
mutableprotected

Definition at line 278 of file strumpack.hpp.

◆ sol_

template<typename STRUMPACKSolverType >
Vector mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::sol_
protected

Definition at line 278 of file strumpack.hpp.

◆ solve_verbose_

template<typename STRUMPACKSolverType >
bool mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::solve_verbose_
protected

Definition at line 275 of file strumpack.hpp.

◆ solver_

template<typename STRUMPACKSolverType >
STRUMPACKSolverType* mfem::STRUMPACKSolverBase< STRUMPACKSolverType >::solver_
protected

Definition at line 272 of file strumpack.hpp.


The documentation for this class was generated from the following files: