Reference documentation for deal.II version 9.2.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\}}\)
Public Member Functions | Private Member Functions | Private Attributes | List of all members
parallel::distributed::ErrorPredictor< dim, spacedim > Class Template Reference

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

Public Member Functions

 ErrorPredictor (const hp::DoFHandler< dim, spacedim > &dof)
 
void prepare_for_coarsening_and_refinement (const std::vector< const Vector< float > * > &all_in, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
 
void prepare_for_coarsening_and_refinement (const Vector< float > &in, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
 
void unpack (std::vector< Vector< float > * > &all_out)
 
void unpack (Vector< float > &out)
 

Private Member Functions

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

Private Attributes

SmartPointer< const hp::DoFHandler< dim, spacedim >, ErrorPredictor< dim, spacedim > > dof_handler
 
std::vector< const Vector< float > * > error_indicators
 
unsigned int handle
 
double gamma_p
 
double gamma_h
 
double gamma_n
 

Detailed Description

template<int dim, int spacedim = dim>
class parallel::distributed::ErrorPredictor< dim, spacedim >

Predict how current error indicators will change after refinement and coarsening were to happen on the provided hp::DoFHandler in context of a parallel::distributed::Triangulation.

This algorithm follows the same logic as the error prediction algorithm presented by [51] and has been implemented in deal.II via the function hp::Refinement::predict_error(). See the documentation of this particular function for details on the algorithm and all featured formulas.

In theory, the combination of hp::Refinement::predict_error() to predict how the error will change with parallel::distributed::CellDataTransfer to transfer data across meshes forms the intended way to apply the prediction algorithm. However, p4est determines the details of grid refinement in the parallel distributed case; consequently, it yields more reliable and trustworthy results when we determine the predicted errors during the adaptation process. This is achieved with this class, similar to the parallel::distributed::SolutionTransfer and parallel::distributed::CellDataTransfer classes.

Before adaptation happens, one or multiple vectors containing error estimates for every cell have to be registered with this class using the prepare_for_coarsening_and_refinement() functions. Performing adaptation on the triangulation initiates the calculation of the error prediction based on how p4est performed grid adaptation. After the refinement process, the predicted errors need to be unpacked on the updated mesh via one of the unpack() functions.

The following code snippet demonstrates how to impose hp-adaptivity based on refinement history in a parallel distributed application:

// [initialization...]
Vector<float> predicted_error_per_cell(triangulation.n_active_cells());
for(unsigned int i = 0; i < triangulation.n_active_cells(); ++i)
predicted_error_per_cell[i] = std::numeric_limits<float>::infinity();
// [during each refinement step...]
// set h-adaptivity flags
Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
parallel::distributed::GridRefinement::
refine_and_coarsen_fixed_{number|fraction}(...);
// set p-adaptivity flags
hp_dof_handler,
estimated_error_per_cell,
predicted_error_per_cell
std::less<float>(),
std::less<float>());
hp::Refinement::{choose|force}_p_over_h(hp_dof_handler);
// perform adaptation and predict error for the subsequent adaptation
predictor.prepare_for_coarsening_and_refinement(estimated_error_per_cell);
triangulation.execute_coarsening_and_refinement();
// unpack predicted errors
predicted_error_per_cell.reinit(triangulation.n_active_cells());
predictor.unpack(predicted_error_per_cell);

For serialization of predicted errors, use the parallel::distributed::CellDataTransfer class.

Author
Marc Fehling, 2019 - 2020

Definition at line 127 of file error_predictor.h.

Constructor & Destructor Documentation

◆ ErrorPredictor()

template<int dim, int spacedim>
parallel::distributed::ErrorPredictor< dim, spacedim >::ErrorPredictor ( const hp::DoFHandler< dim, spacedim > &  dof)

Constructor.

Parameters
[in]dofThe 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 53 of file error_predictor.cc.

Member Function Documentation

◆ prepare_for_coarsening_and_refinement() [1/2]

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::prepare_for_coarsening_and_refinement ( const std::vector< const Vector< float > * > &  all_in,
const double  gamma_p = std::sqrt(0.4),
const double  gamma_h = 2.,
const double  gamma_n = 1. 
)

Prepare the current object for coarsening and refinement. all_in includes all vectors that are to be interpolated onto the new (refined and/or coarsened) grid.

As specified in the error prediction algorithm in hp::Refinement::predict_error(), the control parameters gamma_p, gamma_h, gamma_n can be used to influence the predicted errors in case of p-adaptation, h-adaptation, and no adapation, respectively. Their default values comply to the studies of [51].

Definition at line 72 of file error_predictor.cc.

◆ prepare_for_coarsening_and_refinement() [2/2]

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::prepare_for_coarsening_and_refinement ( const Vector< float > &  in,
const double  gamma_p = std::sqrt(0.4),
const double  gamma_h = 2.,
const double  gamma_n = 1. 
)

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

Definition at line 110 of file error_predictor.cc.

◆ unpack() [1/2]

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::unpack ( std::vector< Vector< float > * > &  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 124 of file error_predictor.cc.

◆ unpack() [2/2]

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::unpack ( Vector< float > &  out)

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

Multiple calling of this function is NOT allowed. Interpolating several functions can be performed in one step by using the above function.

Definition at line 156 of file error_predictor.cc.

◆ pack_callback()

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

Definition at line 166 of file error_predictor.cc.

◆ unpack_callback()

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::unpack_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename Triangulation< dim, spacedim >::CellStatus  status,
const boost::iterator_range< std::vector< char >::const_iterator > &  data_range,
std::vector< Vector< float > * > &  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 302 of file error_predictor.cc.

◆ register_data_attach()

template<int dim, int spacedim>
void parallel::distributed::ErrorPredictor< dim, spacedim >::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 89 of file error_predictor.cc.

Member Data Documentation

◆ dof_handler

template<int dim, int spacedim = dim>
SmartPointer<const hp::DoFHandler<dim, spacedim>, ErrorPredictor<dim, spacedim> > parallel::distributed::ErrorPredictor< dim, spacedim >::dof_handler
private

Pointer to the degree of freedom handler to work with.

Definition at line 197 of file error_predictor.h.

◆ error_indicators

template<int dim, int spacedim = dim>
std::vector<const Vector<float> *> parallel::distributed::ErrorPredictor< dim, spacedim >::error_indicators
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 203 of file error_predictor.h.

◆ handle

template<int dim, int spacedim = dim>
unsigned int parallel::distributed::ErrorPredictor< dim, 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 209 of file error_predictor.h.

◆ gamma_p

template<int dim, int spacedim = dim>
double parallel::distributed::ErrorPredictor< dim, spacedim >::gamma_p
private

Control parameters for the prediction algorithm.

Definition at line 214 of file error_predictor.h.

◆ gamma_h

template<int dim, int spacedim = dim>
double parallel::distributed::ErrorPredictor< dim, spacedim >::gamma_h
private

Definition at line 214 of file error_predictor.h.

◆ gamma_n

template<int dim, int spacedim = dim>
double parallel::distributed::ErrorPredictor< dim, spacedim >::gamma_n
private

Definition at line 214 of file error_predictor.h.


The documentation for this class was generated from the following files:
parallel::distributed::ErrorPredictor
Definition: error_predictor.h:127
KellyErrorEstimator::estimate
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
hp::Refinement
Definition: refinement.h:130
hp::Refinement::p_adaptivity_from_reference
void p_adaptivity_from_reference(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const Vector< Number > &references, const ComparisonFunction< typename identity< Number >::type > &compare_refine, const ComparisonFunction< typename identity< Number >::type > &compare_coarsen)
Definition: refinement.cc:458
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
Vector
Definition: mapping_q1_eulerian.h:32