KASKADE 7 development version
Public Types | Public Member Functions | Related Functions | List of all members
Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation > Class Template Reference

A general multiplicative multigrid preconditioner. More...

#include <multiplicativeMultigrid.hh>

Detailed Description

template<class Entry, class Index, class Smoother, class Prolongation>
class Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >

A general multiplicative multigrid preconditioner.

This realizes a classical V-cycle with a specified number of pre- and post-smoothings. Two (possibly different) preconditioners are used for (i) smoothing on the levels and (ii) approximately solving the coarse grid system.

Note that this is only a preconditioner, no solver. In fact, the computed step is not guaranteed to decrease the energy of the error (depending on the smoother - Gauss-Seidel does, but Jacobi need not).

Template Parameters
Entrythe type of projected Galerkin matrix entries
Indexthe index type of Galerkin matrices (usually size_t)
Smoothera smoother type
Prolongationthe type of prolongations (conceptually a sparse matrix type)

The smoother may rely on the matrix provided on construction remaining available during any method calls besides the smoother destruction. Hence it is perfectly suitable for the smoother just to reference parts of the matrix.

A convenient construction of multiplicative multigrid preconditioners is provided by the functions makeMultiplicativeMultiGrid and makeJacobiMultiGrid for h-multigrid in grid hierarchies, and makeMultiplicativePMultiGrid makeJacobiPMultiGrid for p-h-multigrid exploiting both polynomial order and grid hierarchies.

Multigrid preconditioners can be put together easily from the three main ingredients multigrid stack, coarse solver, and smoother as follows:

auto mgStack = makeGeometricMultiGridStack(gridman,duplicate(A),onlyLowerTriangle);
auto coarsePreconditioner = makeDirectPreconditioner(std::move(mgStack.coarseGridMatrix()));
auto mg = makeMultiplicativeMultiGrid(std::move(mgStack),MakeJacobiSmoother(),moveUnique(std::move(coarsePreconditioner)),nPre,nPost);
mg.setSmootherStepSize(0.5);
Functor for creating Jacobi smoothers.
auto makeDirectPreconditioner(Matrix &&A, DirectType directType=DirectType::MUMPS)
Creates a direct solver for the given matrix.
auto makeMultiplicativeMultiGrid(MultiGridStack< Prolongation, Entry, Index > &&mgStack, MakeSmoother const &makeSmoother, std::unique_ptr< CoarsePreconditioner > &&coarsePreconditioner, int nPre=3, int nPost=3)
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for P1 elements...
MultiGridStack< MGProlongation, Entry, Index > makeGeometricMultiGridStack(GridMan const &gridManager, NumaBCRSMatrix< Entry, Index > &&A, size_t minNodes=10000, bool onlyLowerTriangle=false)
convenience routine for creating multigrid stacks based on geometric coarsening for P1 elements
auto moveUnique(A &&a)
Creates a dynamically allocated object from an r-value reference.
Definition: memory.hh:33
A duplicate(A const &a)
Creates a copy of a given object.
Definition: memory.hh:48

Definition at line 73 of file multiplicativeMultigrid.hh.

Inheritance diagram for Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >:
Kaskade::SymmetricPreconditioner< Smoother::domain_type, Smoother::range_type >

Public Types

using field_type = typename Smoother::field_type
 
using domain_type = typename Smoother::domain_type
 
using range_type = typename Smoother::range_type
 
using CoarsePreconditioner = Dune::Preconditioner< domain_type, range_type >
 

Public Member Functions

 MultiplicativeMultiGrid ()=default
 Default constructor. More...
 
template<class MakeSmoother >
 MultiplicativeMultiGrid (MultiGridStack< Prolongation, Entry, Index > &&mgStack_, MakeSmoother const &makeSmoother, std::unique_ptr< CoarsePreconditioner > &&coarsePreconditioner_, int nPre_=3, int nPost_=3)
 Constructor. More...
 
 MultiplicativeMultiGrid (MultiplicativeMultiGrid &&other)=default
 
MultiplicativeMultiGridoperator= (MultiplicativeMultiGrid &&other)=default
 
virtual void apply (domain_type &x, range_type const &r)
 Application of preconditioner. More...
 
virtual field_type applyDp (domain_type &x, range_type const &r)
 Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \). More...
 
virtual bool requiresInitializedInput () const
 Returns true if the target vector x has to be initialized to zero before calling apply or applyDp. More...
 
void setSmoothings (int nPre_, int nPost_)
 Sets the number of pre- and post-smoothing iterations to perform. More...
 
