Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Precompute weights for parallel::CellWeights. #16537

Open
wants to merge 6 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
5 changes: 5 additions & 0 deletions doc/news/changes/minor/20240124Fehling
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
New: parallel::CellWeights can now work with pre-computed weights.
This requires that the weights are determined only by information
about finite elements.
<br>
(Marc Fehling, 2024/01/24)
18 changes: 12 additions & 6 deletions examples/step-75/step-75.cc
Original file line number Diff line number Diff line change
Expand Up @@ -871,8 +871,8 @@ namespace Step75
LinearAlgebra::distributed::Vector<double> locally_relevant_solution;
LinearAlgebra::distributed::Vector<double> system_rhs;

std::unique_ptr<FESeries::Legendre<dim>> legendre;
std::unique_ptr<parallel::CellWeights<dim>> cell_weights;
std::unique_ptr<FESeries::Legendre<dim>> legendre;
parallel::CellWeights<dim> cell_weights;

Vector<float> estimated_error_per_cell;
Vector<float> hp_decision_indicators;
Expand Down Expand Up @@ -979,15 +979,21 @@ namespace Step75
// the form that $a (n_\text{dofs})^b$ with a provided pair of parameters
// $(a,b)$. We register such a function in the following.
//
// Since we are only using information about the number of degrees of
// freedom per cell, which is a quantity unique to every finite element, we
// can compute the weights in advance.
//
// For load balancing, efficient solvers like the one we use should scale
// linearly with the number of degrees of freedom owned. We set the
// parameters for cell weighting correspondingly: A weighting factor of $1$
// and an exponent of $1$ (see the definitions of the `weighting_factor` and
// `weighting_exponent` above).
cell_weights = std::make_unique<parallel::CellWeights<dim>>(
dof_handler,
parallel::CellWeights<dim>::ndofs_weighting(
{prm.weighting_factor, prm.weighting_exponent}));
const auto weighting_function = parallel::CellWeights<dim>::ndofs_weighting(
{prm.weighting_factor, prm.weighting_exponent});
const auto precomputed_weights =
parallel::CellWeights<dim>::precompute_weights(fe_collection,
weighting_function);
cell_weights.reinit(dof_handler, precomputed_weights);

// In h-adaptive applications, we ensure a 2:1 mesh balance by limiting the
// difference of refinement levels of neighboring cells to one. With the
Expand Down
150 changes: 113 additions & 37 deletions include/deal.II/distributed/cell_weights.h
Original file line number Diff line number Diff line change
Expand Up @@ -71,9 +71,20 @@ namespace parallel
* hp_dof_handler.get_triangulation().signals.weight.connect(
* parallel::CellWeights<dim, spacedim>::make_weighting_callback(
* hp_dof_handler,
* parallel::CellWeights<dim, spacedim>::ndofs_weighting({1, 1}));
* parallel::CellWeights<dim, spacedim>::ndofs_weighting({1, 1})));
* @endcode
*
* Oftentimes, the weighting function only requires information about the
* finite element (for example about the number of degrees of freedom per
* cell), but not the cell we are currently considering (even though the
* weighting function takes a cell iterator as argument). In this case, it is
* wasteful to re-compute the weight on every cell we visit; instead, the
* weights can be computed in advance for each finite element in use on the
* DoFHandler with precompute_weights(). The returned container can then be
* either passed to the corresponding constructor, reinit() or
* make_weighting_callback() function. Of course, in this case no other
* cell-specific information can be used in the computation of the weights.
*
* The use of this class is demonstrated in step-75.
*
* @note See Triangulation::Signals::weight for more information on
Expand Down Expand Up @@ -117,6 +128,58 @@ namespace parallel
unsigned int(const typename DoFHandler<dim, spacedim>::cell_iterator &,
const FiniteElement<dim, spacedim> &)>;

/**
* @name Selection of weighting functions
* @{
*/

/**
* Choose a constant weight @p factor on each cell.
*/
static WeightingFunction
constant_weighting(const unsigned int factor = 1);

/**
* The pair of floating point numbers $(a,b)$ provided via
* @p coefficients determines the weight $w_K$ of each cell $K$ with
* $n_K$ degrees of freedom in the following way: \f[ w_K =
* a \, n_K^b \f]
*
* The right hand side will be rounded to the nearest integer since cell
* weights are required to be integers.
*/
static WeightingFunction
ndofs_weighting(const std::pair<float, float> &coefficients);

/**
* The container @p 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: \f[ w_K =
* \sum_i a_i \, n_K^{b_i} \f]
*
* The right hand side will be rounded to the nearest integer since cell
* weights are required to be integers.
*/
static WeightingFunction
ndofs_weighting(const std::vector<std::pair<float, float>> &coefficients);

