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
Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
parallel::distributed::CellDataTransfer< dim, spacedim, VectorType > Class Template Reference

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

Inheritance diagram for parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >:
[legend]

Public Member Functions

 CellDataTransfer (const parallel::distributed::Triangulation< dim, spacedim > &triangulation, const bool transfer_variable_size_data=false, const std::function< std::vector< value_type >(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)> refinement_strategy=&::AdaptationStrategies::Refinement::preserve< dim, spacedim, value_type >, const std::function< value_type(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)> coarsening_strategy=&::AdaptationStrategies::Coarsening::check_equality< dim, spacedim, value_type >)
 
void prepare_for_coarsening_and_refinement (const VectorType &in)
 
void prepare_for_coarsening_and_refinement (const std::vector< const VectorType * > &all_in)
 
void prepare_for_serialization (const VectorType &in)
 
void prepare_for_serialization (const std::vector< const VectorType * > &all_in)
 
void unpack (VectorType &out)
 
void unpack (std::vector< VectorType * > &all_out)
 
void deserialize (VectorType &out)
 
void deserialize (std::vector< VectorType * > &all_out)
 

Private Types

using value_type = typename VectorType::value_type
 

Private Member Functions

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

Private Attributes

SmartPointer< const parallel::distributed::Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation
 
const bool transfer_variable_size_data
 
const std::function< std::vector< value_type >(const typename Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)> refinement_strategy
 
const std::function< value_type(const typename Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)> coarsening_strategy
 
std::vector< const VectorType * > input_vectors
 
unsigned int handle
 

Detailed Description

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
class parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >

Transfer data that is associated with each active cell (like error indicators) while refining and/or coarsening a distributed triangulation and handle the necessary communication.

This class therefore does for cell-related information what parallel::distributed::SolutionTransfer does for the values of degrees of freedom defined on a parallel::distributed::Triangulation.

This class has been designed to operate on any kind of datatype that is serializable. A non-distributed container (like Vector or std::vector) has to be provided, which holds the cell-wise data in the same order as active cells are traversed. In other words, each entry corresponds to the cell with the same index CellAccessor::active_cell_index(), and the container has to be of size Triangulation::n_active_cells().

Transferring cell-wise data

The following code snippet demonstrates how to transfer cell-related data across refinement/coarsening of the registered triangulation.

// prepare the triangulation,
triangulation.prepare_coarsening_and_refinement();
// prepare the CellDataTransfer object for coarsening and refinement
// and give the cell data vector that we intend to unpack later,
Vector<double> data_to_transfer(triangulation.n_active_cells());
//[fill data_to_transfer with cell-wise values...]
cell_data_trans(triangulation);
cell_data_trans.prepare_for_coarsening_and_refinement(data_to_transfer);
// actually execute the refinement,
triangulation.execute_coarsening_and_refinement();
// unpack transferred data,
Vector<double> transferred_data(triangulation.n_active_cells());
cell_data_trans.unpack(transferred_data);
SmartPointer< const parallel::distributed::Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation

Use for serialization

This class can be used to serialize and later deserialize a distributed mesh with attached data to separate files.

For serialization, the following code snippet saves not only the triangulation itself, but also the cell-wise data attached:

Vector<double> data_to_transfer(triangulation.n_active_cells());
//[fill data_to_transfer with cell-wise values...]
cell_data_trans(triangulation);
cell_data_trans.prepare_for_serialization(data_to_transfer);
triangulation.save(filename);

Later, during deserialization, both triangulation and data can be restored as follows:

//[create coarse mesh...]
triangulation.load(filename);
cell_data_trans(triangulation);
Vector<double> transferred_data(triangulation.n_active_cells());
cell_data_trans.deserialize(transferred_data);
Note
If you use more than one object to transfer data via the parallel::distributed::Triangulation::register_data_attach() and parallel::distributed::Triangulation::notify_ready_for_unpack() interface with the aim of serialization, the calls to the corresponding prepare_for_serialization() and deserialize() functions need to happen in the same order, respectively. Classes relying on this interface are e.g. parallel::distributed::CellDataTransfer, parallel::distributed::SolutionTransfer, and Particles::ParticleHandler.
See the documentation of parallel::distributed::SolutionTransfer for matching code snippets for both transfer and serialization.

Definition at line 127 of file cell_data_transfer.h.

Member Typedef Documentation

◆ value_type

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
using parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::value_type = typename VectorType::value_type
private

An alias that defines the data type of provided container template.

Definition at line 133 of file cell_data_transfer.h.

