Reference documentation for deal.II version 9.5.0
|
#include <deal.II/distributed/cell_weights.h>
Public Types | |
using | WeightingFunction = std::function< unsigned int(const typename DoFHandler< dim, spacedim >::cell_iterator &, const FiniteElement< dim, spacedim > &)> |
Public Member Functions | |
CellWeights ()=default | |
CellWeights (const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function) | |
~CellWeights () | |
void | reinit (const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_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 DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function) |
Selection of weighting functions | |
static WeightingFunction | constant_weighting (const unsigned int factor=1) |
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 DoFHandler< dim, spacedim > &dof_handler, const parallel::TriangulationBase< dim, spacedim > &triangulation, const WeightingFunction &weighting_function) |
Private Attributes | |
boost::signals2::connection | connection |
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 DoFHandler objects with hp-capabilities, 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 DoFHandler. One can choose from predefined weighting algorithms provided by this class or provide a custom one.
If the associated DoFHandler has not been initialized yet, i.e., its hp::FECollection is empty, all cell weights will be evaluated as zero.
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.
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.
The use of this class is demonstrated in step-75.
Definition at line 100 of file cell_weights.h.
using parallel::CellWeights< dim, spacedim >::WeightingFunction = std::function< unsigned int(const typename 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 113 of file cell_weights.h.
|
default |
Constructor.
No weighting function will be connected yet. Please call reinit().
parallel::CellWeights< dim, spacedim >::CellWeights | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
const WeightingFunction & | weighting_function | ||
) |
Constructor.
[in] | dof_handler | The DoFHandler which will be used to determine each cell's finite element. |
[in] | weighting_function | The function that determines each cell's weight during load balancing. |
Definition at line 29 of file cell_weights.cc.
parallel::CellWeights< dim, spacedim >::~CellWeights |
Destructor.
Disconnects the function previously connected to the weighting signal.
Definition at line 39 of file cell_weights.cc.
void parallel::CellWeights< dim, spacedim >::reinit | ( | const 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 48 of file cell_weights.cc.
|
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 132 of file cell_weights.cc.
|
static |
Choose a constant weight factor
on each cell.
Definition at line 66 of file cell_weights.cc.
|
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.
|
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.
|
staticprivate |
A callback function that will be connected to the 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. Returns zero if dof_handler
has not been initialized yet.
Definition at line 165 of file cell_weights.cc.
|
private |
A connection to the corresponding weight
signal of the Triangulation which is attached to the DoFHandler.
Definition at line 210 of file cell_weights.h.