Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Member Functions | Private Attributes | List of all members
parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType > Class Template Reference

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

Public Member Functions

 SolutionTransfer (const DoFHandlerType &dof)
 
 ~SolutionTransfer ()=default
 
void prepare_for_coarsening_and_refinement (const std::vector< const VectorType *> &all_in)
 
void prepare_for_coarsening_and_refinement (const VectorType &in)
 
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 prepare_serialization (const VectorType &in)
 
void prepare_serialization (const std::vector< const VectorType *> &all_in)
 
void deserialize (VectorType &in)
 
void deserialize (std::vector< VectorType *> &all_in)
 

Private Member Functions

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

Private Attributes

SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
 
std::vector< const VectorType * > input_vectors
 
unsigned int handle
 

Detailed Description

template<int dim, typename VectorType, typename DoFHandlerType = DoFHandler<dim>>
class parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType >

Transfer a discrete FE function (like a solution vector) by interpolation while refining and/or coarsening a distributed grid and handles the necessary communication.

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.

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.

Transferring a solution

Here VectorType is your favorite vector type, e.g. 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);
// 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);
// and interpolate the solution
VectorType interpolated_solution;
//create VectorType in the right size here
soltrans.interpolate(interpolated_solution);

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_ghosts();
soltrans.interpolate(interpolated_solution);
interpolated_solution.update_ghost_values();

As 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
IndexSet locally_owned_dofs, locally_relevant_dofs;
locally_owned_dofs = dof_handler.locally_owned_dofs();
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;
...
// Refine grid
// Recreate locally_owned_dofs and locally_relevant_dofs index sets
...
solution.reinit(locally_owned_dofs, mpi_communicator);
soltrans.refine_interpolate(old_solution, solution);

Use for Serialization

This class can be used to serialize and later deserialize a distributed 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(dof_handler);
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);
sol_trans(dof_handler);
sol_trans.deserialize (distributed_vector);

Note on usage with hp::DoFHandler

If an object of the hp::DoFHandler class is registered with an instantiation of this parallel::distributed::SolutionTransfer class, it is necessary to explicitly specify this in the template argument list of this class, i.e.:

parallel::distributed::SolutionsTransfer<dim, VectorType,
hp::DoFHandler<dim, spacedim>> sol_trans(hp_dof_handler);

Since data on hp::DoFHandler objects is associated with many different FiniteElement objects, each cell's data has to be processed with its correpsonding active_fe_index. Further, if refinement is involved, data will be packed on the parent cell with its active_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 dominating finite element of their common subspace, and unpacked on the parent with this particular finite element.

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 parallel::distributed::SolutionTransfer class with hp::DoFHandler objects is provided in the following. Here VectorType is your favorite vector type, e.g. PETScWrappers::MPI::Vector, TrilinosWrappers::MPI::Vector, or corresponding block vectors.

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

parallel::distributed::
SolutionTransfer<dim, VectorType, hp::DoFHandler<dim,spacedim>>
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...]
hp::DoFHandler<dim,spacedim> hp_dof_handler(triangulation);
// We need to introduce our dof_handler to the fe_collection
// before setting all active_fe_indices.
hp_dof_handler.set_fe(fe_collection);
hp_dof_handler.deserialize_active_fe_indices();
hp_dof_handler.distribute_dofs(fe_collection);
parallel::distributed::
SolutionTransfer<dim,VectorType,hp::DoFHandler<dim,spacedim>>
sol_trans(dof_handler);
sol_trans.deserialize(distributed_vector);

Interaction with hanging nodes

In essence, this class implements the same steps as does SolutionTransfer (though the implementation is entirely separate). Consequently, the same issue with hanging nodes and coarsening can happen with this class as happens with SolutionTransfer. See there for an extended discussion.

Author
Timo Heister, 2009-2011

Definition at line 230 of file solution_transfer.h.

Constructor & Destructor Documentation

◆ SolutionTransfer()

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

Constructor.

Parameters
[in]dofThe DoFHandler or hp::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.

Definition at line 50 of file solution_transfer.cc.

◆ ~SolutionTransfer()

template<int dim, typename VectorType, typename DoFHandlerType = DoFHandler<dim>>
parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType >::~SolutionTransfer ( )
default

Destructor.

Member Function Documentation

◆ prepare_for_coarsening_and_refinement() [1/2]

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

Prepare the current object 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 coarsened) grid.

Definition at line 68 of file solution_transfer.cc.

◆ prepare_for_coarsening_and_refinement() [2/2]

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::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 104 of file solution_transfer.cc.

◆ interpolate() [1/2]

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::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 180 of file solution_transfer.cc.

◆ interpolate() [2/2]

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

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

◆ prepare_for_serialization() [1/2]

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::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 115 of file solution_transfer.cc.

◆ prepare_for_serialization() [2/2]

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

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

Definition at line 126 of file solution_transfer.cc.

◆ prepare_serialization() [1/2]

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::prepare_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.

Deprecated:
Use parallel::distributed::SolutionTransfer::prepare_for_serialization() instead.

Definition at line 135 of file solution_transfer.cc.

◆ prepare_serialization() [2/2]

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

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

Deprecated:
Use parallel::distributed::SolutionTransfer::prepare_for_serialization() instead.

Definition at line 145 of file solution_transfer.cc.

◆ deserialize() [1/2]

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::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 155 of file solution_transfer.cc.

◆ deserialize() [2/2]

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

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

Definition at line 166 of file solution_transfer.cc.

◆ pack_callback()

template<int dim, typename VectorType , typename DoFHandlerType >
std::vector< char > SolutionTransfer< dim, VectorType, DoFHandlerType >::pack_callback ( const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &  cell,
const typename Triangulation< dim, DoFHandlerType::space_dimension >::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 229 of file solution_transfer.cc.

◆ unpack_callback()

template<int dim, typename VectorType , typename DoFHandlerType >
void SolutionTransfer< dim, VectorType, DoFHandlerType >::unpack_callback ( const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &  cell,
const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus  status,
const boost::iterator_range< std::vector< char >::const_iterator > &  data_range,
std::vector< VectorType *> &  all_out 
)
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 328 of file solution_transfer.cc.

◆ register_data_attach()

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

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

Definition at line 79 of file solution_transfer.cc.

Member Data Documentation

◆ dof_handler

template<int dim, typename VectorType, typename DoFHandlerType = DoFHandler<dim>>
SmartPointer<const DoFHandlerType, SolutionTransfer<dim, VectorType, DoFHandlerType> > parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType >::dof_handler
private

Pointer to the degree of freedom handler to work with.

Definition at line 352 of file solution_transfer.h.

◆ input_vectors

template<int dim, typename VectorType, typename DoFHandlerType = DoFHandler<dim>>
std::vector<const VectorType *> parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType >::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 358 of file solution_transfer.h.

◆ handle

template<int dim, typename VectorType, typename DoFHandlerType = DoFHandler<dim>>
unsigned int parallel::distributed::SolutionTransfer< dim, VectorType, DoFHandlerType >::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 364 of file solution_transfer.h.


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