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::CellWeights< dim, spacedim > Class Template Reference

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

Public Member Functions

 CellWeights (const ::hp::DoFHandler< dim, spacedim > &dof_handler)
 
 ~CellWeights ()
 
void register_constant_weighting (const unsigned int factor=1000)
 
void register_ndofs_weighting (const unsigned int factor=1000)
 
void register_ndofs_squared_weighting (const unsigned int factor=1000)
 
void register_custom_weighting (const std::function< unsigned int(const FiniteElement< dim, spacedim > &, const typename hp::DoFHandler< dim, spacedim >::cell_iterator &)> custom_function)
 

Private Member Functions

unsigned int weight_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status)
 

Private Attributes

SmartPointer< const ::hp::DoFHandler< dim, spacedim >, CellWeightsdof_handler
 
SmartPointer< const parallel::Triangulation< dim, spacedim >, CellWeightstriangulation
 
std::function< unsigned int(const FiniteElement< dim, spacedim > &, const typename hp::DoFHandler< dim, spacedim >::cell_iterator &)> weighting_function
 
boost::signals2::connection tria_listener
 

Detailed Description

template<int dim, int spacedim = dim>
class parallel::CellWeights< dim, spacedim >

Anytime a parallel::Triangulation is repartitioned, either upon request or by refinement/coarsening, cells will be distributed amongst all subdomains to achieve an equally balanced workload. If the workload per cell varies, which is in general the case for hp::DoFHandler objects, we can take that into account by introducing individual weights for different cells.

This class allows computing these weights for load balancing by consulting the FiniteElement that is associated with each cell of a hp::DoFHandler. One can choose from predefined weighting algorithms provided by this class or provide a custom one. The chosen weighting function will be connected to the corresponding signal of the linked parallel::Triangulation via callback.

An object of this class needs to exist for every DoFHandler associated with the Triangulation we work on to achieve satisfying work balancing results.

A small code snippet follows explaining how to achieve each cell being weighted by its current number of degrees of freedom.

parallel::CellWeights<dim, spacedim> cell_weights(hp_dof_handler);
cell_weights.register_ndofs_weighting();

The weighting function can be changed anytime. Even more ambitious approaches are possible by submitting customized functions, e.g.

cell_weights.register_custom_weighting(
[](const FiniteElement<dim, spacedim> &active_fe,
-> unsigned int {
return 1000 * std::pow(active_fe.dofs_per_cell, 1.6);
});

The returned value has to be an unsigned integer and is thus limited by some large number. It is interpreted as the additional computational load of each cell. See Triangulation::Signals::cell_weight for a discussion on this topic.

Note
Be aware that this class connects the weight function to the Triangulation during its construction. If the Triangulation associated with the DoFHandler changes during the lifetime of the latter, an assertion will be triggered in the weight_callback() function. It is recommended to create a separate object in this case and to destroy the previous one.
Author
Marc Fehling, 2018

Definition at line 81 of file cell_weights.h.

Constructor & Destructor Documentation

◆ CellWeights()

template<int dim, int spacedim = dim>
parallel::CellWeights< dim, spacedim >::CellWeights ( const ::hp::DoFHandler< dim, spacedim > &  dof_handler)

Constructor.

Parameters
[in]dof_handlerThe hp::DoFHandler which will be used to determine each cell's finite element.

Definition at line 27 of file cell_weights.cc.

◆ ~CellWeights()

template<int dim, int spacedim>
parallel::CellWeights< dim, spacedim >::~CellWeights ( )

Destructor.

Definition at line 57 of file cell_weights.cc.

Member Function Documentation

◆ register_constant_weighting()

template<int dim, int spacedim>
void parallel::CellWeights< dim, spacedim >::register_constant_weighting ( const unsigned int  factor = 1000)

Choose a constant weight factor on each cell.

Definition at line 66 of file cell_weights.cc.

◆ register_ndofs_weighting()

template<int dim, int spacedim>
void parallel::CellWeights< dim, spacedim >::register_ndofs_weighting ( const unsigned int  factor = 1000)

Choose a weight for each cell that is proportional to its number of degrees of freedom with a factor factor.

Definition at line 78 of file cell_weights.cc.

◆ register_ndofs_squared_weighting()

template<int dim, int spacedim>
void parallel::CellWeights< dim, spacedim >::register_ndofs_squared_weighting ( const unsigned int  factor = 1000)

Choose a weight for each cell that is proportional to its number of degrees of freedom squared with a factor factor.

Definition at line 90 of file cell_weights.cc.

◆ register_custom_weighting()

template<int dim, int spacedim>
void parallel::CellWeights< dim, spacedim >::register_custom_weighting ( const std::function< unsigned int(const FiniteElement< dim, spacedim > &, const typename hp::DoFHandler< dim, spacedim >::cell_iterator &)>  custom_function)

Register a custom weight for each cell by providing a function as a parameter.

Parameters
[in]custom_functionThe custom weighting function returning the weight of each cell as an unsigned integer. It is required to have two arguments, namely the FiniteElement that will be active on the particular cell, and the cell itself of type hp::DoFHandler::cell_iterator. We require both to make sure to get the right active FiniteElement on each cell in case that we coarsen the Triangulation.

Definition at line 104 of file cell_weights.cc.

◆ weight_callback()

template<int dim, int spacedim>
unsigned int parallel::CellWeights< dim, spacedim >::weight_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename Triangulation< dim, spacedim >::CellStatus  status 
)
private

A callback function that will be attached to the cell_weight signal of the Triangulation, that is a member of the DoFHandler. Ultimately returns the weight for each cell, determined by the weighting_function member.

Definition at line 117 of file cell_weights.cc.

Member Data Documentation

◆ dof_handler

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

Pointer to the degree of freedom handler.

Definition at line 141 of file cell_weights.h.

◆ triangulation

template<int dim, int spacedim = dim>
SmartPointer<const parallel::Triangulation<dim, spacedim>, CellWeights> parallel::CellWeights< dim, spacedim >::triangulation
private

Pointer to the Triangulation associated with the degree of freedom handler.

We store both to make sure to always work on the correct combination of both.

Definition at line 151 of file cell_weights.h.

◆ weighting_function

template<int dim, int spacedim = dim>
std::function<unsigned int( const FiniteElement<dim, spacedim> &, const typename hp::DoFHandler<dim, spacedim>::cell_iterator &)> parallel::CellWeights< dim, spacedim >::weighting_function
private

Function that will determine each cell's weight.

Can be set using the register_constant_weighting(), register_ndofs_weighting(), register_ndofs_squared_weighting(), and register_custom_weighting() member functions.

The function requires the active FiniteElement object on each cell as an argument, as well as the cell itself of type hp::DoFHandler::cell_iterator.

Definition at line 167 of file cell_weights.h.

◆ tria_listener

template<int dim, int spacedim = dim>
boost::signals2::connection parallel::CellWeights< dim, spacedim >::tria_listener
private

A connection to the Triangulation of the DoFHandler.

Definition at line 172 of file cell_weights.h.


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