MFEM
v4.3.0
Finite element discretization library
|
Multigrid solver class. More...
#include <multigrid.hpp>
Public Types | |
enum | CycleType { CycleType::VCYCLE, CycleType::WCYCLE } |
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 } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Public Member Functions | |
Multigrid () | |
Constructs an empty multigrid hierarchy. More... | |
Multigrid (const Array< Operator * > &operators_, const Array< Solver * > &smoothers_, const Array< Operator * > &prolongations_, const Array< bool > &ownedOperators_, const Array< bool > &ownedSmoothers_, const Array< bool > &ownedProlongations_) | |
Constructs a multigrid hierarchy from the given inputs. More... | |
virtual | ~Multigrid () |
Destructor. More... | |
void | AddLevel (Operator *opr, Solver *smoother, bool ownOperator, bool ownSmoother) |
Adds a level to the multigrid operator hierarchy. More... | |
int | NumLevels () const |
Returns the number of levels. More... | |
int | GetFinestLevelIndex () const |
Returns the index of the finest level. More... | |
const Operator * | GetOperatorAtLevel (int level) const |
Returns operator at given level. More... | |
Operator * | GetOperatorAtLevel (int level) |
Returns operator at given level. More... | |
const Operator * | GetOperatorAtFinestLevel () const |
Returns operator at finest level. More... | |
Operator * | GetOperatorAtFinestLevel () |
Returns operator at finest level. More... | |
Solver * | GetSmootherAtLevel (int level) const |
Returns smoother at given level. More... | |
Solver * | GetSmootherAtLevel (int level) |
Returns smoother at given level. More... | |
void | SetCycleType (CycleType cycleType_, int preSmoothingSteps_, int postSmoothingSteps_) |
Set cycle type and number of pre- and post-smoothing steps used by Mult. More... | |
virtual void | Mult (const Vector &x, Vector &y) const override |
Application of the multigrid as a preconditioner. More... | |
virtual void | SetOperator (const Operator &op) override |
Not supported for multigrid. More... | |
Public Member Functions inherited from mfem::Solver | |
Solver (int s=0, bool iter_mode=false) | |
Initialize a square Solver with size s. More... | |
Solver (int h, int w, bool iter_mode=false) | |
Initialize a Solver with height h and width w. More... | |
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. More... | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). More... | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). More... | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). More... | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). More... | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). More... | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). More... | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. More... | |
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. More... | |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More... | |
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. More... | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More... | |
virtual const Operator * | GetOutputRestriction () const |
Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
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. More... | |
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. More... | |
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(). More... | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. More... | |
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). More... | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. More... | |
void | PrintMatlab (std::ostream &out, int n=0, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Type | GetType () const |
Return the type ID of the Operator class. More... | |
Protected Attributes | |
Array< Operator * > | operators |
Array< Solver * > | smoothers |
Array< Operator * > | prolongations |
Array< bool > | ownedOperators |
Array< bool > | ownedSmoothers |
Array< bool > | ownedProlongations |
CycleType | cycleType |
int | preSmoothingSteps |
int | postSmoothingSteps |
Array< Vector * > | X |
Array< Vector * > | Y |
Array< Vector * > | R |
Array< Vector * > | Z |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. More... | |
int | width |
Dimension of the input / number of columns in the matrix. More... | |
Additional Inherited Members | |
Public Attributes inherited from mfem::Solver | |
bool | iterative_mode |
If true, use the second argument of Mult() as an initial guess. More... | |
Protected Member Functions inherited from mfem::Operator | |
void | FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout) |
see FormSystemOperator() More... | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() More... | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". More... | |
Multigrid solver class.
Definition at line 25 of file multigrid.hpp.
|
strong |
Enumerator | |
---|---|
VCYCLE | |
WCYCLE |
Definition at line 28 of file multigrid.hpp.
mfem::Multigrid::Multigrid | ( | ) |
Constructs an empty multigrid hierarchy.
Definition at line 17 of file multigrid.cpp.
mfem::Multigrid::Multigrid | ( | const Array< Operator * > & | operators_, |
const Array< Solver * > & | smoothers_, | ||
const Array< Operator * > & | prolongations_, | ||
const Array< bool > & | ownedOperators_, | ||
const Array< bool > & | ownedSmoothers_, | ||
const Array< bool > & | ownedProlongations_ | ||
) |
Constructs a multigrid hierarchy from the given inputs.
Inputs include operators and smoothers on all levels, prolongation operators that go from coarser to finer levels, and ownership of the given operators, smoothers, and prolongations.
Definition at line 21 of file multigrid.cpp.
|
virtual |
Destructor.
Definition at line 51 of file multigrid.cpp.
void mfem::Multigrid::AddLevel | ( | Operator * | opr, |
Solver * | smoother, | ||
bool | ownOperator, | ||
bool | ownSmoother | ||
) |
Adds a level to the multigrid operator hierarchy.
The ownership of the operators and solvers/smoothers may be transferred to the Multigrid by setting the according boolean variables.
Definition at line 86 of file multigrid.cpp.
int mfem::Multigrid::GetFinestLevelIndex | ( | ) | const |
Returns the index of the finest level.
Definition at line 108 of file multigrid.cpp.
const Operator * mfem::Multigrid::GetOperatorAtFinestLevel | ( | ) | const |
Returns operator at finest level.
Definition at line 120 of file multigrid.cpp.
Operator * mfem::Multigrid::GetOperatorAtFinestLevel | ( | ) |
Returns operator at finest level.
Definition at line 125 of file multigrid.cpp.
const Operator * mfem::Multigrid::GetOperatorAtLevel | ( | int | level | ) | const |
Returns operator at given level.
Definition at line 110 of file multigrid.cpp.
Operator * mfem::Multigrid::GetOperatorAtLevel | ( | int | level | ) |
Returns operator at given level.
Definition at line 115 of file multigrid.cpp.
Solver * mfem::Multigrid::GetSmootherAtLevel | ( | int | level | ) | const |
Returns smoother at given level.
Definition at line 130 of file multigrid.cpp.
Solver * mfem::Multigrid::GetSmootherAtLevel | ( | int | level | ) |
Returns smoother at given level.
Definition at line 135 of file multigrid.cpp.
Application of the multigrid as a preconditioner.
Implements mfem::Operator.
Definition at line 148 of file multigrid.cpp.
int mfem::Multigrid::NumLevels | ( | ) | const |
Returns the number of levels.
Definition at line 106 of file multigrid.cpp.
void mfem::Multigrid::SetCycleType | ( | CycleType | cycleType_, |
int | preSmoothingSteps_, | ||
int | postSmoothingSteps_ | ||
) |
Set cycle type and number of pre- and post-smoothing steps used by Mult.
Definition at line 140 of file multigrid.cpp.
|
overridevirtual |
|
protected |
Definition at line 43 of file multigrid.hpp.
Definition at line 35 of file multigrid.hpp.
|
protected |
Definition at line 39 of file multigrid.hpp.
|
protected |
Definition at line 41 of file multigrid.hpp.
|
protected |
Definition at line 40 of file multigrid.hpp.
|
protected |
Definition at line 45 of file multigrid.hpp.
|
protected |
Definition at line 44 of file multigrid.hpp.
Definition at line 37 of file multigrid.hpp.
Definition at line 49 of file multigrid.hpp.
Definition at line 36 of file multigrid.hpp.
Definition at line 47 of file multigrid.hpp.
Definition at line 48 of file multigrid.hpp.
Definition at line 50 of file multigrid.hpp.