Reference documentation for deal.II version 8.1.0
Public Types | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | Friends | List of all members
Multigrid< VECTOR > Class Template Reference

#include <multigrid.h>

Inheritance diagram for Multigrid< VECTOR >:
[legend]

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 DoFHandler< 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 DEAL_II_DEPRECATED
 
void vmult_add (VECTOR &dst, const VECTOR &src) const DEAL_II_DEPRECATED
 
void Tvmult (VECTOR &dst, const VECTOR &src) const DEAL_II_DEPRECATED
 
void Tvmult_add (VECTOR &dst, const VECTOR &src) const DEAL_II_DEPRECATED
 
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)
 
template<class VECTOR >
DEAL_II_NAMESPACE_OPEN 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, typename Multigrid< VECTOR >::Cycle cycle)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
 DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.")
 
 DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

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

template<int dim, class VECTOR2 , class TRANSFER >
class PreconditionMG
 

Detailed Description

template<class VECTOR>
class Multigrid< VECTOR >

Implementation of the multigrid method.

Warning
multigrid on locally refined meshes only works with discontinuous finite elements right now. It is not clear, whether the paradigm of local smoothing we use is applicable to continuous elements with hanging nodes; in fact, most people you meet on conferences seem to deny this.

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.

Note
The interface of this class is still very clumsy. In particular, you will have to set up quite a few auxiliary objects before you can use it. Unfortunately, it seems that this can be avoided only be restricting the flexibility of this class in an unacceptable way.
Author
Guido Kanschat, 1999 - 2005

Definition at line 68 of file multigrid.h.

Constructor & Destructor Documentation

template<class VECTOR>
template<int dim>
Multigrid< VECTOR >::Multigrid ( const DoFHandler< 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 DoFHandler 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.

template<class VECTOR>
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 DoFHandler is available.

Warning
Not intended for general use.

Member Function Documentation

template<class VECTOR >
void Multigrid< VECTOR >::reinit ( const unsigned int  minlevel,
const unsigned int  maxlevel 
)

Reinit this class according to minlevel and maxlevel.

Definition at line 60 of file multigrid.templates.h.

template<class VECTOR >
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.

Definition at line 379 of file multigrid.templates.h.

template<class VECTOR >
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().

Definition at line 406 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::vmult ( VECTOR &  dst,
const VECTOR &  src 
) const
Deprecated:
This function is purely experimental and will probably never be implemented in a way that it can be released.

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.

Definition at line 425 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::vmult_add ( VECTOR &  dst,
const VECTOR &  src 
) const
Deprecated:
This function is purely experimental and will probably never be implemented in a way that it can be released.

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.

Definition at line 444 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::Tvmult ( VECTOR &  dst,
const VECTOR &  src 
) const
Deprecated:
Even worse than vmult(), this function is not even implemented, but just declared such that certain objects relying on it can be constructed.

Definition at line 455 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::Tvmult_add ( VECTOR &  dst,
const VECTOR &  src 
) const
Deprecated:
Even worse than vmult(), this function is not even implemented, but just declared such that certain objects relying on it can be constructed.

Definition at line 463 of file multigrid.templates.h.

template<class VECTOR >
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.

Note
While 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.

Definition at line 109 of file multigrid.templates.h.

template<class VECTOR >
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.

Note
While 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.

Definition at line 119 of file multigrid.templates.h.

template<class VECTOR>
unsigned int Multigrid< VECTOR >::get_maxlevel ( ) const

Return the finest level for multigrid.

template<class VECTOR>
unsigned int Multigrid< VECTOR >::get_minlevel ( ) const

Return the coarsest level for multigrid.

template<class VECTOR >
void Multigrid< VECTOR >::set_maxlevel ( const unsigned int  l)

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.

Definition at line 71 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::set_minlevel ( const unsigned int  level,
bool  relative = false 
)

Set the coarse level for which the multilevel method is performed. By default, this is zero. Accepted are non-negative values not larger than than the current maxlevel.

If relative ist true, then this function determins the number of levels used, that is, it sets minlevel to maxlevel-level.

Note
The mesh on the coarsest level must cover the whole domain. There may not be hanging nodes on minlevel.
If minlevel is set to a nonzero value, do not forget to adjust your coarse grid solver!

Definition at line 81 of file multigrid.templates.h.

template<class VECTOR>
void Multigrid< VECTOR >::set_cycle ( Cycle  )

Chance cycle_type used in cycle().

Definition at line 93 of file multigrid.templates.h.

template<class VECTOR >
void Multigrid< VECTOR >::set_debug ( const unsigned int  d)

Set the debug level. Higher values will create more debugging output during the multigrid cycles.

Definition at line 101 of file multigrid.templates.h.

template<class VECTOR >
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.

Definition at line 129 of file multigrid.templates.h.

template<class VECTOR >
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.

Definition at line 249 of file multigrid.templates.h.

Member Data Documentation

template<class VECTOR>
Cycle Multigrid< VECTOR >::cycle_type
private

Cycle type performed by the method cycle().

Definition at line 386 of file multigrid.h.

template<class VECTOR>
unsigned int Multigrid< VECTOR >::minlevel
private

Level for coarse grid solution.

Definition at line 391 of file multigrid.h.

template<class VECTOR>
unsigned int Multigrid< VECTOR >::maxlevel
private

Highest level of cells.

Definition at line 396 of file multigrid.h.

template<class VECTOR>
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 405 of file multigrid.h.

template<class VECTOR>
MGLevelObject<VECTOR> Multigrid< VECTOR >::solution

The solution update after the multigrid step.

Definition at line 411 of file multigrid.h.

template<class VECTOR>
MGLevelObject<VECTOR> Multigrid< VECTOR >::t
private

Auxiliary vector.

Definition at line 417 of file multigrid.h.

template<class VECTOR>
MGLevelObject<VECTOR> Multigrid< VECTOR >::defect2
private

Auxiliary vector for W- and F-cycles. Left uninitialized in V-cycle.

Definition at line 424 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::matrix
private

The matrix for each level.

Definition at line 430 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGCoarseGridBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::coarse
private

The matrix for each level.

Definition at line 435 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGTransferBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::transfer
private

Object for grid tranfer.

Definition at line 440 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGSmootherBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::pre_smooth
private

The pre-smoothing object.

Definition at line 445 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGSmootherBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::post_smooth
private

The post-smoothing object.

Definition at line 450 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGMatrixBase<VECTOR> > Multigrid< VECTOR >::edge_out
private

Edge matrix from the interior of the refined part to the refinement edge.

Note
Only vmult is used for these matrices.

Definition at line 460 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGMatrixBase<VECTOR> > Multigrid< VECTOR >::edge_in
private

Transpose edge matrix from the refinement edge to the interior of the refined part.

Note
Only Tvmult is used for these matrices.

Definition at line 470 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::edge_down
private

Edge matrix from fine to coarse.

Note
Only vmult is used for these matrices.

Definition at line 478 of file multigrid.h.

template<class VECTOR>
SmartPointer<const MGMatrixBase<VECTOR>,Multigrid<VECTOR> > Multigrid< VECTOR >::edge_up
private

Transpose edge matrix from coarse to fine.

Note
Only Tvmult is used for these matrices.

Definition at line 486 of file multigrid.h.

template<class VECTOR>
unsigned int Multigrid< VECTOR >::debug
private

Level for debug output. Defaults to zero and can be set by set_debug().

Definition at line 493 of file multigrid.h.


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