deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
\(\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 | Private Member Functions | Private Attributes | List of all members
SolutionTransfer< dim, VectorType, spacedim > Class Template Reference

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

Public Member Functions

 SolutionTransfer (const DoFHandler< dim, spacedim > &dof_handler, const bool average_values=false)
 
 ~SolutionTransfer ()=default
 
void prepare_for_coarsening_and_refinement (const std::vector< const VectorType * > &all_in)
 
void prepare_for_coarsening_and_refinement (const std::vector< VectorType > &all_in)
 
void prepare_for_coarsening_and_refinement (const VectorType &in)
 
void interpolate (std::vector< VectorType * > &all_out)
 
void interpolate (std::vector< VectorType > &all_out)
 
void interpolate (VectorType &out)
 
void prepare_for_serialization (const VectorType &in)
 
void prepare_for_serialization (const std::vector< const VectorType * > &all_in)
 
void deserialize (VectorType &in)
 
void deserialize (std::vector< VectorType * > &all_in)
 
void clear ()
 

Private Member Functions

std::vector< char > pack_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status)
 
void unpack_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range, std::vector< VectorType * > &all_out, VectorType &valence)
 
void register_data_attach ()
 

Private Attributes

ObserverPointer< const DoFHandler< dim, spacedim >, SolutionTransfer< dim, VectorType, spacedim > > dof_handler
 
const bool average_values
 
std::vector< const VectorType * > input_vectors
 
unsigned int handle
 

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 filled with the interpolated values. This class is used in the step-15, step-26, step-31, and step-33 tutorial programs. This class works both for serial and distributed meshes.

Usage

Note
It is important to note, that if you use more than one SolutionTransfer object at the same time, that the calls to prepare_*() and interpolate()/deserialize() need to be in the same order.

Transferring a solution

Here VectorType is your favorite vector type, e.g. Vector, LinearAlgebra::distributed::Vector, PETScWrappers::MPI::Vector, TrilinosWrappers::MPI::Vector, or corresponding block vectors.

// flag some cells for refinement and coarsening, e.g.
error_indicators,
0.3,
0.05);
// Use parallel::distributed::GridRefinement for distributed triangulations.
// prepare the triangulation,
tria.prepare_coarsening_and_refinement();
// prepare the SolutionTransfer object for coarsening and refinement
// and give the solution vector that we intend to interpolate later,
soltrans.prepare_for_coarsening_and_refinement(solution);
// actually execute the refinement,
tria.execute_coarsening_and_refinement();
// redistribute dofs,
dof_handler.distribute_dofs(fe);
// adjust the vector size so that it can hold all dofs
// (for distributed grids, see instructions below),
solution.reinit(...);
// and interpolate the solution
soltrans.interpolate(solution);
ObserverPointer< const DoFHandler< dim, spacedim >, SolutionTransfer< dim, VectorType, spacedim > > dof_handler
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max(), const VectorTools::NormType norm_type=VectorTools::L1_norm)

Usage on distributed grids

If the grid is distributed, it is important to note that the old solution(s) must be copied to one that also provides access to the locally relevant DoF values (these values required for the interpolation process):

// Create initial indexsets pertaining to the grid before refinement
const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
const IndexSet locally_relevant_dofs =
// The solution vector only knows about locally owned DoFs
solution.reinit(locally_owned_dofs,
mpi_communicator);
...
// Transfer solution to vector that provides access to
// locally relevant DoFs
TrilinosWrappers::MPI::Vector old_solution;
old_solution.reinit(locally_owned_dofs,
locally_relevant_dofs,
mpi_communicator);
old_solution = solution;
// Initialize SolutionTransfer object
soltrans.prepare_for_coarsening_and_refinement(old_solution);
...
// Refine grid
// Recreate locally_owned_dofs and locally_relevant_dofs index sets
...
solution.reinit(locally_owned_dofs, mpi_communicator);
soltrans.interpolate(solution);
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)

Note on ghost elements

In a parallel computation PETSc or Trilinos vector may contain ghost elements or not. For reading in information with prepare_for_coarsening_and_refinement() or prepare_for_serialization() you need to supply vectors with ghost elements, so that all locally_active elements can be read. On the other hand, ghosted vectors are generally not writable, so for calls to interpolate() or deserialize() you need to supply distributed vectors without ghost elements. More precisely, during interpolation the current algorithm writes into all locally active degrees of freedom.

Different from PETSc and Trilinos vectors, LinearAlgebra::distributed::Vector allows writing into ghost elements. For a ghosted vector the interpolation step can be accomplished via

interpolated_solution.zero_out_ghost_values();
soltrans.interpolate(interpolated_solution);
interpolated_solution.update_ghost_values();

Use for serialization

This class can be used to serialize and later deserialize a mesh with solution vectors to a file. If you use more than one DoFHandler and therefore more than one SolutionTransfer object, they need to be serialized and deserialized in the same order.

If vector has the locally relevant DoFs, serialization works as follows:

sol_trans.prepare_for_serialization(vector);
triangulation.save(filename);
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation

For deserialization the vector needs to be a distributed vector (without ghost elements):

//[create coarse mesh...]
triangulation.load(filename);
sol_trans.deserialize(distributed_vector);

Note on usage with DoFHandler with hp-capabilities

Since data on DoFHandler objects with hp-capabilities is associated with many different FiniteElement objects, each cell's data has to be processed with its corresponding future_fe_index. Further, if refinement is involved, data will be packed on the parent cell with its future_fe_index and unpacked later with the same index on its children. If cells get coarsened into one, data will be packed on the children with the least dominant finite element of their common subspace, and unpacked on the parent with this particular finite element (consult hp::FECollection::find_dominated_fe_extended() for more information).

