Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Attributes | Private Attributes | List of all members
MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType > Class Template Reference

#include <deal.II/multigrid/mg_smoother.h>

Inheritance diagram for MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >:
[legend]

Public Member Functions

 MGSmootherPrecondition (const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
 
template<typename MatrixType2 >
void initialize (const MGLevelObject< MatrixType2 > &matrices, const typename PreconditionerType::AdditionalData &additional_data=typename PreconditionerType::AdditionalData())
 
template<typename MatrixType2 >
void initialize_matrices (const MGLevelObject< MatrixType2 > &matrices)
 
template<typename MatrixType2 , class DATA >
void initialize (const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DATA > &additional_data)
 
template<typename MatrixType2 , class DATA >
void initialize (const MGLevelObject< MatrixType2 > &matrices, const DATA &additional_data, const unsigned int block_row, const unsigned int block_col)
 
template<typename MatrixType2 , class DATA >
void initialize (const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DATA > &additional_data, const unsigned int block_row, const unsigned int block_col)
 
void clear () override
 
virtual void smooth (const unsigned int level, VectorType &u, const VectorType &rhs) const override
 
virtual void apply (const unsigned int level, VectorType &u, const VectorType &rhs) const override
 
std::size_t memory_consumption () const
 
void set_steps (const unsigned int)
 
void set_variable (const bool)
 
void set_symmetric (const bool)
 
void set_transpose (const bool)
 
void set_debug (const unsigned int level)
 

Public Attributes

MGLevelObject< PreconditionerType > smoothers
 

Protected Attributes

GrowingVectorMemory< VectorType > vector_memory
 
unsigned int steps
 
bool variable
 
bool symmetric
 
bool transpose
 
unsigned int debug
 

Private Attributes

MGLevelObject< LinearOperator< VectorType > > matrices
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

Detailed Description

template<typename MatrixType, typename PreconditionerType, typename VectorType>
class MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >

Smoother using preconditioner classes.

This class performs smoothing on each level. The operation can be controlled by several parameters. First, the relaxation parameter omega is used in the underlying relaxation method. steps is the number of relaxation steps on the finest level (on all levels if variable is off). If variable is true, the number of smoothing steps is doubled on each coarser level. This results in a method having the complexity of the W-cycle, but saving grid transfers. This is the method proposed by Bramble at al.

The option symmetric switches on alternating between the smoother and its transpose in each step as proposed by Bramble.

transpose uses the transposed smoothing operation using Tvmult instead of the regular vmult of the relaxation scheme.

If you are using block matrices, the second initialize function offers the possibility to extract a single block for smoothing. In this case, the multigrid method must be used only with the vector associated to that single block.

The library contains instantiation for SparseMatrix<.> and Vector<.>, where the template arguments are all combinations of float and double. Additional instantiations may be created by including the file mg_smoother.templates.h.

Definition at line 444 of file mg_smoother.h.

Constructor & Destructor Documentation

◆ MGSmootherPrecondition()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::MGSmootherPrecondition ( const unsigned int  steps = 1,
const bool  variable = false,
const bool  symmetric = false,
const bool  transpose = false 
)

Constructor. Sets smoothing parameters.

Member Function Documentation

◆ initialize() [1/4]

template<typename MatrixType , typename PreconditionerType , typename VectorType >
template<typename MatrixType2 >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::initialize ( const MGLevelObject< MatrixType2 > &  matrices,
const typename PreconditionerType::AdditionalData &  additional_data = typename PreconditionerType::AdditionalData() 
)

Initialize for matrices. This function stores pointers to the level matrices and initializes the smoothing operator with the same smoother for each level.

additional_data is an object of type PreconditionerType::AdditionalData and is handed to the initialization function of the relaxation method.

◆ initialize_matrices()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
template<typename MatrixType2 >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::initialize_matrices ( const MGLevelObject< MatrixType2 > &  matrices)

In contrast to the function above, only initialize the matrices. The smoothers need to be set up manually by the user as a subsequent step in the code. For this purpose, the public field smoothers can be directly modified. This is useful if one wants full flexibility in the choice of smoothers, e.g., use different smoothers on the levels.

◆ initialize() [2/4]

template<typename MatrixType , typename PreconditionerType , typename VectorType >
template<typename MatrixType2 , class DATA >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::initialize ( const MGLevelObject< MatrixType2 > &  matrices,
const MGLevelObject< DATA > &  additional_data 
)

Initialize for matrices. This function stores pointers to the level matrices and initializes the smoothing operator with the according smoother for each level.

additional_data is an object of type PreconditionerType::AdditionalData and is handed to the initialization function of the relaxation method.

◆ initialize() [3/4]

template<typename MatrixType , typename PreconditionerType , typename VectorType >
template<typename MatrixType2 , class DATA >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::initialize ( const MGLevelObject< MatrixType2 > &  matrices,
const DATA &  additional_data,
const unsigned int  block_row,
const unsigned int  block_col 
)

Initialize for single blocks of matrices. Of this block matrix, the block indicated by block_row and block_col is selected on each level. This function stores pointers to the level matrices and initializes the smoothing operator with the same smoother for each level.

additional_data is an object of type PreconditionerType::AdditionalData and is handed to the initialization function of the relaxation method.

◆ initialize() [4/4]

template<typename MatrixType , typename PreconditionerType , typename VectorType >
template<typename MatrixType2 , class DATA >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::initialize ( const MGLevelObject< MatrixType2 > &  matrices,
const MGLevelObject< DATA > &  additional_data,
const unsigned int  block_row,
const unsigned int  block_col 
)

Initialize for single blocks of matrices. Of this block matrix, the block indicated by block_row and block_col is selected on each level. This function stores pointers to the level matrices and initializes the smoothing operator with the according smoother for each level.

additional_data is an object of type PreconditionerType::AdditionalData and is handed to the initialization function of the relaxation method.

◆ clear()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::clear ( )
overridevirtual

Empty all vectors.

Implements MGSmootherBase< VectorType >.

◆ smooth()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
virtual void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::smooth ( const unsigned int  level,
VectorType &  u,
const VectorType &  rhs 
) const
overridevirtual

The actual smoothing method.

Implements MGSmootherBase< VectorType >.

◆ apply()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
virtual void MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::apply ( const unsigned int  level,
VectorType &  u,
const VectorType &  rhs 
) const
overridevirtual

The apply variant of smoothing, setting the vector u to zero before calling the smooth function. This function is equivalent to the following code

u = 0;
smooth(level, u, rhs);
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const override
unsigned int level
Definition: grid_out.cc:4606

In the multigrid preconditioner interfaces, the apply() method is used for the pre-smoothing operation because the previous content in the solution vector needs to be overwritten for a new incoming residual. On the other hand, all subsequent operations need to smooth the content already present in the vector u given the right hand side, which is done by smooth().

Reimplemented from MGSmootherBase< VectorType >.

◆ memory_consumption()

template<typename MatrixType , typename PreconditionerType , typename VectorType >
std::size_t MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::memory_consumption ( ) const

Memory used by this object.

◆ set_steps()

template<typename VectorType >
void MGSmoother< VectorType >::set_steps ( const unsigned int  )
inherited

Modify the number of smoothing steps on finest level.

◆ set_variable()

template<typename VectorType >
void MGSmoother< VectorType >::set_variable ( const bool  )
inherited

Switch on/off variable smoothing.

◆ set_symmetric()

template<typename VectorType >
void MGSmoother< VectorType >::set_symmetric ( const bool  )
inherited

Switch on/off symmetric smoothing.

◆ set_transpose()

template<typename VectorType >
void MGSmoother< VectorType >::set_transpose ( const bool  )
inherited

Switch on/off transposed smoothing. The effect is overridden by set_symmetric().

◆ set_debug()

template<typename VectorType >
void MGSmoother< VectorType >::set_debug ( const unsigned int  level)
inherited

Set debug to a nonzero value to get debug information logged to deallog. Increase to get more information

Member Data Documentation

◆ smoothers

template<typename MatrixType , typename PreconditionerType , typename VectorType >
MGLevelObject<PreconditionerType> MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::smoothers

Object containing relaxation methods.

Definition at line 566 of file mg_smoother.h.

◆ matrices

template<typename MatrixType , typename PreconditionerType , typename VectorType >
MGLevelObject<LinearOperator<VectorType> > MGSmootherPrecondition< MatrixType, PreconditionerType, VectorType >::matrices
private

Pointer to the matrices.

Definition at line 579 of file mg_smoother.h.

◆ vector_memory

template<typename VectorType >
GrowingVectorMemory<VectorType> MGSmoother< VectorType >::vector_memory
mutableprotectedinherited

A memory object to be used for temporary vectors.

The object is marked as mutable since we will need to use it to allocate temporary vectors also in functions that are const.

Definition at line 98 of file mg_smoother.h.

◆ steps

template<typename VectorType >
unsigned int MGSmoother< VectorType >::steps
protectedinherited

Number of smoothing steps on the finest level. If no variable smoothing is chosen, this is the number of steps on all levels.

Definition at line 104 of file mg_smoother.h.

◆ variable

template<typename VectorType >
bool MGSmoother< VectorType >::variable
protectedinherited

Variable smoothing: double the number of smoothing steps whenever going to the next coarser level

Definition at line 110 of file mg_smoother.h.

◆ symmetric

template<typename VectorType >
bool MGSmoother< VectorType >::symmetric
protectedinherited

Symmetric smoothing: in the smoothing iteration, alternate between the relaxation method and its transpose.

Definition at line 116 of file mg_smoother.h.

◆ transpose

template<typename VectorType >
bool MGSmoother< VectorType >::transpose
protectedinherited

Use the transpose of the relaxation method instead of the method itself. This has no effect if symmetric smoothing is chosen.

Definition at line 122 of file mg_smoother.h.

◆ debug

template<typename VectorType >
unsigned int MGSmoother< VectorType >::debug
protectedinherited

Output debugging information to deallog if this is nonzero.

Definition at line 127 of file mg_smoother.h.


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