Skip to content

Commit

Permalink
Merge pull request #13561 from marcfehling/distributed-weight
Browse files Browse the repository at this point in the history
cell_weight: suggestions from code review.
  • Loading branch information
marcfehling committed Mar 22, 2022
2 parents b997598 + 78e07ee commit b0e7df8
Show file tree
Hide file tree
Showing 6 changed files with 76 additions and 50 deletions.
6 changes: 3 additions & 3 deletions include/deal.II/distributed/cell_weights.h
Original file line number Diff line number Diff line change
Expand Up @@ -84,7 +84,7 @@ namespace parallel
* weighting function on the old Triangulation and connect it to the new one.
*
* @note A hp::FECollection needs to be attached to your DoFHandler object via
* DoFHandler::distribute_dofs() <em>once before</em> the
* DoFHandler::distribute_dofs() <em>before</em> the
* 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
Expand Down Expand Up @@ -197,13 +197,13 @@ namespace parallel

private:
/**
* A connection to the corresponding 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 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
Expand Down
6 changes: 3 additions & 3 deletions include/deal.II/distributed/tria.h
Original file line number Diff line number Diff line change
Expand Up @@ -499,7 +499,7 @@ namespace parallel
* @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 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
Expand All @@ -508,7 +508,7 @@ namespace parallel
* 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 weight 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.
*/
Expand Down Expand Up @@ -542,7 +542,7 @@ namespace parallel
* 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 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
Expand Down
4 changes: 2 additions & 2 deletions include/deal.II/grid/grid_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -2166,7 +2166,7 @@ namespace GridTools
* single-partition case without packages installed, and only requires them
* installed when multiple partitions are required.
*
* @note If the @p weight signal has been attached to the @p triangulation,
* @note If the `weight` signal has been attached to the @p triangulation,
* then this will be used and passed to the partitioner.
*/
template <int dim, int spacedim>
Expand Down Expand Up @@ -2236,7 +2236,7 @@ namespace GridTools
* case like this, partitioning algorithm may sometimes make bad decisions and
* you may want to build your own connectivity graph.
*
* @note If the @p weight signal has been attached to the @p triangulation,
* @note If the `weight` signal has been attached to the @p triangulation,
* then this will be used and passed to the partitioner.
*/
template <int dim, int spacedim>
Expand Down
67 changes: 54 additions & 13 deletions include/deal.II/grid/tria.h
Original file line number Diff line number Diff line change
Expand Up @@ -2044,7 +2044,7 @@ class Triangulation : public Subscriptor
};

/**
* A structure used to accumulate the results of the weight signal 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.
*/
Expand Down Expand Up @@ -2173,8 +2173,7 @@ class Triangulation : public Subscriptor

/**
* This signal is triggered for each cell during every automatic or manual
* repartitioning. This signal will only be triggered if functions
* are connected to it. It is intended to allow a weighted repartitioning
* repartitioning. It is intended to allow a weighted repartitioning
* of the domain to balance the computational load across processes in a
* different way than balancing the number of cells. Any connected
* function is expected to take an iterator to a cell, and a CellStatus
Expand All @@ -2191,23 +2190,31 @@ class Triangulation : public Subscriptor
* refinement. If this cell is going to be coarsened, the signal is called
* for the parent cell and you need to provide the weight of the future
* parent cell. If this cell is going to be refined, the function is called
* on all children and you should ideally return the same weight for all
* children.
* on all children while `cell_iterator` refers to their parent cell. In
* this case, you need to pick a weight for each individual child based on
* information given by the parent cell.
*
* If several functions are connected to this signal, their return values
* will be summed to calculate the final weight via `CellWeightSum`.
* will be summed to calculate the final weight of a cell. This allows
* different parts of a larger code base to have their own functions
* computing the weight of a cell; for example in a code that does both
* finite element and particle computations on each cell, the code could
* separate the computation of a cell's weight into two functions, each
* implemented in their respective files, that provide the finite
* element-based and the particle-based weights.
*
* This function is used in step-68 and implicitely in step-75 using the
* This function is used in step-68 and implicitly in step-75 using the
* parallel::CellWeights class.
*/
boost::signals2::signal<unsigned int(const cell_iterator &,
const CellStatus),
CellWeightSum<unsigned int>>
weight;

#ifndef DOXYGEN
/**
* Legacy constructor that connects deprecated signals to their new ones.
* Constructor.
*
* Connects a deprecated signal to its successor.
*/
Signals()
: cell_weight(weight)
Expand All @@ -2227,16 +2234,28 @@ class Triangulation : public Subscriptor
using slot_type =
boost::signals2::slot<signature_type, slot_function_type>;

/**
* Constructor.
*/
LegacySignal(
boost::signals2::signal<signature_type, combiner_type> &new_signal)
: new_signal(new_signal)
{}

