MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::HypreSmoother Class Reference

Parallel smoothers in hypre. More...

#include <hypre.hpp>

Inheritance diagram for mfem::HypreSmoother:
[legend]
Collaboration diagram for mfem::HypreSmoother:
[legend]

Public Types

enum  Type {
  Jacobi = 0 , l1Jacobi = 1 , l1GS = 2 , l1GStr = 4 ,
  lumpedJacobi = 5 , GS = 6 , OPFS = 10 , Chebyshev = 16 ,
  Taubin = 1001 , FIR = 1002
}
 HYPRE smoother types. More...
 
- 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 Member Functions

 HypreSmoother ()
 
 HypreSmoother (const HypreParMatrix &A_, int type=DefaultType(), int relax_times=1, real_t relax_weight=1.0, real_t omega=1.0, int poly_order=2, real_t poly_fraction=.3, int eig_est_cg_iter=10)
 
void SetType (HypreSmoother::Type type, int relax_times=1)
 Set the relaxation type and number of sweeps.
 
void SetSOROptions (real_t relax_weight, real_t omega)
 Set SOR-related parameters.
 
void SetPolyOptions (int poly_order, real_t poly_fraction, int eig_est_cg_iter=10)
 Set parameters for polynomial smoothing.
 
void SetTaubinOptions (real_t lambda, real_t mu, int iter)
 Set parameters for Taubin's lambda-mu method.
 
void SetWindowByName (const char *window_name)
 Convenience function for setting canonical windowing parameters.
 
void SetWindowParameters (real_t a, real_t b, real_t c)
 Set parameters for windowing function for FIR smoother.
 
void SetFIRCoefficients (real_t max_eig)
 Compute window and Chebyshev coefficients for given polynomial order.
 
void SetPositiveDiagonal (bool pos=true)
 After computing l1-norms, replace them with their absolute values.
 
void SetOperatorSymmetry (bool is_sym)
 
virtual void SetOperator (const Operator &op)
 
virtual void Mult (const HypreParVector &b, HypreParVector &x) const
 Relax the linear system Ax=b.
 
virtual void Mult (const Vector &b, Vector &x) const
 Operator application: y=A(x).
 
virtual void MultTranspose (const Vector &b, Vector &x) const
 Apply transpose of the smoother to relax the linear system Ax=b.
 
virtual ~HypreSmoother ()
 
- 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.
 
- 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 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.
 

Static Public Member Functions

static Type DefaultType ()
 Default value for the smoother type used by the constructors: Type::l1GS when HYPRE is running on CPU and Type::l1Jacobi when HYPRE is running on GPU.
 

Static Public Attributes

static MFEM_DEPRECATED constexpr Type default_type = l1GS
 

Protected Attributes

HypreParMatrixA
 The linear system matrix.
 
HypreParVectorB
 Right-hand side and solution vectors.
 
HypreParVectorX
 
Memory< real_tauxB
 Auxiliary buffers for the case when the input or output arrays in methods like Mult(const Vector &, Vector &) need to be deep copied in order to be used by hypre.
 
Memory< real_tauxX
 
HypreParVectorV
 Temporary vectors.
 
HypreParVectorZ
 
HypreParVectorX0
 FIR Filter Temporary Vectors.
 
HypreParVectorX1
 
int type
 
int relax_times
 Number of relaxation sweeps.
 
real_t relax_weight
 Damping coefficient (usually <= 1)
 
real_t omega
 SOR parameter (usually in (0,2))
 
int poly_order
 Order of the smoothing polynomial.
 
real_t poly_fraction
 Fraction of spectrum to smooth for polynomial relaxation.
 
int poly_scale
 Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
 
real_t lambda
 Taubin's lambda-mu method parameters.
 
real_t mu
 
int taubin_iter
 
real_tl1_norms
 l1 norms of the rows of A
 
bool pos_l1_norms
 If set, take absolute values of the computed l1_norms.
 
int eig_est_cg_iter
 Number of CG iterations to determine eigenvalue estimates.
 
real_t max_eig_est
 Maximal eigenvalue estimate for polynomial smoothing.
 
real_t min_eig_est
 Minimal eigenvalue estimate for polynomial smoothing.
 
real_t window_params [3]
 Parameters for windowing function of FIR filter.
 
real_tfir_coeffs
 Combined coefficients for windowing and Chebyshev polynomials.
 
bool A_is_symmetric
 A flag that indicates whether the linear system matrix A is symmetric.
 
- 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 Attributes inherited from mfem::Solver
bool iterative_mode
 If true, use the second argument of Mult() as an initial guess.
 
- 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".
 

Detailed Description

Parallel smoothers in hypre.

Definition at line 1019 of file hypre.hpp.

Member Enumeration Documentation

◆ Type

HYPRE smoother types.

Enumerator
Jacobi 

Jacobi.

l1Jacobi 

l1-scaled Jacobi

l1GS 

l1-scaled block Gauss-Seidel/SSOR

l1GStr 