void setSmootherStepSize (double w)
 Define the smoother step size. More...
 
double getSmootherStepSize () const
 
void setLinesearch (bool ls)
 Enable or disable the line search option. More...
 
CoarsePreconditionergetCoarsePreconditioner ()
 Provides access to the coarse grid preconditioner. More...
 
std::pair< int, int > getSmoothings ()
 Returns a pair of number of pre- and post-smoothing iterations to perform. More...
 
virtual void pre (Smoother::domain_type &, Smoother::range_type &)
 Preconditioner preparation. More...
 
virtual void post (Smoother::domain_type &x)
 Preconditioner cleanup. More...
 
virtual Dune::SolverCategory::Category category () const override
 returns the category of the operator More...
 

Related Functions

(Note that these are not member functions.)

template<class Entry , class Index , class Prolongation , class MakeSmoother , class CoarsePreconditioner >
auto makeMultiplicativeMultiGrid (MultiGridStack< Prolongation, Entry, Index > &&mgStack, MakeSmoother const &makeSmoother, std::unique_ptr< CoarsePreconditioner > &&coarsePreconditioner, int nPre=3, int nPost=3)
 Convenience function creating multiplicative multigrid preconditioner of V-cycle type for P1 elements. More...
 
template<class Entry , class Index , class GridMan >
auto makeJacobiMultiGrid (NumaBCRSMatrix< Entry, Index > const &A, GridMan const &gridman, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS)
 Convenience function creating multiplicative Jacobi multigrid for P1 elements. More...
 
template<class Entry , class Index , class Space , class MakeSmoother >
auto makeMultiplicativePMultiGrid (NumaBCRSMatrix< Entry, Index > &&A, Space const &space, MakeSmoother const &makeSmoother, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS)
 Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More...
 
template<class Entry , class Index , class Space , class CoarseSpace >
auto makeBlockJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, Space const &space, NumaBCRSMatrix< Entry, Index > coarseA, CoarseSpace const &coarseSpace, int nPre=3, int nPost=3, DirectType directType=DirectType::MUMPS)
 Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More...
 
template<class Entry , class Index , class Space >
auto makeBlockJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, Space const &space, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS)
 Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More...
 
template<class Entry , class Index , class FineSpace >
auto makeJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, FineSpace const &space, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS)
 Convenience function creating multiplicative multigrid preconditioner of V-cycle type with Jacobi smoother for higher order elements. More...
 

Member Typedef Documentation

◆ CoarsePreconditioner

template<class Entry , class Index , class Smoother , class Prolongation >
using Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::CoarsePreconditioner = Dune::Preconditioner<domain_type,range_type>

Definition at line 79 of file multiplicativeMultigrid.hh.

◆ domain_type

template<class Entry , class Index , class Smoother , class Prolongation >
using Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::domain_type = typename Smoother::domain_type

Definition at line 77 of file multiplicativeMultigrid.hh.

◆ field_type

template<class Entry , class Index , class Smoother , class Prolongation >
using Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::field_type = typename Smoother::field_type

Definition at line 76 of file multiplicativeMultigrid.hh.

◆ range_type

template<class Entry , class Index , class Smoother , class Prolongation >
using Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::range_type = typename Smoother::range_type

Definition at line 78 of file multiplicativeMultigrid.hh.

Constructor & Destructor Documentation

◆ MultiplicativeMultiGrid() [1/3]

template<class Entry , class Index , class Smoother , class Prolongation >
Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::MultiplicativeMultiGrid ( )
default

Default constructor.

◆ MultiplicativeMultiGrid() [2/3]

template<class Entry , class Index , class Smoother , class Prolongation >
template<class MakeSmoother >
Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::MultiplicativeMultiGrid ( MultiGridStack< Prolongation, Entry, Index > &&  mgStack_,
MakeSmoother const &  makeSmoother,
std::unique_ptr< CoarsePreconditioner > &&  coarsePreconditioner_,
int  nPre_ = 3,
int  nPost_ = 3 
)
inline

Constructor.

Parameters
Athe Galerkin matrix to be preconditioned
Psthe stack of multigrid prolongations
makeSmoothera callable object with arguments (Matrix const&, int level) that creates a smoother of type Smoother
See also
makeMultiplicativeMultiGrid
makeJacobiMultiGrid

Definition at line 96 of file multiplicativeMultigrid.hh.

◆ MultiplicativeMultiGrid() [3/3]

template<class Entry , class Index , class Smoother , class Prolongation >
Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::MultiplicativeMultiGrid ( MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation > &&  other)
default

Member Function Documentation

◆ apply()

