Reference documentation for deal.II version 9.5.0
\(\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
Classes | Public Member Functions | Static Public Member Functions | Private Types | Private Attributes | List of all members
SolutionTransfer< dim, VectorType, spacedim > Class Template Reference

#include <deal.II/numerics/solution_transfer.h>

Classes

struct  Pointerstruct
 

Public Member Functions

 SolutionTransfer (const DoFHandler< dim, spacedim > &dof)
 
 ~SolutionTransfer ()
 
void clear ()
 
void prepare_for_pure_refinement ()
 
void prepare_for_coarsening_and_refinement (const std::vector< VectorType > &all_in)
 
void prepare_for_coarsening_and_refinement (const VectorType &in)
 
void refine_interpolate (const VectorType &in, VectorType &out) const
 
void interpolate (const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
 
void interpolate (const VectorType &in, VectorType &out) const
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static ::ExceptionBaseExcNotPrepared ()
 
static ::ExceptionBaseExcAlreadyPrepForRef ()
 
static ::ExceptionBaseExcAlreadyPrepForCoarseAndRef ()
 

Private Types

enum  PreparationState { none , pure_refinement , coarsening_and_refinement }
 

Private Attributes

SmartPointer< const DoFHandler< dim, spacedim >, SolutionTransfer< dim, VectorType, spacedim > > dof_handler
 
types::global_dof_index n_dofs_old
 
PreparationState prepared_for
 
std::vector< std::vector< types::global_dof_index > > indices_on_cell
 
std::map< std::pair< unsigned int, unsigned int >, Pointerstructcell_map
 
std::vector< std::vector< Vector< typename VectorType::value_type > > > dof_values_on_cell
 

Detailed Description

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
class SolutionTransfer< dim, VectorType, spacedim >

This class implements the transfer of a discrete FE function (e.g. a solution vector) from one mesh to another that is obtained from the first by a single refinement and/or coarsening step. During interpolation the vector is reinitialized to the new size and filled with the interpolated values. This class is used in the step-15, step-26, step-31, and step-33 tutorial programs. A version of this class that works on parallel triangulations is available as parallel::distributed::SolutionTransfer.

Usage

This class implements the algorithms in two different ways:

For deleting all stored data in SolutionTransfer and reinitializing it use the clear() function.

The template argument VectorType denotes the type of data container you want to transfer.

Interpolating in the presence of hanging nodes and boundary values

The interpolation onto the new mesh is a local operation, i.e., it interpolates onto the new mesh only. If that new mesh has hanging nodes, you will therefore get a solution that does not satisfy hanging node constraints. The same is true with boundary values: the interpolated solution will just be the interpolation of the old solution at the boundary, and this may or may not satisfy boundary values at newly introduced boundary nodes.

Consequently, you may have to apply hanging node or boundary value constraints after interpolation. step-15 and step-26 have examples of dealing with this.

Implementation

Interaction with hanging nodes

This class does its best to represent on the new mesh the finite element function that existed on the old mesh, but this may lead to situations where the function on the new mesh is no longer conforming at hanging nodes. To this end, consider a situation of a twice refined mesh that started with a single square cell (i.e., we now have 16 cells). Consider also that we coarsen 4 of the cells back to the first refinement level. In this case, we end up with a mesh that will look as follows if we were to use a \(Q_1\) element:

The process of interpolating from the old to the new mesh would imply that the values of the finite element function will not change on all of the cells that remained as they are (i.e., the fine cells) but that on the coarse cell at the top right, the four values at the vertices are obtained by interpolating down from its former children. If the original function was not linear, this implies that the marked hanging nodes will retain their old values which, in general, will not lead to a continuous function along the corresponding edges. In other words, the solution vector obtained after SolutionTransfer::interpolate() does not satisfy hanging node constraints: it corresponds to the pointwise interpolation, but not to the interpolation onto the new finite element space that contains constraints from hanging nodes.

Whether this is a problem you need to worry about or not depends on your application. The situation is easily corrected, of course, by applying AffineConstraints::distribute() to your solution vector after transfer, using a constraints object computed on the new DoFHandler object (you probably need to create this object anyway if you have hanging nodes). This is also what is done, for example, in step-15.

Note
This situation can only happen if you do coarsening. If all cells remain as they are or are refined, then SolutionTransfer::interpolate() computes a new vector of nodal values, but the function represented is of course exactly the same because the old finite element space is a subspace of the new one. Thus, if the old function was conforming (i.e., satisfied hanging node constraints), then so does the new one, and it is not necessary to call AffineConstraints::distribute().

Implementation in the context of hp-finite elements

In the case of DoFHandlers with hp-capabilities, nothing defines which of the finite elements that are part of the hp::FECollection associated with the DoFHandler, should be considered on cells that are not active (i.e., that have children). This is because degrees of freedom are only allocated for active cells and, in fact, it is not allowed to set an active FE index on non-active cells using DoFAccessor::set_active_fe_index().

It is, thus, not entirely natural what should happen if, for example, a few cells are coarsened away. This class then implements the following algorithm:

Note
In the context of hp-refinement, if cells are coarsened or the polynomial degree is lowered on some cells, then the old finite element space is not a subspace of the new space and you may run into the same situation as discussed above with hanging nodes. You may want to consider calling AffineConstraints::distribute() on the vector obtained by transferring the solution.

Definition at line 336 of file solution_transfer.h.

Member Enumeration Documentation

◆ PreparationState

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
enum SolutionTransfer::PreparationState
private

Declaration of PreparationState that denotes the three possible states of the SolutionTransfer: being prepared for 'pure refinement', prepared for 'coarsening and refinement' or not prepared.

Enumerator
none 

The SolutionTransfer is not yet prepared.

pure_refinement 

The SolutionTransfer is prepared for purely refinement.

coarsening_and_refinement 

The SolutionTransfer is prepared for coarsening and refinement.

Definition at line 485 of file solution_transfer.h.

Constructor & Destructor Documentation

◆ SolutionTransfer()

template<int dim, typename VectorType , int spacedim>
SolutionTransfer< dim, VectorType, spacedim >::SolutionTransfer ( const DoFHandler< dim, spacedim > &  dof)

Constructor, takes the current DoFHandler as argument.

Definition at line 46 of file solution_transfer.cc.

◆ ~SolutionTransfer()

template<int dim, typename VectorType , int spacedim>
SolutionTransfer< dim, VectorType, spacedim >::~SolutionTransfer

Destructor

Definition at line 65 of file solution_transfer.cc.

Member Function Documentation

◆ clear()

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::clear

Reinit this class to the state that it has directly after calling the Constructor

Definition at line 74 of file solution_transfer.cc.

◆ prepare_for_pure_refinement()

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::prepare_for_pure_refinement

Prepares the SolutionTransfer for pure refinement. It stores the dof indices of each cell. After calling this function only calling the refine_interpolate functions is allowed.

Definition at line 87 of file solution_transfer.cc.

◆ prepare_for_coarsening_and_refinement() [1/2]

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::prepare_for_coarsening_and_refinement ( const std::vector< VectorType > &  all_in)

Prepares the SolutionTransfer for coarsening and refinement. It stores the dof indices of each cell and stores the dof values of the vectors in all_in in each cell that'll be coarsened. all_in includes all vectors that are to be interpolated onto the new (refined and/or coarsenend) grid.

Definition at line 280 of file solution_transfer.cc.

◆ prepare_for_coarsening_and_refinement() [2/2]

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::prepare_for_coarsening_and_refinement ( const VectorType &  in)

Same as previous function but for only one discrete function to be interpolated.

Definition at line 436 of file solution_transfer.cc.

◆ refine_interpolate()

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::refine_interpolate ( const VectorType &  in,
VectorType &  out 
) const

This function interpolates the discrete function in, which is a vector on the grid before the refinement, to the function out which then is a vector on the refined grid. It assumes the vectors having the right sizes (i.e. in.size()==n_dofs_old, out.size()==n_dofs_refined)

Calling this function is allowed only if prepare_for_pure_refinement is called and the refinement is executed before. Multiple calling of this function is allowed. e.g. for interpolating several functions.

Definition at line 136 of file solution_transfer.cc.

◆ interpolate() [1/2]

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::interpolate ( const std::vector< VectorType > &  all_in,
std::vector< VectorType > &  all_out 
) const

This function interpolates the discrete functions that are stored in all_in onto the refined and/or coarsenend grid. It assumes the vectors in all_in denote the same vectors as in all_in as parameter of prepare_for_refinement_and_coarsening(all_in). However, there is no way of verifying this internally, so be careful here.

Calling this function is allowed only if first Triangulation::prepare_coarsening_and_refinement, second SolutionTransfer::prepare_for_coarsening_and_refinement, an then third Triangulation::execute_coarsening_and_refinement are called before. Multiple calling of this function is NOT allowed. Interpolating several functions can be performed in one step.

The number of output vectors is assumed to be the same as the number of input vectors. Also, the sizes of the output vectors are assumed to be of the right size (n_dofs_refined). Otherwise an assertion will be thrown.

Definition at line 447 of file solution_transfer.cc.

◆ interpolate() [2/2]

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::interpolate ( const VectorType &  in,
VectorType &  out 
) const

Same as the previous function. It interpolates only one function. It assumes the vectors having the right sizes (i.e. in.size()==n_dofs_old, out.size()==n_dofs_refined)

Multiple calling of this function is NOT allowed. Interpolating several functions can be performed in one step by using interpolate (all_in, all_out)

Definition at line 592 of file solution_transfer.cc.

◆ memory_consumption()

template<int dim, typename VectorType , int spacedim>
std::size_t SolutionTransfer< dim, VectorType, spacedim >::memory_consumption

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 611 of file solution_transfer.cc.

Member Data Documentation

◆ dof_handler

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SmartPointer<const DoFHandler<dim, spacedim>, SolutionTransfer<dim, VectorType, spacedim> > SolutionTransfer< dim, VectorType, spacedim >::dof_handler
private

Pointer to the degree of freedom handler to work with.

Definition at line 473 of file solution_transfer.h.

◆ n_dofs_old

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
types::global_dof_index SolutionTransfer< dim, VectorType, spacedim >::n_dofs_old
private

Stores the number of DoFs before the refinement and/or coarsening.

Definition at line 478 of file solution_transfer.h.

◆ prepared_for

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
PreparationState SolutionTransfer< dim, VectorType, spacedim >::prepared_for
private

Definition of the respective variable.

Definition at line 504 of file solution_transfer.h.

◆ indices_on_cell

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
std::vector<std::vector<types::global_dof_index> > SolutionTransfer< dim, VectorType, spacedim >::indices_on_cell
private

Is used for prepare_for_refining (of course also for repare_for_refining_and_coarsening) and stores all dof indices of the cells that'll be refined

Definition at line 512 of file solution_transfer.h.

◆ cell_map

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
std::map<std::pair<unsigned int, unsigned int>, Pointerstruct> SolutionTransfer< dim, VectorType, spacedim >::cell_map
private

Map mapping from level and index of cell to the Pointerstructs (cf. there). This map makes it possible to keep all the information needed to transfer the solution inside this object rather than using user pointers of the Triangulation for this purpose.

Definition at line 559 of file solution_transfer.h.

◆ dof_values_on_cell

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
std::vector<std::vector<Vector<typename VectorType::value_type> > > SolutionTransfer< dim, VectorType, spacedim >::dof_values_on_cell
private

Is used for prepare_for_refining_and_coarsening The interpolated dof values of all cells that'll be coarsened will be stored in this vector.

Definition at line 566 of file solution_transfer.h.


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