/**
* Destructor.
*/
~LegacySignal()
{
base_weight.disconnect();
}

/**
* Connects a function to the signal.
*
* Connects an additional base weight function if signal was previously
* empty.
*/
DEAL_II_DEPRECATED_EARLY
boost::signals2::connection
connect(
Expand All @@ -2256,6 +2275,10 @@ class Triangulation : public Subscriptor
return new_signal.connect(slot, position);
}

/**
* Returns the number of connected functions <em>without</em> the base
* weight.
*/
DEAL_II_DEPRECATED_EARLY
std::size_t
num_slots() const
Expand All @@ -2264,6 +2287,9 @@ class Triangulation : public Subscriptor
static_cast<std::size_t>(base_weight.connected());
}

/**
* Checks if there are any connected functions to the signal.
*/
DEAL_II_DEPRECATED_EARLY
bool
empty() const
Expand All @@ -2276,6 +2302,12 @@ class Triangulation : public Subscriptor
return false;
}

/**
* Disconnects a function from the signal.
*
* Also disconnects the base weight function if it is the last connected
* function.
*/
template <typename S>
DEAL_II_DEPRECATED_EARLY void
disconnect(const S &connection)
Expand All @@ -2290,6 +2322,9 @@ class Triangulation : public Subscriptor
}
}

/**
* Triggers the signal.
*/
DEAL_II_DEPRECATED_EARLY
unsigned int
operator()(const cell_iterator &iterator, const CellStatus status)
Expand All @@ -2298,23 +2333,29 @@ class Triangulation : public Subscriptor
}

private:
/**
* Monitors the connection of the base weight function.
*/
boost::signals2::connection base_weight;

/**
* Reference to the successor signal.
*/
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
* 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.
* @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:
* that returns the base weight as follows. This function should be added
* <em>in addition</em> to the one that actually computes the weight.
* @code{.cc}
* triangulation.signals.weight.connect(
* [](const typename Triangulation<dim>::cell_iterator &,
Expand Down
25 changes: 8 additions & 17 deletions source/distributed/tria.cc
Original file line number Diff line number Diff line change
Expand Up @@ -3361,19 +3361,14 @@ namespace parallel
// get cell weights for a weighted repartitioning.
const std::vector<unsigned int> cell_weights = get_cell_weights();

# ifdef DEBUG
// verify that the global sum of weights is larger than 0
const auto local_weights_sum =
std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0));
const auto global_weights_sum =
Utilities::MPI::sum(local_weights_sum, this->mpi_communicator);
Assert(global_weights_sum > 0,
Assert(Utilities::MPI::sum(std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0)),
this->mpi_communicator) > 0,
ExcMessage(
"The global sum of weights over all active cells "
"is zero. Please verify how you generate weights."));
# endif

PartitionWeights<dim, spacedim> partition_weights(cell_weights);

Expand Down Expand Up @@ -3543,18 +3538,14 @@ namespace parallel
// get cell weights for a weighted repartitioning.
const std::vector<unsigned int> cell_weights = get_cell_weights();

# ifdef DEBUG
// verify that the global sum of weights is larger than 0
const auto local_weights_sum = std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0));
const auto global_weights_sum =
Utilities::MPI::sum(local_weights_sum, this->mpi_communicator);
Assert(global_weights_sum > 0,
Assert(Utilities::MPI::sum(std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0)),
this->mpi_communicator) > 0,
ExcMessage(
"The global sum of weights over all active cells "
"is zero. Please verify how you generate weights."));
# endif

PartitionWeights<dim, spacedim> partition_weights(cell_weights);

Expand Down
18 changes: 6 additions & 12 deletions source/grid/grid_tools.cc
Original file line number Diff line number Diff line change
Expand Up @@ -3842,15 +3842,12 @@ namespace GridTools
shared_tria->get_communicator(),
cell_weights);

#ifdef DEBUG
// verify that the global sum of weights is larger than 0
const auto global_cell_sum = std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0));
Assert(global_cell_sum > 0,
Assert(std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0)) > 0,
ExcMessage("The global sum of weights over all active cells "
"is zero. Please verify how you generate weights."));
#endif
}

// Call the other more general function
Expand Down Expand Up @@ -3946,15 +3943,12 @@ namespace GridTools
shared_tria->get_communicator(),
cell_weights);

#ifdef DEBUG
// verify that the global sum of weights is larger than 0
const auto global_cell_sum = std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0));
Assert(global_cell_sum > 0,
Assert(std::accumulate(cell_weights.begin(),
cell_weights.end(),
std::uint64_t(0)) > 0,
ExcMessage("The global sum of weights over all active cells "
"is zero. Please verify how you generate weights."));
#endif
}

// Call the other more general function
Expand Down

0 comments on commit b0e7df8

Please sign in to comment.