Transferring a solution across refinement works exactly like in the non-hp-case. However, when considering serialization, we also have to store the active FE indices in an additional step. A code snippet demonstrating serialization with the SolutionTransfer class with DoFHandler objects with hp-capabilities is provided in the following. Here VectorType is your favorite vector type, e.g. Vector, LinearAlgebra::distributed::Vector, PETScWrappers::MPI::Vector, TrilinosWrappers::MPI::Vector, or corresponding block vectors.

If vector has the locally relevant DoFs, serialization works as follows:

sol_trans(hp_dof_handler);
hp_dof_handler.prepare_for_serialization_of_active_fe_indices();
sol_trans.prepare_for_serialization(vector);
triangulation.save(filename);

For deserialization the vector needs to be a distributed vector (without ghost elements):

//[create coarse mesh...]
triangulation.load(filename);
//[prepare identical fe_collection...]
// We need to introduce our dof_handler to the fe_collection
// before setting all active FE indices.
hp_dof_handler.deserialize_active_fe_indices();
hp_dof_handler.distribute_dofs(fe_collection);
sol_trans(hp_dof_handler);
sol_trans.deserialize(distributed_vector);

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 302 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_handler,
const bool  average_values = false 
)

Constructor.

Parameters
[in]dof_handlerThe DoFHandler on which all operations will happen. At the time when this constructor is called, the DoFHandler still points to the Triangulation before the refinement in question happens.
[in]average_valuesAverage the contribututions to the same DoF coming from different cells. Note: averaging requires an additional communication step, since the valence of the DoF has to be determined.

Definition at line 116 of file solution_transfer.cc.

◆ ~SolutionTransfer()

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SolutionTransfer< dim, VectorType, spacedim >::~SolutionTransfer ( )
default

Destructor.

Member Function Documentation

◆ prepare_for_coarsening_and_refinement() [1/3]

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

Prepare the current object for coarsening and refinement. It stores pointers to the vectors in all_in to be interpolated onto the new (refined and/or coarsened) grid, and registers this object for data transfer on the grid.

Definition at line 128 of file solution_transfer.cc.

◆ prepare_for_coarsening_and_refinement() [2/3]

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

Same as above but without pointers.

Definition at line 148 of file solution_transfer.cc.

◆ prepare_for_coarsening_and_refinement() [3/3]

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

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

Definition at line 186 of file solution_transfer.cc.

◆ interpolate() [1/3]

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

Interpolate the data previously stored in this object before the mesh was refined or coarsened onto the current set of cells. Do so for each of the vectors provided to prepare_for_coarsening_and_refinement() and write the result into the given set of vectors.

Definition at line 242 of file solution_transfer.cc.

◆ interpolate() [2/3]

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

Same as above but without pointers.

Definition at line 316 of file solution_transfer.cc.

◆ interpolate() [3/3]

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

Same as the previous function. It interpolates only one function. It assumes the vector has the right size (i.e. 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_out)

Definition at line 331 of file solution_transfer.cc.

◆ prepare_for_serialization() [1/2]

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

Prepare the serialization of the given vector. The serialization is done by Triangulation::save(). The given vector needs all information on the locally active DoFs (it must be ghosted). See documentation of this class for more information.

Definition at line 197 of file solution_transfer.cc.

◆ prepare_for_serialization() [2/2]

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

Same as the function above, only for a list of vectors.

Definition at line 208 of file solution_transfer.cc.

◆ deserialize() [1/2]

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

Execute the deserialization of the given vector. This needs to be done after calling Triangulation::load(). The given vector must be a fully distributed vector without ghost elements. See documentation of this class for more information.

Definition at line 218 of file solution_transfer.cc.

◆ deserialize() [2/2]

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

Same as the function above, only for a list of vectors.

Definition at line 228 of file solution_transfer.cc.

◆ 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 505 of file solution_transfer.cc.

◆ pack_callback()

template<int dim, typename VectorType , int spacedim>
std::vector< char > SolutionTransfer< dim, VectorType, spacedim >::pack_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellStatus  status 
)
private

A callback function used to pack the data on the current mesh into objects that can later be retrieved after refinement, coarsening and repartitioning.

Definition at line 341 of file solution_transfer.cc.

◆ unpack_callback()

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::unpack_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellStatus  status,
const boost::iterator_range< std::vector< char >::const_iterator > &  data_range,
std::vector< VectorType * > &  all_out,
VectorType &  valence 
)
private

A callback function used to unpack the data on the current mesh that has been packed up previously on the mesh before refinement, coarsening and repartitioning.

Definition at line 408 of file solution_transfer.cc.

◆ register_data_attach()

template<int dim, typename VectorType , int spacedim>
void SolutionTransfer< dim, VectorType, spacedim >::register_data_attach ( )
private

Registers the pack_callback() function to the Triangulation that has been assigned to the DoFHandler class member and stores the returning handle.

Definition at line 163 of file solution_transfer.cc.

Member Data Documentation

◆ dof_handler

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
ObserverPointer<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 420 of file solution_transfer.h.

◆ average_values

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
const bool SolutionTransfer< dim, VectorType, spacedim >::average_values
private

Flag indicating if averaging should be performed.

Definition at line 425 of file solution_transfer.h.

◆ input_vectors

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
std::vector<const VectorType *> SolutionTransfer< dim, VectorType, spacedim >::input_vectors
private

A vector that stores pointers to all the vectors we are supposed to copy over from the old to the new mesh.

Definition at line 431 of file solution_transfer.h.

◆ handle

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
unsigned int SolutionTransfer< dim, VectorType, spacedim >::handle
private

The handle that the Triangulation has assigned to this object with which we can access our memory offset and our pack function.

Definition at line 437 of file solution_transfer.h.


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