Reference documentation for deal.II version 9.1.1
|
#include <deal.II/multigrid/multigrid.h>
Public Types | |
enum | Cycle { v_cycle, w_cycle, f_cycle } |
Public Member Functions | |
template<int dim> | |
Multigrid (const DoFHandler< dim > &mg_dof_handler, const MGMatrixBase< VectorType > &matrix, const MGCoarseGridBase< VectorType > &coarse, const MGTransferBase< VectorType > &transfer, const MGSmootherBase< VectorType > &pre_smooth, const MGSmootherBase< VectorType > &post_smooth, const unsigned int minlevel=0, const unsigned int maxlevel=numbers::invalid_unsigned_int, Cycle cycle=v_cycle) | |
Multigrid (const MGMatrixBase< VectorType > &matrix, const MGCoarseGridBase< VectorType > &coarse, const MGTransferBase< VectorType > &transfer, const MGSmootherBase< VectorType > &pre_smooth, const MGSmootherBase< VectorType > &post_smooth, const unsigned int minlevel=0, const unsigned int maxlevel=numbers::invalid_unsigned_int, Cycle cycle=v_cycle) | |
void | reinit (const unsigned int minlevel, const unsigned int maxlevel) |
void | cycle () |
void | vcycle () |
void | set_edge_matrices (const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in) |
void | set_edge_flux_matrices (const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up) |
unsigned int | get_maxlevel () const |
unsigned int | get_minlevel () const |
void | set_maxlevel (const unsigned int) |
void | set_minlevel (const unsigned int level, bool relative=false) |
void | set_cycle (Cycle) |
void | set_debug (const unsigned int) |
boost::signals2::connection | connect_coarse_solve (const std::function< void(const bool, const unsigned int)> &slot) |
boost::signals2::connection | connect_restriction (const std::function< void(const bool, const unsigned int)> &slot) |
boost::signals2::connection | connect_prolongation (const std::function< void(const bool, const unsigned int)> &slot) |
boost::signals2::connection | connect_pre_smoother_step (const std::function< void(const bool, const unsigned int)> &slot) |
boost::signals2::connection | connect_post_smoother_step (const std::function< void(const bool, const unsigned int)> &slot) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
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) |
Public Attributes | |
MGLevelObject< VectorType > | defect |
MGLevelObject< VectorType > | solution |
Private Member Functions | |
void | level_v_step (const unsigned int level) |
void | level_step (const unsigned int level, Cycle cycle) |
Private Attributes | |
mg::Signals | signals |
Cycle | cycle_type |
unsigned int | minlevel |
unsigned int | maxlevel |
MGLevelObject< VectorType > | t |
MGLevelObject< VectorType > | defect2 |
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > | matrix |
SmartPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > | coarse |
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > | transfer |
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > | pre_smooth |
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > | post_smooth |
SmartPointer< const MGMatrixBase< VectorType > > | edge_out |
SmartPointer< const MGMatrixBase< VectorType > > | edge_in |
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > | edge_down |
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > | edge_up |
unsigned int | debug |
Friends | |
template<int dim, class OtherVectorType , class TRANSFER > | |
class | PreconditionMG |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Implementation of the multigrid method. The implementation supports both continuous and discontinuous elements and follows the procedure described in the multigrid paper by Janssen and Kanschat.
The function which starts a multigrid cycle on the finest level is cycle(). Depending on the cycle type chosen with the constructor (see enum Cycle), this function triggers one of the cycles level_v_step() or level_step(), where the latter one can do different types of cycles.
Using this class, it is expected that the right hand side has been converted from a vector living on the locally finest level to a multilevel vector. This is a nontrivial operation, usually initiated automatically by the class PreconditionMG and performed by the classes derived from MGTransferBase.
Definition at line 137 of file multigrid.h.
enum Multigrid::Cycle |
List of implemented cycle types.
Enumerator | |
---|---|
v_cycle | The V-cycle. |
w_cycle | The W-cycle. |
f_cycle | The F-cycle. |
Definition at line 143 of file multigrid.h.
Multigrid< VectorType >::Multigrid | ( | const DoFHandler< dim > & | mg_dof_handler, |
const MGMatrixBase< VectorType > & | matrix, | ||
const MGCoarseGridBase< VectorType > & | coarse, | ||
const MGTransferBase< VectorType > & | transfer, | ||
const MGSmootherBase< VectorType > & | pre_smooth, | ||
const MGSmootherBase< VectorType > & | post_smooth, | ||
const unsigned int | minlevel = 0 , |
||
const unsigned int | maxlevel = numbers::invalid_unsigned_int , |
||
Cycle | cycle = v_cycle |
||
) |
Constructor. The DoFHandler is used to check whether the provided minlevel and maxlevel are in the range of valid levels. If maxlevel is set to the default value, the highest valid level is used. transfer
is an object performing prolongation and restriction.
This function already initializes the vectors which will be used later in the course of the computations. You should therefore create objects of this type as late as possible.
Multigrid< VectorType >::Multigrid | ( | const MGMatrixBase< VectorType > & | matrix, |
const MGCoarseGridBase< VectorType > & | coarse, | ||
const MGTransferBase< VectorType > & | transfer, | ||
const MGSmootherBase< VectorType > & | pre_smooth, | ||
const MGSmootherBase< VectorType > & | post_smooth, | ||
const unsigned int | minlevel = 0 , |
||
const unsigned int | maxlevel = numbers::invalid_unsigned_int , |
||
Cycle | cycle = v_cycle |
||
) |
Constructor. transfer
is an object performing prolongation and restriction. For levels in [minlevel, maxlevel] matrix has to contain valid matrices. By default the maxlevel is set to the maximal valid level.
This function already initializes the vectors which will be used later in the course of the computations. You should therefore create objects of this type as late as possible.
void Multigrid< VectorType >::reinit | ( | const unsigned int | minlevel, |
const unsigned int | maxlevel | ||
) |
void Multigrid< VectorType >::cycle | ( | ) |
Execute one multigrid cycle. The type of cycle is selected by the constructor argument cycle. See the enum Cycle for available types.
void Multigrid< VectorType >::vcycle | ( | ) |
Execute one step of the V-cycle algorithm. This function assumes, that the multilevel vector defect is filled with the residual of an outer defect correction scheme. This is usually taken care of by PreconditionMG). After vcycle(), the result is in the multilevel vector solution. See copy_*_mg
in the MGTools namespace if you want to use these vectors yourself.
The actual work for this function is done in level_v_step().
void Multigrid< VectorType >::set_edge_matrices | ( | const MGMatrixBase< VectorType > & | edge_out, |
const MGMatrixBase< VectorType > & | edge_in | ||
) |
Set additional matrices to correct residual computation at refinement edges. Since we only smoothen in the interior of the refined part of the mesh, the coupling across the refinement edge is missing. This coupling is provided by these two matrices.
edge_out.vmult
is used, for the second argument, we use edge_in.Tvmult
. Thus, edge_in
should be assembled in transposed form. This saves a second sparsity pattern for edge_in
. In particular, for symmetric operators, both arguments can refer to the same matrix, saving assembling of one of them. void Multigrid< VectorType >::set_edge_flux_matrices | ( | const MGMatrixBase< VectorType > & | edge_down, |
const MGMatrixBase< VectorType > & | edge_up | ||
) |
Set additional matrices to correct residual computation at refinement edges. These matrices originate from discontinuous Galerkin methods (see FE_DGQ etc.), where they correspond to the edge fluxes at the refinement edge between two levels.
edge_down.vmult
is used, for the second argument, we use edge_up.Tvmult
. Thus, edge_up
should be assembled in transposed form. This saves a second sparsity pattern for edge_up
. In particular, for symmetric operators, both arguments can refer to the same matrix, saving assembling of one of them. unsigned int Multigrid< VectorType >::get_maxlevel | ( | ) | const |
Return the finest level for multigrid.
unsigned int Multigrid< VectorType >::get_minlevel | ( | ) | const |
Return the coarsest level for multigrid.
void Multigrid< VectorType >::set_maxlevel | ( | const unsigned | int | ) |
Set the highest level for which the multilevel method is performed. By default, this is the finest level of the Triangulation. Accepted are values not smaller than the current minlevel.
void Multigrid< VectorType >::set_minlevel | ( | const unsigned int | level, |
bool | relative = false |
||
) |
Chance cycle_type used in cycle().
void Multigrid< VectorType >::set_debug | ( | const unsigned | int | ) |
Set the debug level. Higher values will create more debugging output during the multigrid cycles.
boost::signals2::connection Multigrid< VectorType >::connect_coarse_solve | ( | const std::function< void(const bool, const unsigned int)> & | slot | ) |
Connect a function to mg::Signals::coarse_solve.
boost::signals2::connection Multigrid< VectorType >::connect_restriction | ( | const std::function< void(const bool, const unsigned int)> & | slot | ) |
Connect a function to mg::Signals::restriction.
boost::signals2::connection Multigrid< VectorType >::connect_prolongation | ( | const std::function< void(const bool, const unsigned int)> & | slot | ) |
Connect a function to mg::Signals::prolongation.
boost::signals2::connection Multigrid< VectorType >::connect_pre_smoother_step | ( | const std::function< void(const bool, const unsigned int)> & | slot | ) |
Connect a function to mg::Signals::pre_smoother_step.
boost::signals2::connection Multigrid< VectorType >::connect_post_smoother_step | ( | const std::function< void(const bool, const unsigned int)> & | slot | ) |
Connect a function to mg::Signals::post_smoother_step.
|
private |
The V-cycle multigrid method. level
is the level the function starts on. It will usually be called for the highest level from outside, but will then call itself recursively for level-1
, unless we are on minlevel where the coarse grid solver solves the problem exactly.
|
private |
The actual W-cycle or F-cycle multigrid method. level
is the level the function starts on. It will usually be called for the highest level from outside, but will then call itself recursively for level-1
, unless we are on minlevel where the coarse grid solver solves the problem exactly.
|
private |
Signals for the various actions that the Multigrid algorithm uses.
Definition at line 353 of file multigrid.h.
Cycle type performed by the method cycle().
Definition at line 377 of file multigrid.h.
|
private |
Level for coarse grid solution.
Definition at line 382 of file multigrid.h.
|
private |
Highest level of cells.
Definition at line 387 of file multigrid.h.
MGLevelObject<VectorType> Multigrid< VectorType >::defect |
Input vector for the cycle. Contains the defect of the outer method projected to the multilevel vectors.
Definition at line 394 of file multigrid.h.
MGLevelObject<VectorType> Multigrid< VectorType >::solution |
The solution update after the multigrid step.
Definition at line 399 of file multigrid.h.
|
private |
Auxiliary vector.
Definition at line 405 of file multigrid.h.
|
private |
Auxiliary vector for W- and F-cycles. Left uninitialized in V-cycle.
Definition at line 410 of file multigrid.h.
|
private |
The matrix for each level.
Definition at line 416 of file multigrid.h.
|
private |
The matrix for each level.
Definition at line 422 of file multigrid.h.
|
private |
Object for grid transfer.
Definition at line 428 of file multigrid.h.
|
private |
The pre-smoothing object.
Definition at line 434 of file multigrid.h.
|
private |
The post-smoothing object.
Definition at line 440 of file multigrid.h.
|
private |
Edge matrix from the interior of the refined part to the refinement edge.
vmult
is used for these matrices. Definition at line 447 of file multigrid.h.
|
private |
Transpose edge matrix from the refinement edge to the interior of the refined part.
Tvmult
is used for these matrices. Definition at line 455 of file multigrid.h.
|
private |
Edge matrix from fine to coarse.
vmult
is used for these matrices. Definition at line 462 of file multigrid.h.
|
private |
Transpose edge matrix from coarse to fine.
Tvmult
is used for these matrices. Definition at line 469 of file multigrid.h.
|
private |
Level for debug output. Defaults to zero and can be set by set_debug().
Definition at line 474 of file multigrid.h.