/**
* If your weights only require information about the finite element on each
* cell, consider pre-computing them with this function. All other
* cell-specific characteristics in determining weights will be omitted.
*
* Returns pre-computed weights for each finite element. Pass the container
* to the corresponding constructor, reinit() or make_weighting_callback()
* functions.
*/
static std::vector<unsigned int>
precompute_weights(const hp::FECollection<dim, spacedim> &fe_collection,
const WeightingFunction &weighting_function);

/**
* @}
*/

/**
* Constructor.
*
Expand All @@ -135,6 +198,17 @@ namespace parallel
CellWeights(const DoFHandler<dim, spacedim> &dof_handler,
const WeightingFunction &weighting_function);

/**
* Constructor.
*
* @param[in] dof_handler The DoFHandler which will be used to
* determine each cell's finite element.
* @param[in] precomputed_weights Weights for each finite element that
* determine each cell's weight during load balancing.
*/
CellWeights(const DoFHandler<dim, spacedim> &dof_handler,
const std::vector<unsigned int> &precomputed_weights);

/**
* Destructor.
*
Expand All @@ -152,6 +226,17 @@ namespace parallel
reinit(const DoFHandler<dim, spacedim> &dof_handler,
const WeightingFunction &weighting_function);

/**
* Connect a different weighting mechanism to the Triangulation
* associated with the @p dof_handler. Values in @p precomputed_weights will
* be used as weights for each finite element.
*
* Disconnects the function previously connected to the weighting signal.
*/
void
reinit(const DoFHandler<dim, spacedim> &dof_handler,
const std::vector<unsigned int> &precomputed_weights);

/**
* Converts a @p weighting_function to a different type that qualifies as
* a callback function, which can be connected to a weighting signal of a
Expand All @@ -167,43 +252,19 @@ namespace parallel
const WeightingFunction &weighting_function);

/**
* @name Selection of weighting functions
* @{
*/

/**
* Choose a constant weight @p factor on each cell.
*/
static WeightingFunction
constant_weighting(const unsigned int factor = 1);

/**
* The pair of floating point numbers $(a,b)$ provided via
* @p coefficients determines the weight $w_K$ of each cell $K$ with
* $n_K$ degrees of freedom in the following way: \f[ w_K =
* a \, n_K^b \f]
* Use @p precomputed_weights as weights for each finite element of
* @p dof_handler. Returns a callback function, which can be connected to a
* weighting signal of a Triangulation.
*
* The right hand side will be rounded to the nearest integer since cell
* weights are required to be integers.
*/
static WeightingFunction
ndofs_weighting(const std::pair<float, float> &coefficients);

/**
* The container @p 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: \f[ w_K =
* \sum_i a_i \, n_K^{b_i} \f]
*
* The right hand side will be rounded to the nearest integer since cell
* weights are required to be integers.
*/
static WeightingFunction
ndofs_weighting(const std::vector<std::pair<float, float>> &coefficients);

/**
* @}
* This function does <b>not</b> connect the callback function to the
* Triangulation associated with the @p dof_handler.
*/
static std::function<unsigned int(
const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
const CellStatus status)>
make_weighting_callback(
const DoFHandler<dim, spacedim> &dof_handler,
const std::vector<unsigned int> &precomputed_weights);

private:
/**
Expand All @@ -213,7 +274,7 @@ namespace parallel
boost::signals2::connection connection;

/**
* A callback function that will be connected to the `weight` signal of
* A callback function that can be connected to the `weight` signal of
* the @p triangulation, to which the @p dof_handler is attached. Ultimately
* returns the weight for each cell, determined by the @p weighting_function
* provided as a parameter. Returns zero if @p dof_handler has not been
Expand All @@ -226,6 +287,21 @@ namespace parallel
const DoFHandler<dim, spacedim> &dof_handler,
const parallel::TriangulationBase<dim, spacedim> &triangulation,
const WeightingFunction &weighting_function);

/**
* A callback function that can be connected to the `weight` signal of
* the @p triangulation, to which the @p dof_handler is attached. Ultimately
* returns the weight for each cell, determined by the @p precomputed_weights
* provided as a parameter. Returns zero if @p dof_handler has not been
* initialized yet.
*/
static unsigned int
weighting_callback(
const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
const CellStatus status,
const DoFHandler<dim, spacedim> &dof_handler,
const parallel::TriangulationBase<dim, spacedim> &triangulation,
const std::vector<unsigned int> &precomputed_weights);
};
} // namespace parallel

Expand Down