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

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

Public Types

using WeightingFunction = std::function< unsigned int(const typename hp::DoFHandler< dim, spacedim >::cell_iterator &, const FiniteElement< dim, spacedim > &)>
 

Public Member Functions

 CellWeights (const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
 
 ~CellWeights ()
 
void reinit (const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
 
Deprecated functions
 CellWeights (const hp::DoFHandler< dim, spacedim > &dof_handler)
 
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)
 

Static Public Member Functions

static std::function< unsigned int(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status)> make_weighting_callback (const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
 
Selection of weighting functions
static WeightingFunction constant_weighting (const unsigned int factor=1000)
 
static WeightingFunction ndofs_weighting (const std::pair< float, float > &coefficients)
 
static WeightingFunction ndofs_weighting (const std::vector< std::pair< float, float >> &coefficients)
 

Static Private Member Functions

static unsigned int weighting_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const hp::DoFHandler< dim, spacedim > &dof_handler, const parallel::TriangulationBase< dim, spacedim > &triangulation, const WeightingFunction &weighting_function)
 

Private Attributes

boost::signals2::connection connection
 
Deprecated members
SmartPointer< const ::hp::DoFHandler< dim, spacedim >, CellWeightsdof_handler
 
SmartPointer< const parallel::TriangulationBase< dim, spacedim >, CellWeightstriangulation
 

Detailed Description

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

Anytime a parallel::TriangulationBase 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.

This class offers two different ways of connecting the chosen weighting function to the corresponding signal of the linked parallel::TriangulationBase. The recommended way involves creating an object of this class which will automatically take care of registering the weighting function upon creation and de-registering it once destroyed. An object of this class needs to exist for every DoFHandler associated with the Triangulation we work on to achieve satisfying work balancing results. The connected weighting function may be changed anytime using the CellWeights::reinit() function. The following code snippet demonstrates how to achieve each cell being weighted by its current number of degrees of freedom. We chose a factor of 1000 that corresponds to the initial weight each cell is assigned to upon creation.

On the other hand, you are also able to take care of handling the signal connection manually by using the static member function of this class. In this case, an analogous code example looks as follows.

boost::signals2::connection connection =
hp_dof_handler.get_triangulation().signals.cell_weight.connect(
hp_dof_handler,
{1000, 1}));
Note
See Triangulation::Signals::cell_weight for more information on weighting and load balancing.
Be aware that this class connects the weight function to the Triangulation during this class's constructor. If the Triangulation associated with the DoFHandler changes during the lifetime of the latter via hp::DoFHandler::initialize(), an assertion will be triggered in the weight_callback() function. Use CellWeights::reinit() to deregister the weighting function on the old Triangulation and connect it to the new one.
Author
Marc Fehling, 2018, 2019

Definition at line 87 of file cell_weights.h.

Member Typedef Documentation

◆ WeightingFunction

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

An alias that defines the characteristics of a function that can be used for weighting cells during load balancing.

Such weighting functions take as arguments an iterator to a cell and the future finite element that will be assigned to it after repartitioning. They return an unsigned integer which is interpreted as the cell's weight or, in other words, the additional computational load associated with it.

Definition at line 102 of file cell_weights.h.

Constructor & Destructor Documentation

◆ CellWeights() [1/2]

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

Constructor.

Parameters
[in]dof_handlerThe hp::DoFHandler which will be used to determine each cell's finite element.
[in]weighting_functionThe function that determines each cell's weight during load balancing.

Definition at line 28 of file cell_weights.cc.

◆ ~CellWeights()

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

Destructor.

Disconnects the function previously connected to the weighting signal.

Definition at line 38 of file cell_weights.cc.

◆ CellWeights() [2/2]

template<int dim, int spacedim>
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 231 of file cell_weights.cc.

Member Function Documentation

◆ reinit()

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

Connect a different weighting_function to the Triangulation associated with the dof_handler.

Disconnects the function previously connected to the weighting signal.

Definition at line 47 of file cell_weights.cc.

◆ make_weighting_callback()

template<int dim, int spacedim = dim>
std::function< unsigned int(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status)> parallel::CellWeights< dim, spacedim >::make_weighting_callback ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const WeightingFunction weighting_function 
)
static

Converts a weighting_function to a different type that qualifies as a callback function, which can be connected to a weighting signal of a Triangulation.

This function does not connect the converted function to the Triangulation associated with the dof_handler.

Definition at line 131 of file cell_weights.cc.

◆ constant_weighting()

template<int dim, int spacedim>
CellWeights< dim, spacedim >::WeightingFunction parallel::CellWeights< dim, spacedim >::constant_weighting ( const unsigned int  factor = 1000)
static

Choose a constant weight factor on each cell.

Definition at line 65 of file cell_weights.cc.

◆ ndofs_weighting() [1/2]

template<int dim, int spacedim>
CellWeights< dim, spacedim >::WeightingFunction parallel::CellWeights< dim, spacedim >::ndofs_weighting ( const std::pair< float, float > &  coefficients)
static

The pair of floating point numbers \((a,b)\) provided via coefficients determines the weight \(w_K\) of each cell \(K\) with \(n_K\) degrees of freedom in the following way:

\[ w_K = a \, n_K^b \]

The right hand side will be rounded to the nearest integer since cell weights are required to be integers.

Definition at line 78 of file cell_weights.cc.

◆ ndofs_weighting() [2/2]

template<int dim, int spacedim>
CellWeights< dim, spacedim >::WeightingFunction parallel::CellWeights< dim, spacedim >::ndofs_weighting ( const std::vector< std::pair< float, float >> &  coefficients)
static

The container coefficients provides pairs of floating point numbers \((a_i, b_i)\) that determine the weight \(w_K\) of each cell \(K\) with \(n_K\) degrees of freedom in the following way:

\[ w_K = \sum_i a_i \, n_K^{b_i} \]

The right hand side will be rounded to the nearest integer since cell weights are required to be integers.

Definition at line 102 of file cell_weights.cc.

◆ 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 250 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 272 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 294 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 316 of file cell_weights.cc.

◆ weighting_callback()

template<int dim, int spacedim = dim>
unsigned int parallel::CellWeights< dim, spacedim >::weighting_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename Triangulation< dim, spacedim >::CellStatus  status,
const hp::DoFHandler< dim, spacedim > &  dof_handler,
const parallel::TriangulationBase< dim, spacedim > &  triangulation,
const WeightingFunction weighting_function 
)
staticprivate

A callback function that will be connected to the cell_weight signal of the triangulation, to which the dof_handler is attached. Ultimately returns the weight for each cell, determined by the weighting_function provided as a parameter.

Definition at line 162 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 259 of file cell_weights.h.

◆ triangulation

template<int dim, int spacedim = dim>
SmartPointer<const parallel::TriangulationBase<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 270 of file cell_weights.h.

◆ connection

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

A connection to the corresponding cell_weight signal of the Triangulation which is attached to the DoFHandler.

Definition at line 280 of file cell_weights.h.


The documentation for this class was generated from the following files:
parallel::CellWeights
Definition: cell_weights.h:87
parallel::CellWeights::connection
boost::signals2::connection connection
Definition: cell_weights.h:280