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 Attributes | List of all members
CellDataTransfer< dim, spacedim, VectorType > Class Template Reference

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

Public Member Functions

 CellDataTransfer (const Triangulation< dim, spacedim > &triangulation, 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 ()
 
void unpack (const VectorType &in, VectorType &out)
 

Private Types

using value_type = typename VectorType::value_type
 

Private Attributes

SmartPointer< const Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation
 
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_indices)> coarsening_strategy
 
std::map< const typename Triangulation< dim, spacedim >::cell_iterator, const unsigned intpersisting_cells_active_index
 
std::map< const typename Triangulation< dim, spacedim >::cell_iterator, const unsigned intrefined_cells_active_index
 
std::map< const typename Triangulation< dim, spacedim >::cell_iterator, const std::set< unsigned int > > coarsened_cells_active_index
 
unsigned int n_active_cells_pre
 

Detailed Description

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

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

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

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();
// actually execute the refinement,
triangulation.execute_coarsening_and_refinement();
// unpack transferred data,
Vector<double> transferred_data(triangulation.n_active_cells());
cell_data_trans.unpack(data_to_transfer, transferred_data);
SmartPointer< const Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation

When using a parallel::shared::Triangulation, we need to ensure that we have the global data available in our local vector before refinement happened. We can achieve this as follows:

Vector<double> data_to_transfer(triangulation.n_active_cells());
//[fill data_to_transfer with cell-wise values...]
distributed_data_to_transfer(mpi_communicator,
triangulation.n_active_cells(),
triangulation.n_locally_owned_active_cells());
for (const auto &cell : triangulation.active_cell_iterators())
if (cell->is_locally_owned())
{
const unsigned int index = cell->active_cell_index();
distributed_data_to_transfer(index) = data_to_transfer(index);
}
distributed_data_to_transfer.compress(VectorOperation::insert);
data_to_transfer = distributed_data_to_transfer;

For the parallel distributed case, a designated class parallel::distributed::CellDataTransfer is available. Please refer to this particular class when using a parallel::distributed::Triangulation.

Note
See the documentation of SolutionTransfer for matching code snippets for transfer.

Definition at line 107 of file cell_data_transfer.h.

Member Typedef Documentation

◆ value_type

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
using 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 113 of file cell_data_transfer.h.

Constructor & Destructor Documentation

◆ CellDataTransfer()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
CellDataTransfer< dim, spacedim, VectorType >::CellDataTransfer ( const Triangulation< dim, spacedim > &  triangulation,
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]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()

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
void CellDataTransfer< dim, spacedim, VectorType >::prepare_for_coarsening_and_refinement ( )

Prepare the current object for coarsening and refinement.

Stores the active_cell_indices of all active cells on the associated triangulation and attribute them to either persisting, refined or coarsened cells.

◆ unpack()

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

Transfer the information from the previous mesh to the updated one.

Data from the previous mesh supplied by in will be transferred to the updated mesh and stored in out. out has to provide enough space to hold the transferred data, i.e. has to be of size triangulation.n_active_cells().

Member Data Documentation

◆ triangulation

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

Pointer to the triangulation to work with.

Definition at line 165 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)> CellDataTransfer< dim, spacedim, VectorType >::refinement_strategy
private

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

Definition at line 174 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_indices)> 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 183 of file cell_data_transfer.h.

◆ persisting_cells_active_index

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
std::map<const typename Triangulation<dim, spacedim>::cell_iterator, const unsigned int> CellDataTransfer< dim, spacedim, VectorType >::persisting_cells_active_index
private

Container to temporarily store the iterator and active cell index of cells that persist.

Definition at line 191 of file cell_data_transfer.h.

◆ refined_cells_active_index

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
std::map<const typename Triangulation<dim, spacedim>::cell_iterator, const unsigned int> CellDataTransfer< dim, spacedim, VectorType >::refined_cells_active_index
private

Container to temporarily store the iterator and active cell index of cells that will be refined.

Definition at line 199 of file cell_data_transfer.h.

◆ coarsened_cells_active_index

template<int dim, int spacedim = dim, typename VectorType = Vector<double>>
std::map<const typename Triangulation<dim, spacedim>::cell_iterator, const std::set<unsigned int> > CellDataTransfer< dim, spacedim, VectorType >::coarsened_cells_active_index
private

Container to temporarily store the iterator of parent cells that will remain after coarsening along with the active cell indices of the corresponding children cells.

Definition at line 208 of file cell_data_transfer.h.

◆ n_active_cells_pre

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

Number of active cells on the initial triangulation that has not been refined yet.

It will be set in prepare_for_coarsening_and_refinement() and used to validate user inputs after refinement happened (only in debug mode).

Definition at line 217 of file cell_data_transfer.h.


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