template<class Entry , class Index , class Smoother , class Prolongation >
virtual void Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::apply ( domain_type x,
range_type const &  r 
)
inlinevirtual

Application of preconditioner.

Precondition: x = 0

Definition at line 119 of file multiplicativeMultigrid.hh.

Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::applyDp().

◆ applyDp()

template<class Entry , class Index , class Smoother , class Prolongation >
virtual field_type Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::applyDp ( domain_type x,
range_type const &  y 
)
inlinevirtual

Computes \( x \leftarrow By \) and returns \( \langle By, y \rangle \).

Implements Kaskade::SymmetricPreconditioner< Smoother::domain_type, Smoother::range_type >.

Definition at line 126 of file multiplicativeMultigrid.hh.

◆ category()

virtual Dune::SolverCategory::Category Kaskade::SymmetricPreconditioner< Smoother::domain_type , Smoother::range_type >::category ( ) const
inlineoverridevirtualinherited

returns the category of the operator

From the Dune doxygen documentation it is unclear what this is supposed to mean. We return a dummy here.

Definition at line 242 of file symmetricOperators.hh.

◆ getCoarsePreconditioner()

template<class Entry , class Index , class Smoother , class Prolongation >
CoarsePreconditioner & Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::getCoarsePreconditioner ( )
inline

Provides access to the coarse grid preconditioner.

Definition at line 185 of file multiplicativeMultigrid.hh.

◆ getSmootherStepSize()

template<class Entry , class Index , class Smoother , class Prolongation >
double Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::getSmootherStepSize ( ) const
inline

Definition at line 163 of file multiplicativeMultigrid.hh.

◆ getSmoothings()

template<class Entry , class Index , class Smoother , class Prolongation >
std::pair< int, int > Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::getSmoothings ( )
inline

Returns a pair of number of pre- and post-smoothing iterations to perform.

Definition at line 193 of file multiplicativeMultigrid.hh.

◆ operator=()

template<class Entry , class Index , class Smoother , class Prolongation >
MultiplicativeMultiGrid & Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::operator= ( MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation > &&  other)
default

◆ post()

virtual void Kaskade::SymmetricPreconditioner< Smoother::domain_type , Smoother::range_type >::post ( Smoother::domain_type &  x)
inlinevirtualinherited

Preconditioner cleanup.

The provided default implementation does nothing.

Definition at line 225 of file symmetricOperators.hh.

◆ pre()

virtual void Kaskade::SymmetricPreconditioner< Smoother::domain_type , Smoother::range_type >::pre ( Smoother::domain_type &  ,
Smoother::range_type &   
)
inlinevirtualinherited

Preconditioner preparation.

The provided default implementation does nothing.

Definition at line 218 of file symmetricOperators.hh.

◆ requiresInitializedInput()

template<class Entry , class Index , class Smoother , class Prolongation >
virtual bool Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::requiresInitializedInput ( ) const
inlinevirtual

Returns true if the target vector x has to be initialized to zero before calling apply or applyDp.

Implements Kaskade::SymmetricPreconditioner< Smoother::domain_type, Smoother::range_type >.

Definition at line 132 of file multiplicativeMultigrid.hh.

◆ setLinesearch()

template<class Entry , class Index , class Smoother , class Prolongation >
void Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::setLinesearch ( bool  ls)
inline

Enable or disable the line search option.

With line search option, at the end of each level in the V-cycle a line search is performed for an optimal scaling of the correction. This may improve the contraction (but need not, for simple Laplace type problems it does not) and guarantee convergence of the multigrid as a fixed point iteration (removing the need for an outer stepsize loop), but incurs one more matrix-vector product (an overhead of up to 30%) and renders the V-cycle a nonlinear scheme, not suited as preconditioner in CG.

The default value on construction is off.

Definition at line 177 of file multiplicativeMultigrid.hh.

◆ setSmootherStepSize()

template<class Entry , class Index , class Smoother , class Prolongation >
void Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::setSmootherStepSize ( double  w)
inline

Define the smoother step size.

The default step length is 1.

Definition at line 157 of file multiplicativeMultigrid.hh.

◆ setSmoothings()

template<class Entry , class Index , class Smoother , class Prolongation >
void Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::setSmoothings ( int  nPre_,
int  nPost_ 
)
inline

Sets the number of pre- and post-smoothing iterations to perform.

Note that the sum of pre- and postsmoothings must be positive.

Parameters
nPrenumber of pre-smoothings (nonnegative)
nPostnumber of post-smoothings (nonnegative)

Definition at line 145 of file multiplicativeMultigrid.hh.


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