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

Class for multigrid stacks. More...

#include <prolongation.hh>

Detailed Description

template<class Prolongation, class Entry, class Index>
class Kaskade::MultiGridStack< Prolongation, Entry, Index >

Class for multigrid stacks.

This provides the storage of and access to prolongations and projected Galerkin matrices, but leaves the construction of these to derived classes.

Definition at line 747 of file prolongation.hh.

Public Member Functions

 MultiGridStack (std::vector< Prolongation > &&ps, std::vector< NumaBCRSMatrix< Entry, Index > > &&as)
 Constructor. More...
 
 MultiGridStack (std::vector< Prolongation > &&ps, NumaBCRSMatrix< Entry, Index > &&A, bool onlyLowerTriangle)
 Constructor. More...
 
 MultiGridStack (MultiGridStack &&other)=default
 
int levels () const
 The number of grid levels. More...
 
Prolongation const & p (int level) const
 Returns the prolongation from given level to next higher one. More...
 
NumaBCRSMatrix< Entry, Index > const & a (int level) const
 Returns the projected Galerkin matrix on the given level. More...
 
NumaBCRSMatrix< Entry, Index > & coarseGridMatrix ()
 Returns the projected Galerkin matrix on the coarsest level. More...
 
void report (std::ostream &out) const
 

Related Functions

(Note that these are not member functions.)

template<class Prolongation , class Entry , class Index >
MultiGridStack< Prolongation, Entry, Index > makeMultiGridStack (std::vector< Prolongation > &&ps, NumaBCRSMatrix< Entry, Index > &&A, bool onlyLowerTriangle)
 Convenience routine for creating multigrid stacks. More...
 

Constructor & Destructor Documentation

◆ MultiGridStack() [1/3]

template<class Prolongation , class Entry , class Index >
Kaskade::MultiGridStack< Prolongation, Entry, Index >::MultiGridStack ( std::vector< Prolongation > &&  ps,
std::vector< NumaBCRSMatrix< Entry, Index > > &&  as 
)
inline

Constructor.

This takes both a stack of prolongations and a matching stack of projected.

Parameters
psa stack of prolongations
asa stack of stiffness matrices. as.size()==ps.size()+1 has to hold.

Usually, it holds that \( A_i = P_i^T A_{i+1} P_i \), i.e. as[0] contains the coarse grid matrix.

Definition at line 760 of file prolongation.hh.

◆ MultiGridStack() [2/3]

template<class Prolongation , class Entry , class Index >
Kaskade::MultiGridStack< Prolongation, Entry, Index >::MultiGridStack ( std::vector< Prolongation > &&  ps,
NumaBCRSMatrix< Entry, Index > &&  A,
bool  onlyLowerTriangle 
)
inline

Constructor.

This takes a stack of prolongations and creates the stack of projected Galerkin matrices from the given fine grid matrix. Most useful for geometric multigrid, where the prolongations are defined solely in terms of the grid.

See also
makeMultiGridStack

Definition at line 774 of file prolongation.hh.

◆ MultiGridStack() [3/3]

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

Member Function Documentation

◆ a()

template<class Prolongation , class Entry , class Index >
NumaBCRSMatrix< Entry, Index > const & Kaskade::MultiGridStack< Prolongation, Entry, Index >::a ( int  level) const
inline

Returns the projected Galerkin matrix on the given level.

Parameters
levelprecondition 0 <= level < levels()

For level==0, the returned matrix can be in an undefined state if it has been modified previously via the non-const coarseGridMatrix() reference.

Definition at line 819 of file prolongation.hh.

Referenced by Kaskade::MultiGridStack< Prolongation, Entry, Index >::report().

◆ coarseGridMatrix()

template<class Prolongation , class Entry , class Index >
NumaBCRSMatrix< Entry, Index > & Kaskade::MultiGridStack< Prolongation, Entry, Index >::coarseGridMatrix ( )
inline

Returns the projected Galerkin matrix on the coarsest level.

This is explicitly a mutable reference, such that the coarse grid matrix can be moved from. This is useful, as in the multigrid, the coarsest level matrix is not referenced (the coarse grid preconditioner has its own copy, maybe obtained by moving from here...):

auto coarseSolver = makeDirectPreconditioner(std::move(mgStack.coarseGridMatrix()));
auto makeDirectPreconditioner(Matrix &&A, DirectType directType=DirectType::MUMPS)
Creates a direct solver for the given matrix.

Definition at line 835 of file prolongation.hh.

◆ levels()

template<class Prolongation , class Entry , class Index >
int Kaskade::MultiGridStack< Prolongation, Entry, Index >::levels ( ) const
inline

◆ p()

template<class Prolongation , class Entry , class Index >
Prolongation const & Kaskade::MultiGridStack< Prolongation, Entry, Index >::p ( int  level) const
inline

Returns the prolongation from given level to next higher one.

Parameters
levelprecondition 0 <= level < levels()-1

Definition at line 805 of file prolongation.hh.

Referenced by Kaskade::MultiGridStack< Prolongation, Entry, Index >::report().

◆ report()

template<class Prolongation , class Entry , class Index >
void Kaskade::MultiGridStack< Prolongation, Entry, Index >::report ( std::ostream &  out) const
inline

Definition at line 840 of file prolongation.hh.

Referenced by Kaskade::operator<<().


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