truncated l1-scaled block Gauss-Seidel/SSOR

lumpedJacobi 

lumped Jacobi

GS 

Gauss-Seidel.

OPFS 

On-processor forward solve for matrix w/ triangular structure

Chebyshev 

Chebyshev.

Taubin 

Taubin polynomial smoother.

FIR 

FIR polynomial smoother.

Definition at line 1077 of file hypre.hpp.

Constructor & Destructor Documentation

◆ HypreSmoother() [1/2]

mfem::HypreSmoother::HypreSmoother ( )

Definition at line 3504 of file hypre.cpp.

◆ HypreSmoother() [2/2]

mfem::HypreSmoother::HypreSmoother ( const HypreParMatrix & A_,
int type = DefaultType(),
int relax_times = 1,
real_t relax_weight = 1.0,
real_t omega = 1.0,
int poly_order = 2,
real_t poly_fraction = .3,
int eig_est_cg_iter = 10 )

Definition at line 3526 of file hypre.cpp.

◆ ~HypreSmoother()

mfem::HypreSmoother::~HypreSmoother ( )
virtual

Definition at line 3928 of file hypre.cpp.

Member Function Documentation

◆ DefaultType()

static Type mfem::HypreSmoother::DefaultType ( )
inlinestatic

Default value for the smoother type used by the constructors: Type::l1GS when HYPRE is running on CPU and Type::l1Jacobi when HYPRE is running on GPU.

Definition at line 1102 of file hypre.hpp.

◆ Mult() [1/2]

void mfem::HypreSmoother::Mult ( const HypreParVector & b,
HypreParVector & x ) const
virtual

Relax the linear system Ax=b.

Definition at line 3777 of file hypre.cpp.

◆ Mult() [2/2]

void mfem::HypreSmoother::Mult ( const Vector & x,
Vector & y ) const
virtual

Operator application: y=A(x).

Implements mfem::Operator.

Definition at line 3857 of file hypre.cpp.

◆ MultTranspose()

void mfem::HypreSmoother::MultTranspose ( const Vector & b,
Vector & x ) const
virtual

Apply transpose of the smoother to relax the linear system Ax=b.

Reimplemented from mfem::Operator.

Definition at line 3918 of file hypre.cpp.

◆ SetFIRCoefficients()

void mfem::HypreSmoother::SetFIRCoefficients ( real_t max_eig)

Compute window and Chebyshev coefficients for given polynomial order.

Definition at line 3737 of file hypre.cpp.

◆ SetOperator()

void mfem::HypreSmoother::SetOperator ( const Operator & op)
virtual

Set/update the associated operator. Must be called after setting the HypreSmoother type and options.

Implements mfem::Solver.

Definition at line 3600 of file hypre.cpp.

◆ SetOperatorSymmetry()

void mfem::HypreSmoother::SetOperatorSymmetry ( bool is_sym)
inline

Explicitly indicate whether the linear system matrix A is symmetric. If A is symmetric, the smoother will also be symmetric. In this case, calling MultTranspose will be redirected to Mult. (This is also done if the smoother is diagonal.) By default, A is assumed to be nonsymmetric.

Definition at line 1142 of file hypre.hpp.

◆ SetPolyOptions()

void mfem::HypreSmoother::SetPolyOptions ( int poly_order,
real_t poly_fraction,
int eig_est_cg_iter = 10 )

Set parameters for polynomial smoothing.

By default, 10 iterations of CG are used to estimate the eigenvalues. Setting eig_est_cg_iter = 0 uses hypre's hypre_ParCSRMaxEigEstimate() instead.

Definition at line 3562 of file hypre.cpp.

◆ SetPositiveDiagonal()

void mfem::HypreSmoother::SetPositiveDiagonal ( bool pos = true)
inline

After computing l1-norms, replace them with their absolute values.

By default, the l1-norms take their sign from the corresponding diagonal entries in the associated matrix.

Definition at line 1136 of file hypre.hpp.

◆ SetSOROptions()

void mfem::HypreSmoother::SetSOROptions ( real_t relax_weight,
real_t omega )

Set SOR-related parameters.

Definition at line 3556 of file hypre.cpp.

◆ SetTaubinOptions()

void mfem::HypreSmoother::SetTaubinOptions ( real_t lambda,
real_t mu,
int iter )

Set parameters for Taubin's lambda-mu method.

Definition at line 3570 of file hypre.cpp.

◆ SetType()

void mfem::HypreSmoother::SetType ( HypreSmoother::Type type,
int relax_times = 1 )

Set the relaxation type and number of sweeps.

Definition at line 3550 of file hypre.cpp.

◆ SetWindowByName()

void mfem::HypreSmoother::SetWindowByName ( const char * window_name)

Convenience function for setting canonical windowing parameters.

Definition at line 3578 of file hypre.cpp.

◆ SetWindowParameters()

void mfem::HypreSmoother::SetWindowParameters ( real_t a,
real_t b,
real_t c )

