Changed: For weighted load balancing with parallel::distributed::Triangulation
objects, an intial weight of `1000` will no longer be assigned to each cell.
<br>
-You can invoke the old behavior by connecting a function to the signal
-that returns the base weight like this:
+The old signal Triangulation::Signals::cell_weight has been deprecated.
+Please use the new signal Triangulation::Signals::weight instead.
+<br>
+You can invoke the old behavior by connecting a function to the new
+signal that returns the base weight like this:
@code{.cc}
-triangulation.signals.cell_weight.connect(
+triangulation.signals.weight.connect(
[](const typename parallel::distributed::Triangulation<dim>::cell_iterator &,
const typename parallel::distributed::Triangulation<dim>::CellStatus)
-> unsigned int { return 1000; });
// return value of this function (representing "work for this cell") is
// calculated based on the number of particles in the current cell.
// The function is
- // connected to the cell_weight() signal inside the triangulation, and will be
+ // connected to the `weight` signal inside the triangulation, and will be
// called once per cell, whenever the triangulation repartitions the domain
// between ranks (the connection is created inside the
// generate_particles() function of this class).
// be created once, so we might as well have set it up in the constructor
// of this class, but for the purpose of this example we want to group the
// particle related instructions.
- background_triangulation.signals.cell_weight.connect(
+ background_triangulation.signals.weight.connect(
[&](
const typename parallel::distributed::Triangulation<dim>::cell_iterator
&cell,
* this case, an analogous code example looks as follows.
* @code
* boost::signals2::connection connection =
- * hp_dof_handler.get_triangulation().signals.cell_weight.connect(
+ * 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}));
*
* The use of this class is demonstrated in step-75.
*
- * @note See Triangulation::Signals::cell_weight for more information on
+ * @note See Triangulation::Signals::weight for more information on
* weighting and load balancing.
*
* @note Be aware that this class connects the weight function to the
*
* @note A hp::FECollection needs to be attached to your DoFHandler object via
* DoFHandler::distribute_dofs() <em>once before</em> the
- * Triangulation::Signals::cell_weight signal will be triggered. Otherwise,
+ * Triangulation::Signals::weight signal will be triggered. Otherwise,
* your DoFHandler does not know many degrees of freedom your cells have. In
* other words, you need to call DoFHandler::distribute_dofs() once before you
* call
private:
/**
- * A connection to the corresponding cell_weight signal of the Triangulation
+ * A connection to the corresponding weight signal of the Triangulation
* which is attached to the DoFHandler.
*/
boost::signals2::connection connection;
/**
- * A callback function that will be connected to the cell_weight signal of
+ * A callback function that will 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
* @note This function by default partitions the mesh in such a way that
* the number of cells on all processors is roughly equal. If you want
* to set weights for partitioning, e.g. because some cells are more
- * expensive to compute than others, you can use the signal cell_weight
+ * expensive to compute than others, you can use the signal weight
* as documented in the dealii::Triangulation class. This function will
* check whether a function is connected to the signal and if so use it.
* If you prefer to repartition the mesh yourself at user-defined
* flag to the constructor, which ensures that calling the current
* function only refines and coarsens the triangulation, but doesn't
* partition it. You can then call the repartition() function manually.
- * The usage of the cell_weights signal is identical in both cases, if a
+ * The usage of the weight signal is identical in both cases, if a
* function is connected to the signal it will be used to balance the
* calculated weights, otherwise the number of cells is balanced.
*/
* the same way as execute_coarsening_and_refinement() with respect to
* dealing with data movement (SolutionTransfer, etc.).
*
- * @note If no function is connected to the cell_weight signal described
+ * @note If no function is connected to the weight signal described
* in the dealii::Triangulation class, this function will balance the
* number of cells on each processor. If one or more functions are
* connected, it will calculate the sum of the weights and balance the
* single-partition case without packages installed, and only requires them
* installed when multiple partitions are required.
*
- * @note If the @p cell_weight signal has been attached to the @p triangulation,
+ * @note If the @p weight signal has been attached to the @p triangulation,
* then this will be used and passed to the partitioner.
*/
template <int dim, int spacedim>
* case like this, partitioning algorithm may sometimes make bad decisions and
* you may want to build your own connectivity graph.
*
- * @note If the @p cell_weight signal has been attached to the @p triangulation,
+ * @note If the @p weight signal has been attached to the @p triangulation,
* then this will be used and passed to the partitioner.
*/
template <int dim, int spacedim>
};
/**
- * A structure used to accumulate the results of the cell_weights slot
+ * A structure used to accumulate the results of the weight signal slot
* functions below. It takes an iterator range and returns the sum of
* values.
*/
boost::signals2::signal<unsigned int(const cell_iterator &,
const CellStatus),
CellWeightSum<unsigned int>>
- cell_weight;
+ weight;
+
+#ifndef DOXYGEN
+ /**
+ * Legacy constructor that connects deprecated signals to their new ones.
+ */
+ Signals()
+ : cell_weight(weight)
+ {}
+
+ /**
+ * Legacy signal emulation to deprecate the old signal.
+ */
+ class LegacySignal
+ {
+ public:
+ using signature_type = unsigned int(const cell_iterator &,
+ const CellStatus);
+ using combiner_type = CellWeightSum<unsigned int>;
+
+ using slot_function_type = boost::function<signature_type>;
+ using slot_type =
+ boost::signals2::slot<signature_type, slot_function_type>;
+
+ LegacySignal(
+ boost::signals2::signal<signature_type, combiner_type> &new_signal)
+ : new_signal(new_signal)
+ {}
+
+ ~LegacySignal()
+ {
+ base_weight.disconnect();
+ }
+
+ DEAL_II_DEPRECATED_EARLY
+ boost::signals2::connection
+ connect(
+ const slot_type & slot,
+ boost::signals2::connect_position position = boost::signals2::at_back)
+ {
+ if (base_weight.connected() == false)
+ {
+ base_weight = new_signal.connect(
+ [](const cell_iterator &, const CellStatus) -> unsigned int {
+ return 1000;
+ });
+ Assert(base_weight.connected() && new_signal.num_slots() == 1,
+ ExcInternalError());
+ }
+
+ return new_signal.connect(slot, position);
+ }
+
+ DEAL_II_DEPRECATED_EARLY
+ std::size_t
+ num_slots() const
+ {
+ return new_signal.num_slots() -
+ static_cast<std::size_t>(base_weight.connected());
+ }
+
+ DEAL_II_DEPRECATED_EARLY
+ bool
+ empty() const
+ {
+ if (num_slots() == 0)
+ {
+ Assert(new_signal.num_slots() == 0, ExcInternalError());
+ return true;
+ }
+ return false;
+ }
+
+ template <typename S>
+ DEAL_II_DEPRECATED_EARLY void
+ disconnect(const S &connection)
+ {
+ new_signal.disconnect(connection);
+
+ if (num_slots() == 0)
+ {
+ Assert(base_weight.connected() && new_signal.num_slots() == 1,
+ ExcInternalError());
+ new_signal.disconnect(base_weight);
+ }
+ }
+
+ DEAL_II_DEPRECATED_EARLY
+ unsigned int
+ operator()(const cell_iterator &iterator, const CellStatus status)
+ {
+ return new_signal(iterator, status);
+ }
+
+ private:
+ boost::signals2::connection base_weight;
+
+ boost::signals2::signal<signature_type, combiner_type> &new_signal;
+ };
+#endif
+
+ /**
+ * @copydoc weight
+ *
+ * As a reference a value of 1000 is added for every cell to the total
+ * weight. This means a signal return value of 1000 (resulting in a weight
+ * of 2000) means that it is twice as expensive for a process to handle this
+ * particular cell.
+ *
+ * @deprecated Use the weight signal instead which omits the base weight.
+ * You can invoke the old behavior by connecting a function to the signal
+ * that returns the base weight like this:
+ * @code{.cc}
+ * triangulation.signals.weight.connect(
+ * [](const typename Triangulation<dim>::cell_iterator &,
+ * const typename Triangulation<dim>::CellStatus)
+ * -> unsigned int { return 1000; });
+ * @endcode
+ */
+ LegacySignal cell_weight;
/**
* This signal is triggered at the beginning of execution of the
&weighting_function)
{
connection.disconnect();
- connection = dof_handler.get_triangulation().signals.cell_weight.connect(
+ connection = dof_handler.get_triangulation().signals.weight.connect(
make_weighting_callback(dof_handler, weighting_function));
}
{
public:
/**
- * This constructor assumes the cell_weights are already sorted in the
+ * This constructor assumes the @p cell_weights are already sorted in the
* order that p4est will encounter the cells, and they do not contain
* ghost cells or artificial cells.
*/
{
// partition the new mesh between all processors. If cell weights
// have not been given balance the number of cells.
- if (this->signals.cell_weight.num_slots() == 0)
+ if (this->signals.weight.empty())
dealii::internal::p4est::functions<dim>::partition(
parallel_forest,
/* prepare coarsening */ 1,
(parallel_forest->mpisize + 1));
}
- if (this->signals.cell_weight.num_slots() == 0)
+ if (this->signals.weight.empty())
{
// no cell weights given -- call p4est's 'partition' without a
// callback for cell weights
// Iterate over p4est and Triangulation relations
// to find refined/coarsened/kept
- // cells. Then append cell_weight.
+ // cells. Then append weight.
// Note that we need to follow the p4est ordering
- // instead of the deal.II ordering to get the cell_weights
+ // instead of the deal.II ordering to get the weights
// in the same order p4est will encounter them during repartitioning.
for (const auto &cell_rel : this->local_cell_relations)
{
const auto &cell_it = cell_rel.first;
const auto &cell_status = cell_rel.second;
- weights.push_back(this->signals.cell_weight(cell_it, cell_status));
+ weights.push_back(this->signals.weight(cell_it, cell_status));
}
return weights;
std::vector<unsigned int> cell_weights;
// Get cell weighting if a signal has been attached to the triangulation
- if (!triangulation.signals.cell_weight.empty())
+ if (!triangulation.signals.weight.empty())
{
cell_weights.resize(triangulation.n_active_cells(), 0U);
for (const auto &cell : triangulation.active_cell_iterators())
if (cell->is_locally_owned())
cell_weights[cell->active_cell_index()] =
- triangulation.signals.cell_weight(
+ triangulation.signals.weight(
cell, Triangulation<dim, spacedim>::CellStatus::CELL_PERSIST);
// If this is a parallel triangulation, we then need to also
std::vector<unsigned int> cell_weights;
// Get cell weighting if a signal has been attached to the triangulation
- if (!triangulation.signals.cell_weight.empty())
+ if (!triangulation.signals.weight.empty())
{
cell_weights.resize(triangulation.n_active_cells(), 0U);
for (const auto &cell : triangulation.active_cell_iterators() |
IteratorFilters::LocallyOwnedCell())
cell_weights[cell->active_cell_index()] =
- triangulation.signals.cell_weight(
+ triangulation.signals.weight(
cell, Triangulation<dim, spacedim>::CellStatus::CELL_PERSIST);
// If this is a parallel triangulation, we then need to also
"are already partitioned implicitly and can not be "
"partitioned again explicitly."));
Assert(n_partitions > 0, ExcInvalidNumberOfPartitions(n_partitions));
- Assert(triangulation.signals.cell_weight.empty(), ExcNotImplemented());
+ Assert(triangulation.signals.weight.empty(), ExcNotImplemented());
// signal that partitioning is going to happen
triangulation.signals.pre_partition();
// repartition the mesh as described above, first in some arbitrary
// way, and then with all equal weights
- tr.signals.cell_weight.connect(std::bind(&cell_weight_1<dim>,
- std::placeholders::_1,
- std::placeholders::_2));
+ tr.signals.weight.connect(&cell_weight_1<dim>);
tr.repartition();
- tr.signals.cell_weight.disconnect_all_slots();
+ tr.signals.weight.disconnect_all_slots();
- tr.signals.cell_weight.connect(std::bind(&cell_weight_2<dim>,
- std::placeholders::_1,
- std::placeholders::_2));
+ tr.signals.weight.connect(&cell_weight_2<dim>);
tr.repartition();
const auto n_locally_owned_active_cells_per_processor =
// repartition the mesh as described above, first in some arbitrary
// way, and then with no weights
- tr.signals.cell_weight.connect(std::bind(&cell_weight_1<dim>,
- std::placeholders::_1,
- std::placeholders::_2));
+ tr.signals.weight.connect(&cell_weight_1<dim>);
tr.repartition();
- tr.signals.cell_weight.disconnect_all_slots();
+ tr.signals.weight.disconnect_all_slots();
tr.repartition();
const auto n_locally_owned_active_cells_per_processor =
GridGenerator::subdivided_hyper_cube(tr, 16);
- tr.signals.cell_weight.connect(&cell_weight<dim>);
+ tr.signals.weight.connect(&cell_weight<dim>);
tr.refine_global(1);
const auto n_locally_owned_active_cells_per_processor =
tr.refine_global(1);
- tr.signals.cell_weight.connect(&cell_weight<dim>);
+ tr.signals.weight.connect(&cell_weight<dim>);
tr.repartition();
GridGenerator::subdivided_hyper_cube(tr, 16);
- tr.signals.cell_weight.connect(&cell_weight<dim>);
+ tr.signals.weight.connect(&cell_weight<dim>);
tr.refine_global(1);
// repartition the mesh; attach different weights to all cells
- tr.signals.cell_weight.connect(&cell_weight<dim>);
+ tr.signals.weight.connect(&cell_weight<dim>);
tr.repartition();
// repartition the mesh; attach different weights to all cells
n_global_active_cells = tr.n_global_active_cells();
- tr.signals.cell_weight.connect(&cell_weight<dim>);
+ tr.signals.weight.connect(&cell_weight<dim>);
tr.repartition();
const auto n_locally_owned_active_cells_per_processor =
// return value of this function (representing "work for this cell") is
// calculated based on the number of particles in the current cell.
// The function is
- // connected to the cell_weight() signal inside the triangulation, and will be
+ // connected to the `weight` signal inside the triangulation, and will be
// called once per cell, whenever the triangulation repartitions the domain
// between ranks (the connection is created inside the
// generate_particles() function of this class).
// be created once, so we might as well have set it up in the constructor
// of this class, but for the purpose of this example we want to group the
// particle related instructions.
- background_triangulation.signals.cell_weight.connect(
+ background_triangulation.signals.weight.connect(
[&](
const typename parallel::distributed::Triangulation<dim>::cell_iterator
&cell,