Constructor & Destructor Documentation

◆ CellDataTransfer()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::CellDataTransfer ( const parallel::distributed::Triangulation< dim, spacedim > &  triangulation,
const bool  transfer_variable_size_data = false,
const std::function< std::vector< value_type >(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)>  refinement_strategy = &::AdaptationStrategies::Refinement::preserve< dim, spacedim, value_type >,
const std::function< value_type(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)>  coarsening_strategy = &::AdaptationStrategies::Coarsening::check_equality< dim, spacedim, value_type > 
)

Constructor.

Parameters
[in]triangulationThe triangulation on which all operations will happen. At the time when this constructor is called, the refinement in question has not happened yet.
[in]transfer_variable_size_dataSpecify whether your VectorType container stores values that differ in size. A varying amount of data may be packed per cell, if for example the underlying ValueType of the VectorType container is a container itself.
[in]refinement_strategyFunction deciding how data will be stored on refined cells from its parent cell.
[in]coarsening_strategyFunction deciding which data to store on a cell whose children will get coarsened into.

Member Function Documentation

◆ prepare_for_coarsening_and_refinement() [1/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::prepare_for_coarsening_and_refinement ( const VectorType &  in)

Prepare the current object for coarsening and refinement.

It registers the data transfer of in on the underlying triangulation. in includes data to be interpolated onto the new (refined and/or coarsened) grid. See documentation of this class for more information on how to use this functionality.

This function can be called only once for the specified container until data transfer has been completed. If multiple vectors shall be transferred via this class, use the function below.

◆ prepare_for_coarsening_and_refinement() [2/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::prepare_for_coarsening_and_refinement ( const std::vector< const VectorType * > &  all_in)

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

◆ prepare_for_serialization() [1/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::prepare_for_serialization ( const VectorType &  in)

Prepare the serialization of the given vector.

The serialization is done by Triangulation::save(). See documentation of this class for more information on how to use this functionality.

This function can be called only once for the specified container until data transfer has been completed. If multiple vectors shall be transferred via this class, use the function below.

◆ prepare_for_serialization() [2/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::prepare_for_serialization ( const std::vector< const VectorType * > &  all_in)

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

◆ unpack() [1/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::unpack ( VectorType &  out)

Unpack the information previously stored in this object before the mesh was refined or coarsened onto the current set of cells.

◆ unpack() [2/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::unpack ( std::vector< VectorType * > &  all_out)

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

◆ deserialize() [1/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::deserialize ( VectorType &  out)

Execute the deserialization of the stored information. This needs to be done after calling Triangulation::load().

◆ deserialize() [2/2]

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::deserialize ( std::vector< VectorType * > &  all_out)

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

◆ register_data_attach()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::register_data_attach ( )
private

Registers the pack_callback() function to the triangulation and stores the returning handle.

◆ pack_callback()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
std::vector< char > parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::pack_callback ( const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename parallel::distributed::Triangulation< dim, spacedim >::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.

◆ unpack_callback()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::unpack_callback ( const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename parallel::distributed::Triangulation< dim, spacedim >::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.

Member Data Documentation

◆ triangulation

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
SmartPointer<const parallel::distributed::Triangulation<dim, spacedim>, CellDataTransfer<dim, spacedim, VectorType> > parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::triangulation
private

Pointer to the triangulation to work with.

Definition at line 241 of file cell_data_transfer.h.

◆ transfer_variable_size_data

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
const bool parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::transfer_variable_size_data
private

Specifies if size of data to transfer varies from cell to cell.

Definition at line 246 of file cell_data_transfer.h.

◆ refinement_strategy

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
const std::function<std::vector<value_type>( const typename Triangulation<dim, spacedim>::cell_iterator &parent, const value_type parent_value)> parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::refinement_strategy
private

Function deciding how data will be stored on refined cells from its parent cell.

Definition at line 255 of file cell_data_transfer.h.

◆ coarsening_strategy

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
const std::function<value_type( const typename Triangulation<dim, spacedim>::cell_iterator &parent, const std::vector<value_type> &children_values)> parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::coarsening_strategy
private

Function deciding on how to process data from children to be stored on the parent cell.

Definition at line 264 of file cell_data_transfer.h.

◆ input_vectors

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
std::vector<const VectorType *> parallel::distributed::CellDataTransfer< dim, spacedim, VectorType >::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 270 of file cell_data_transfer.h.

◆ handle

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

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

Definition at line 276 of file cell_data_transfer.h.


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