Set parameters for windowing function for FIR smoother.

Definition at line 3593 of file hypre.cpp.

Member Data Documentation

◆ A

HypreParMatrix* mfem::HypreSmoother::A
protected

The linear system matrix.

Definition at line 1023 of file hypre.hpp.

◆ A_is_symmetric

bool mfem::HypreSmoother::A_is_symmetric
protected

A flag that indicates whether the linear system matrix A is symmetric.

Definition at line 1073 of file hypre.hpp.

◆ auxB

Memory<real_t> mfem::HypreSmoother::auxB
mutableprotected

Auxiliary buffers for the case when the input or output arrays in methods like Mult(const Vector &, Vector &) need to be deep copied in order to be used by hypre.

Definition at line 1029 of file hypre.hpp.

◆ auxX

Memory<real_t> mfem::HypreSmoother::auxX
protected

Definition at line 1029 of file hypre.hpp.

◆ B

HypreParVector* mfem::HypreSmoother::B
mutableprotected

Right-hand side and solution vectors.

Definition at line 1025 of file hypre.hpp.

◆ default_type

static MFEM_DEPRECATED constexpr Type mfem::HypreSmoother::default_type = l1GS
staticconstexpr
Deprecated
Use DefaultType() instead

Definition at line 1094 of file hypre.hpp.

◆ eig_est_cg_iter

int mfem::HypreSmoother::eig_est_cg_iter
protected

Number of CG iterations to determine eigenvalue estimates.

Definition at line 1061 of file hypre.hpp.

◆ fir_coeffs

real_t* mfem::HypreSmoother::fir_coeffs
protected

Combined coefficients for windowing and Chebyshev polynomials.

Definition at line 1070 of file hypre.hpp.

◆ l1_norms

real_t* mfem::HypreSmoother::l1_norms
protected

l1 norms of the rows of A

Definition at line 1057 of file hypre.hpp.

◆ lambda

real_t mfem::HypreSmoother::lambda
protected

Taubin's lambda-mu method parameters.

Definition at line 1052 of file hypre.hpp.

◆ max_eig_est

real_t mfem::HypreSmoother::max_eig_est
protected

Maximal eigenvalue estimate for polynomial smoothing.

Definition at line 1063 of file hypre.hpp.

◆ min_eig_est

real_t mfem::HypreSmoother::min_eig_est
protected

Minimal eigenvalue estimate for polynomial smoothing.

Definition at line 1065 of file hypre.hpp.

◆ mu

real_t mfem::HypreSmoother::mu
protected

Definition at line 1053 of file hypre.hpp.

◆ omega

real_t mfem::HypreSmoother::omega
protected

SOR parameter (usually in (0,2))

Definition at line 1043 of file hypre.hpp.

◆ poly_fraction

real_t mfem::HypreSmoother::poly_fraction
protected

Fraction of spectrum to smooth for polynomial relaxation.

Definition at line 1047 of file hypre.hpp.

◆ poly_order

int mfem::HypreSmoother::poly_order
protected

Order of the smoothing polynomial.

Definition at line 1045 of file hypre.hpp.

◆ poly_scale

int mfem::HypreSmoother::poly_scale
protected

Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.

Definition at line 1049 of file hypre.hpp.

◆ pos_l1_norms

bool mfem::HypreSmoother::pos_l1_norms
protected

If set, take absolute values of the computed l1_norms.

Definition at line 1059 of file hypre.hpp.

◆ relax_times

int mfem::HypreSmoother::relax_times
protected

Number of relaxation sweeps.

Definition at line 1039 of file hypre.hpp.

◆ relax_weight

real_t mfem::HypreSmoother::relax_weight
protected

Damping coefficient (usually <= 1)

Definition at line 1041 of file hypre.hpp.

◆ taubin_iter

int mfem::HypreSmoother::taubin_iter
protected

Definition at line 1054 of file hypre.hpp.

◆ type

int mfem::HypreSmoother::type
protected

Smoother type from hypre_ParCSRRelax() in ams.c plus extensions, see the enumeration Type below.

Definition at line 1037 of file hypre.hpp.

◆ V

HypreParVector* mfem::HypreSmoother::V
mutableprotected

Temporary vectors.

Definition at line 1031 of file hypre.hpp.

◆ window_params

real_t mfem::HypreSmoother::window_params[3]
protected

Parameters for windowing function of FIR filter.

Definition at line 1067 of file hypre.hpp.

◆ X

HypreParVector * mfem::HypreSmoother::X
protected

Definition at line 1025 of file hypre.hpp.

◆ X0

HypreParVector* mfem::HypreSmoother::X0
mutableprotected

FIR Filter Temporary Vectors.

Definition at line 1033 of file hypre.hpp.

◆ X1

HypreParVector * mfem::HypreSmoother::X1
protected

Definition at line 1033 of file hypre.hpp.

◆ Z

HypreParVector * mfem::HypreSmoother::Z
protected

Definition at line 1031 of file hypre.hpp.


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