G+Smo
24.08.0
Geometry + Simulation Modules
|
Multigrid preconditioner.
This class implements an generic geometric multigrid framework. This class implements a preconditioner which can be used, e.g., in an iterative solver, like gsConjugateGradient. If a simple multigrid iteration is required, use this class as a preconditioner in gsGradientMethod with step size 1.
The class expects transfer matrices, which can be constructed with the class gsGridHierarchy.
The intergrid transfers (restriction and prolongation) and the stiffness matrix are stored as gsLinearOperator, which means that also matrix-free variants are possible. If linear operators are provided, operators representing the siffness matrix have to be provided for all grid levles. If the integrid transfers and the stiffness matrix are provided as a sparse matrix, the matrices for the coarser levels are computed automatically based on the Galerkin principle.
For all levels, a smoother has to be provided by defining setSmoother().
The solver for the coarsest grid level is created automatically if not provided by the caller.
The solver can be configured to use V- or W-cycles using setNumCycles(). The number of pre- and post-smoothing steps can be configured with setNumPreSmooth() and setNumPostSmooth(), respectively.
Also full multigrid and cascadic multigrid algorithms for computing an initial guess are provided.
Public Types | |
typedef gsPreconditionerOp< T > | Base |
Direct base class. | |
typedef gsLinearOperator< T >::Ptr | BasePtr |
Base class. | |
typedef gsMatrix< T > | Matrix |
Matrix type. | |
typedef gsLinearOperator< T >::Ptr | OpPtr |
Shared pointer to gsLinearOperator. | |
typedef gsPreconditionerOp< T > ::Ptr | PrecondPtr |
Shared pointer to gsPreconditionerOp. | |
typedef memory::shared_ptr < gsMultiGridOp > | Ptr |
Shared pointer for gsMultiGridOp. | |
typedef gsSparseMatrix< T > | SpMatrix |
Sparse matrix type. | |
typedef memory::shared_ptr < SpMatrix > | SpMatrixPtr |
Smart pointer to sparse matrix type. | |
typedef gsSparseMatrix< T, RowMajor > | SpMatrixRowMajor |
Matrix type for transfers. | |
typedef memory::shared_ptr < SpMatrixRowMajor > | SpMatrixRowMajorPtr |
Smart pointer to matrix type for transfers. | |
typedef memory::unique_ptr < gsMultiGridOp > | uPtr |
Unique pointer for gsMultiGridOp. | |
Public Member Functions | |
void | apply (const gsMatrix< T > &input, gsMatrix< T > &x) const |
apply the operator on the input vector and store the result in x More... | |
void | cascadicMultiGrid (const std::vector< Matrix > &rhs, const std::vector< Matrix > &fixedValues, Matrix &result) const |
Perform one cascadic multigrid cycle and store the resulting solution vector in result. More... | |
T | coarseGridCorrectionDamping () const |
Get the damping of for the coarse-grid correction. | |
OpPtr | coarseSolver () const |
Get the solver for the coarsest problem (level 0) | |
index_t | cols () const override |
Returns the number of columns of the operator. | |
T | estimateLargestEigenvalueOfPreconditionedSystem (index_t steps=10) const |
Estimates the largest eigenvalue of \( PA \). More... | |
index_t | finestLevel () const |
The index of the finest level (0-based) | |
void | fullMultiGrid (const std::vector< Matrix > &rhs, const std::vector< Matrix > &dirichletIntp, Matrix &result) const |
Perform one full multigrid cycle and store the resulting solution vector in result. More... | |
gsMultiGridOp (SpMatrix fineMatrix, std::vector< SpMatrixRowMajor > transferMatrices, OpPtr coarseSolver=OpPtr()) | |
Constructor. More... | |
gsMultiGridOp (SpMatrixPtr fineMatrix, std::vector< SpMatrixRowMajorPtr > transferMatrices, OpPtr coarseSolver=OpPtr()) | |
Constructor. More... | |
gsMultiGridOp (const std::vector< OpPtr > &ops, const std::vector< OpPtr > &prolongation, const std::vector< OpPtr > &restriction, OpPtr coarseSolver=OpPtr()) | |
Constructor for a matix-free variant. More... | |
const SpMatrix & | matrix (index_t lvl) const |
Stiffness matrix for given level. | |
const SpMatrix & | matrix () const |
Stiffness matrix for finest level. | |
void | multiGridStep (index_t level, const Matrix &rhs, Matrix &x) const |
Perform one multigrid cycle, starting from the given level. More... | |
index_t | nDofs (index_t lvl) const |
Number of dofs for the given level. | |
index_t | nDofs () const |
Number of dofs for the finest level. | |
index_t | numCycles (index_t lvl) const |
Get number of coarse grid steps to be applied. | |
index_t | numLevels () const |
Number of levels in the multigrid construction. | |
index_t | numOfSweeps () |
Get the number of sweeps to be applied in the member function apply. | |
void | prolongVector (index_t lc, const Matrix &coarse, Matrix &fine) const |
Prolong a vector of coefficients from the level lc, given by coarse, to the next finer level. | |
void | restrictVector (index_t lf, const Matrix &fine, Matrix &coarse) const |
Restrict a vector of coefficients from the level lf, given by fine, to the next coarser level. | |
index_t | rows () const override |
Returns the number of rows of the operator. | |
void | setCoarseGridCorrectionDamping (T damping) |
Set the damping of for the coarse-grid correction. | |
void | setCoarseSolver (const OpPtr &sol) |
Set the solver for the coarsest problem (level 0) | |
void | setNumCycles (index_t n) |
void | setNumCycles (index_t lvl, index_t n) |
void | setNumOfSweeps (index_t n) |
Set the number of sweeps to be applied in the member function apply. | |
void | setNumPostSmooth (index_t n) |
Set number of post-smoothing steps to perform. | |
index_t | setNumPostSmooth () const |
Get number of post-smoothing steps to perform. | |
void | setNumPreSmooth (index_t n) |
Set number of pre-smoothing steps to perform. | |
index_t | setNumPreSmooth () const |
Get number of pre-smoothing steps to perform. | |
void | setOptions (const gsOptionList &opt) override |
Set the options based on a gsOptionList. | |
void | setSmoother (index_t lvl, const PrecondPtr &sm) |
void | setSymmSmooth (bool s) |
void | setUnderlyingOp (index_t lvl, OpPtr op) |
Set underlying operator (=stiffness matrix) for certain level. | |
PrecondPtr | smoother (index_t lvl) const |
void | smoothingStep (index_t level, const Matrix &rhs, Matrix &x) const |
Apply smoothing steps on the corresponding level. More... | |
void | smoothingStep (const Matrix &rhs, Matrix &x) const |
Apply smoothing steps on finest grid. More... | |
void | solveCoarse (const Matrix &rhs, Matrix &result) const |
Solve the problem with direct solver on coarsest level. | |
void | step (const Matrix &rhs, Matrix &x) const override |
Perform one multigrid cycle. More... | |
void | stepT (const Matrix &rhs, Matrix &x) const override |
Perform one transposed multigrid cycle. More... | |
bool | symmSmooth () const |
Get symmetric post-smoothing option. | |
OpPtr | underlyingOp (index_t lvl) const |
Underlying operator (=stiffness matrix) for given level. | |
OpPtr | underlyingOp () const override |
Underlying operator (=stiffness matrix) for finest level. | |
Static Public Member Functions | |
static gsOptionList | defaultOptions () |
Returns a list of default options. | |
static gsIdentityOp< T > | Identity (const index_t dim) |
Identity operator. | |
static uPtr | make (SpMatrix fineMatrix, std::vector< SpMatrixRowMajor > transferMatrices, OpPtr coarseSolver=OpPtr()) |
Make function returning smart pointer. More... | |
static uPtr | make (SpMatrixPtr fineMatrix, std::vector< SpMatrixRowMajorPtr > transferMatrices, OpPtr coarseSolver=OpPtr()) |
Make function returning smart pointer. More... | |
static uPtr | make (const std::vector< OpPtr > &ops, const std::vector< OpPtr > &prolongation, const std::vector< OpPtr > &restriction, OpPtr coarseSolver=OpPtr()) |
Make function returning a shared pointer for a matix-free variant. More... | |
Private Member Functions | |
void | init (SpMatrixPtr fineMatrix, std::vector< SpMatrixRowMajorPtr > transferMatrices) |
Init function that is used by matrix based constructors. | |
void | initCoarseSolver () const |
Init solver on coarsest grid level. | |
Private Attributes | |
OpPtr | m_coarseSolver |
Solver for the coarsest-grid problem. | |
T | m_damping |
Damping for the coarse-grid correction. | |
index_t | m_nLevels |
Number of levels. | |
gsVector< index_t > | m_numCycles |
Number of coarse-grid steps (1=V-cycle, 2=W-cycle) | |
index_t | m_numPostSmooth |
Number of post-smoothing steps. | |
index_t | m_numPreSmooth |
Number of pre-smoothing steps. | |
std::vector< OpPtr > | m_ops |
Underlying operators (=stiffness matrix) on each level. | |
std::vector< OpPtr > | m_prolong |
Prolongation matrix for each grid transition. | |
std::vector< OpPtr > | m_restrict |
Restriction matrix for each grid transition. | |
std::vector< PrecondPtr > | m_smoother |
Smoothers on each level. | |
bool | m_symmSmooth |
Symmmetric smoothing. | |
gsMultiGridOp | ( | SpMatrix | fineMatrix, |
std::vector< SpMatrixRowMajor > | transferMatrices, | ||
OpPtr | coarseSolver = OpPtr() |
||
) |
Constructor.
fineMatrix | Stiffness matrix on the finest grid |
transferMatrices | Intergrid transfer matrices representing restriction and prolongation operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
gsMultiGridOp | ( | SpMatrixPtr | fineMatrix, |
std::vector< SpMatrixRowMajorPtr > | transferMatrices, | ||
OpPtr | coarseSolver = OpPtr() |
||
) |
Constructor.
fineMatrix | Stiffness matrix (as smart pointers) on the finest grid |
transferMatrices | Intergrid transfer matrices representing restriction and prolongation operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
gsMultiGridOp | ( | const std::vector< OpPtr > & | ops, |
const std::vector< OpPtr > & | prolongation, | ||
const std::vector< OpPtr > & | restriction, | ||
OpPtr | coarseSolver = OpPtr() |
||
) |
Constructor for a matix-free variant.
ops | Linear operators representing the stiffness matrix on all levels |
prolongation | Linear operators representing the prolongation operators |
restriction | Linear operators representing the restriction operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
apply the operator on the input vector and store the result in x
input | Input vector |
x | result vector |
Implements gsLinearOperator< T >.
Reimplemented in gsPreconditionerFromOp< T >.
void cascadicMultiGrid | ( | const std::vector< Matrix > & | rhs, |
const std::vector< Matrix > & | fixedValues, | ||
Matrix & | result | ||
) | const |
Perform one cascadic multigrid cycle and store the resulting solution vector in result.
[in] | rhs | The right-hand-sides for all grid levels |
[in] | fixedValues | The values for the fixed dofs (like eliminated Dirichlet values) |
[out] | result | The destination for the result |
|
inlineinherited |
Estimates the largest eigenvalue of \( PA \).
steps | Number of steps to be performed. |
void fullMultiGrid | ( | const std::vector< Matrix > & | rhs, |
const std::vector< Matrix > & | dirichletIntp, | ||
Matrix & | result | ||
) | const |
Perform one full multigrid cycle and store the resulting solution vector in result.
[in] | rhs | The right-hand-sides for all grid levels |
[in] | fixedValues | The values for the fixed dofs (like eliminated Dirichlet values) |
[out] | result | The destination for the result |
|
inlinestatic |
Make function returning smart pointer.
fineMatrix | Stiffness matrix on the finest grid |
transferMatrices | Intergrid transfer matrices representing restriction and prolongation operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
|
inlinestatic |
Make function returning smart pointer.
fineMatrix | Stiffness matrix (as smart pointers) on the finest grid |
transferMatrices | Intergrid transfer matrices representing restriction and prolongation operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
|
inlinestatic |
Make function returning a shared pointer for a matix-free variant.
ops | Linear operators representing the stiffness matrix on all levels |
prolongation | Linear operators representing the prolongation operators |
restriction | Linear operators representing the restriction operators |
coarseSolver | Linear operator representing the exact solver on the coarsest grid level, defaulted to a direct solver (PartialPivLUSolver) |
Perform one multigrid cycle, starting from the given level.
[in] | level | Level |
[in] | rhs | Right-hand side |
[in,out] | x | Current iterate vector |
|
inline |
Set number of coarse grid steps to be applied
n | Number of levels (1 for V-cycle, 2 for W-cycle) |
Set number of coarse grid steps to be applied
lvl | Level of coarser grid (multiGridStep on level lvl+1 calls n times multiGridStep for level lvl) |
n | Number of calls (1 for V-cycle, 2 for W-cycle) |
void setSmoother | ( | index_t | lvl, |
const PrecondPtr & | sm | ||
) |
Set the smoother
lvl | The corresponding level |
sm | The smoother |
|
inline |
Set symmetric post-smoothing option
s | Iff true (default), stepT is called for post-smoothing |
|
inline |
Get the smoother
lvl | The corresponding level |
Apply smoothing steps on the corresponding level.
[in] | level | Level |
[in] | rhs | Right hand side |
[in,out] | x | Current iterate vector |
Apply smoothing steps on finest grid.
[in] | rhs | Right hand side |
[in,out] | x | Current iterate vector |
Perform one multigrid cycle.
[in] | rhs | Right-hand side |
[in,out] | x | Current iterate vector |
Implements gsPreconditionerOp< T >.
Perform one transposed multigrid cycle.
[in] | rhs | Right-hand side |
[in,out] | x | Current iterate vector |
This assumes that the stiffness matrix is symmetric, that SymmSmooth is true and that the smoothers implement stepT correctly.
Reimplemented from gsPreconditionerOp< T >.