
Public Types | |
| enum | Cycle { v_cycle, w_cycle, f_cycle } |
| typedef VECTOR | vector_type |
| typedef const VECTOR | const_vector_type |
Public Member Functions | |
| template<int dim> | |
| Multigrid (const MGDoFHandler< dim > &mg_dof_handler, const MGMatrixBase< VECTOR > &matrix, const MGCoarseGridBase< VECTOR > &coarse, const MGTransferBase< VECTOR > &transfer, const MGSmootherBase< VECTOR > &pre_smooth, const MGSmootherBase< VECTOR > &post_smooth, Cycle cycle=v_cycle) | |
| Multigrid (const unsigned int minlevel, const unsigned int maxlevel, const MGMatrixBase< VECTOR > &matrix, const MGCoarseGridBase< VECTOR > &coarse, const MGTransferBase< VECTOR > &transfer, const MGSmootherBase< VECTOR > &pre_smooth, const MGSmootherBase< VECTOR > &post_smooth, Cycle cycle=v_cycle) | |
| void | reinit (const unsigned int minlevel, const unsigned int maxlevel) |
| void | cycle () |
| void | vcycle () |
| void | vmult (VECTOR &dst, const VECTOR &src) const |
| void | vmult_add (VECTOR &dst, const VECTOR &src) const |
| void | Tvmult (VECTOR &dst, const VECTOR &src) const |
| void | Tvmult_add (VECTOR &dst, const VECTOR &src) const |
| void | set_edge_matrices (const MGMatrixBase< VECTOR > &edge_out, const MGMatrixBase< VECTOR > &edge_in) |
| void | set_edge_flux_matrices (const MGMatrixBase< VECTOR > &edge_down, const MGMatrixBase< VECTOR > &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) |
Public Attributes | |
| MGLevelObject< VECTOR > | defect |
| MGLevelObject< VECTOR > | solution |
Private Member Functions | |
| void | level_v_step (const unsigned int level) |
| void | level_step (const unsigned int level, Cycle cycle) |
Private Attributes | |
| Cycle | cycle_type |
| unsigned int | minlevel |
| unsigned int | maxlevel |
| MGLevelObject< VECTOR > | t |
| MGLevelObject< VECTOR > | defect2 |
| SmartPointer< const MGMatrixBase< VECTOR > , Multigrid< VECTOR > > | matrix |
| SmartPointer< const MGCoarseGridBase< VECTOR > , Multigrid< VECTOR > > | coarse |
| SmartPointer< const MGTransferBase< VECTOR > , Multigrid< VECTOR > > | transfer |
| SmartPointer< const MGSmootherBase< VECTOR > , Multigrid< VECTOR > > | pre_smooth |
| SmartPointer< const MGSmootherBase< VECTOR > , Multigrid< VECTOR > > | post_smooth |
| SmartPointer< const MGMatrixBase< VECTOR > > | edge_out |
| SmartPointer< const MGMatrixBase< VECTOR > > | edge_in |
| SmartPointer< const MGMatrixBase< VECTOR > , Multigrid< VECTOR > > | edge_down |
| SmartPointer< const MGMatrixBase< VECTOR > , Multigrid< VECTOR > > | edge_up |
| unsigned int | debug |
Friends | |
| class | PreconditionMG |
Implementation of the multigrid method.
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 63 of file multigrid.h.
| typedef VECTOR Multigrid< VECTOR >::vector_type |
Definition at line 79 of file multigrid.h.
| typedef const VECTOR Multigrid< VECTOR >::const_vector_type |
Definition at line 80 of file multigrid.h.
| enum Multigrid::Cycle |
List of implemented cycle types.
Definition at line 69 of file multigrid.h.
| Multigrid< VECTOR >::Multigrid | ( | const MGDoFHandler< dim > & | mg_dof_handler, |
| const MGMatrixBase< VECTOR > & | matrix, | ||
| const MGCoarseGridBase< VECTOR > & | coarse, | ||
| const MGTransferBase< VECTOR > & | transfer, | ||
| const MGSmootherBase< VECTOR > & | pre_smooth, | ||
| const MGSmootherBase< VECTOR > & | post_smooth, | ||
| Cycle | cycle = v_cycle |
||
| ) |
Constructor. The MGDoFHandler is used to determine the highest possible level. 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< VECTOR >::Multigrid | ( | const unsigned int | minlevel, |
| const unsigned int | maxlevel, | ||
| const MGMatrixBase< VECTOR > & | matrix, | ||
| const MGCoarseGridBase< VECTOR > & | coarse, | ||
| const MGTransferBase< VECTOR > & | transfer, | ||
| const MGSmootherBase< VECTOR > & | pre_smooth, | ||
| const MGSmootherBase< VECTOR > & | post_smooth, | ||
| Cycle | cycle = v_cycle |
||
| ) |
Experimental constructor for cases in which no MGDoFHandler is available.
| void Multigrid< VECTOR >::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< VECTOR >::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 class MGTools if you want to use these vectors yourself.
The actual work for this function is done in level_v_step().
| void Multigrid< VECTOR >::vmult | ( | VECTOR & | dst, |
| const VECTOR & | src | ||
| ) | const |
Perform a multigrid cycle with a vector which is already a level vector. Use of this function assumes that there is NO local refinement and that both vectors are on the finest level of this Multigrid object.
| void Multigrid< VECTOR >::vmult_add | ( | VECTOR & | dst, |
| const VECTOR & | src | ||
| ) | const |
Perform a multigrid cycle with a vector which is already a level vector. Use of this function assumes that there is NO local refinement and that both vectors are on the finest level of this Multigrid object.
| void Multigrid< VECTOR >::Tvmult | ( | VECTOR & | dst, |
| const VECTOR & | src | ||
| ) | const |
| void Multigrid< VECTOR >::Tvmult_add | ( | VECTOR & | dst, |
| const VECTOR & | src | ||
| ) | const |
| void Multigrid< VECTOR >::set_edge_matrices | ( | const MGMatrixBase< VECTOR > & | edge_out, |
| const MGMatrixBase< VECTOR > & | 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< VECTOR >::set_edge_flux_matrices | ( | const MGMatrixBase< VECTOR > & | edge_down, |
| const MGMatrixBase< VECTOR > & | 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. Return the finest level for multigrid.
Return the coarsest level for multigrid.
Set the highest level for which the multilevel method is performed. By default, this is the finest level of the Triangulation; therefore, this function will only accept arguments smaller than the current maxlevel and not smaller than the current minlevel.
| void Multigrid< VECTOR >::set_minlevel | ( | const unsigned int | level, |
| bool | relative = false |
||
| ) |
Chance cycle_type used in cycle().
Set the debug level. Higher values will create more debugging output during the multigrid cycles.
| void Multigrid< VECTOR >::level_v_step | ( | const unsigned int | level ) | [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.
| void Multigrid< VECTOR >::level_step | ( | const unsigned int | level, |
| Cycle | cycle | ||
| ) | [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.
friend class PreconditionMG [friend] |
Definition at line 490 of file multigrid.h.
Cycle Multigrid< VECTOR >::cycle_type [private] |
Cycle type performed by the method cycle().
Definition at line 381 of file multigrid.h.
Level for coarse grid solution.
Definition at line 386 of file multigrid.h.
Highest level of cells.
Definition at line 391 of file multigrid.h.
| MGLevelObject<VECTOR> Multigrid< VECTOR >::defect |
Input vector for the cycle. Contains the defect of the outer method projected to the multilevel vectors.
Definition at line 400 of file multigrid.h.
| MGLevelObject<VECTOR> Multigrid< VECTOR >::solution |
The solution update after the multigrid step.
Definition at line 406 of file multigrid.h.
MGLevelObject<VECTOR> Multigrid< VECTOR >::t [private] |
Auxiliary vector.
Definition at line 412 of file multigrid.h.
MGLevelObject<VECTOR> Multigrid< VECTOR >::defect2 [private] |
Auxiliary vector for W- and F-cycles. Left uninitialized in V-cycle.
Definition at line 419 of file multigrid.h.
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::matrix [private] |
The matrix for each level.
Definition at line 425 of file multigrid.h.
SmartPointer<const MGCoarseGridBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::coarse [private] |
The matrix for each level.
Definition at line 430 of file multigrid.h.
SmartPointer<const MGTransferBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::transfer [private] |
Object for grid tranfer.
Definition at line 435 of file multigrid.h.
SmartPointer<const MGSmootherBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::pre_smooth [private] |
The pre-smoothing object.
Definition at line 440 of file multigrid.h.
SmartPointer<const MGSmootherBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::post_smooth [private] |
The post-smoothing object.
Definition at line 445 of file multigrid.h.
SmartPointer<const MGMatrixBase<VECTOR> > Multigrid< VECTOR >::edge_out [private] |
Edge matrix from the interior of the refined part to the refinement edge.
vmult is used for these matrices. Definition at line 455 of file multigrid.h.
SmartPointer<const MGMatrixBase<VECTOR> > Multigrid< VECTOR >::edge_in [private] |
Transpose edge matrix from the refinement edge to the interior of the refined part.
Tvmult is used for these matrices. Definition at line 465 of file multigrid.h.
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::edge_down [private] |
Edge matrix from fine to coarse.
vmult is used for these matrices. Definition at line 473 of file multigrid.h.
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::edge_up [private] |
Transpose edge matrix from coarse to fine.
Tvmult is used for these matrices. Definition at line 481 of file multigrid.h.
Level for debug output. Defaults to zero and can be set by set_debug().
Definition at line 488 of file multigrid.h.
documentation generated on Fri Feb 3 2012 06:04:09 by
doxygen
1.7.2