1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_distributed_cell_weights_h
17 #define dealii_distributed_cell_weights_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/distributed/tria_base.h>
22 
23 #include <deal.II/hp/dof_handler.h>
24 
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 namespace parallel
29 {
30   /**
31    * Anytime a parallel::TriangulationBase is repartitioned, either upon request
32    * or by refinement/coarsening, cells will be distributed amongst all
33    * subdomains to achieve an equally balanced workload. If the workload per
34    * cell varies, which is in general the case for hp::DoFHandler objects, we
35    * can take that into account by introducing individual weights for
36    * different cells.
37    *
38    * This class allows computing these weights for load balancing by
39    * consulting the FiniteElement that is associated with each cell of
40    * a hp::DoFHandler. One can choose from predefined weighting
41    * algorithms provided by this class or provide a custom one.
42    *
43    * This class offers two different ways of connecting the chosen weighting
44    * function to the corresponding signal of the linked
45    * parallel::TriangulationBase. The recommended way involves creating an
46    * object of this class which will automatically take care of registering the
47    * weighting function upon creation and de-registering it once destroyed. An
48    * object of this class needs to exist for every DoFHandler associated with
49    * the Triangulation we work on to achieve satisfying work balancing results.
50    * The connected weighting function may be changed anytime using the
51    * CellWeights::reinit() function. The following code snippet demonstrates how
52    * to achieve each cell being weighted by its current number of degrees of
53    * freedom. We chose a factor of `1000` that corresponds to the initial weight
54    * each cell is assigned to upon creation.
55    * @code
56    * parallel::CellWeights<dim, spacedim> cell_weights(
57    *   hp_dof_handler,
58    *   parallel::CellWeights<dim, spacedim>::ndofs_weighting({1000, 1}));
59    * @endcode
60    *
61    * On the other hand, you are also able to take care of handling the signal
62    * connection manually by using the static member function of this class. In
63    * this case, an analogous code example looks as follows.
64    * @code
65    * boost::signals2::connection connection =
66    *   hp_dof_handler.get_triangulation().signals.cell_weight.connect(
67    *     parallel::CellWeights<dim, spacedim>::make_weighting_callback(
68    *       hp_dof_handler,
69    *       parallel::CellWeights<dim, spacedim>::ndofs_weighting(
70    *         {1000, 1}));
71    * @endcode
72    *
73    * @note See Triangulation::Signals::cell_weight for more information on
74    * weighting and load balancing.
75    *
76    * @note Be aware that this class connects the weight function to the
77    * Triangulation during this class's constructor. If the Triangulation
78    * associated with the DoFHandler changes during the lifetime of the
79    * latter via hp::DoFHandler::initialize(), an assertion will be triggered in
80    * the weight_callback() function. Use CellWeights::reinit() to deregister the
81    * weighting function on the old Triangulation and connect it to the new one.
82    *
83    * @ingroup distributed
84    */
85   template <int dim, int spacedim = dim>
86   class CellWeights
87   {
88   public:
89     /**
90      * An alias that defines the characteristics of a function that can be used
91      * for weighting cells during load balancing.
92      *
93      * Such weighting functions take as arguments an iterator to a cell and the
94      * future finite element that will be assigned to it after repartitioning.
95      * They return an unsigned integer which is interpreted as the cell's
96      * weight or, in other words, the additional computational load associated
97      * with it.
98      */
99     using WeightingFunction = std::function<
100       unsigned int(const typename DoFHandler<dim, spacedim>::cell_iterator &,
101                    const FiniteElement<dim, spacedim> &)>;
102 
103     /**
104      * Constructor.
105      *
106      * @param[in] dof_handler The hp::DoFHandler which will be used to
107      *    determine each cell's finite element.
108      * @param[in] weighting_function The function that determines each
109      *    cell's weight during load balancing.
110      */
111     CellWeights(const dealii::DoFHandler<dim, spacedim> &dof_handler,
112                 const WeightingFunction &                weighting_function);
113 
114     /**
115      * Destructor.
116      *
117      * Disconnects the function previously connected to the weighting signal.
118      */
119     ~CellWeights();
120 
121     /**
122      * Connect a different @p weighting_function to the Triangulation
123      * associated with the @p dof_handler.
124      *
125      * Disconnects the function previously connected to the weighting signal.
126      */
127     void
128     reinit(const DoFHandler<dim, spacedim> &dof_handler,
129            const WeightingFunction &        weighting_function);
130 
131     /**
132      * Converts a @p weighting_function to a different type that qualifies as
133      * a callback function, which can be connected to a weighting signal of a
134      * Triangulation.
135      *
136      * This function does <b>not</b> connect the converted function to the
137      * Triangulation associated with the @p dof_handler.
138      */
139     static std::function<unsigned int(
140       const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
141       const typename dealii::Triangulation<dim, spacedim>::CellStatus status)>
142     make_weighting_callback(const DoFHandler<dim, spacedim> &dof_handler,
143                             const WeightingFunction &weighting_function);
144 
145     /**
146      * @name Selection of weighting functions
147      * @{
148      */
149 
150     /**
151      * Choose a constant weight @p factor on each cell.
152      */
153     static WeightingFunction
154     constant_weighting(const unsigned int factor = 1000);
155 
156     /**
157      * The pair of floating point numbers $(a,b)$ provided via
158      * @p coefficients determines the weight $w_K$ of each cell $K$ with
159      * $n_K$ degrees of freedom in the following way: \f[ w_K =
160      * a \, n_K^b \f]
161      *
162      * The right hand side will be rounded to the nearest integer since cell
163      * weights are required to be integers.
164      */
165     static WeightingFunction
166     ndofs_weighting(const std::pair<float, float> &coefficients);
167 
168     /**
169      * The container @p coefficients provides pairs of floating point numbers
170      * $(a_i, b_i)$ that determine the weight $w_K$ of each cell
171      * $K$ with $n_K$ degrees of freedom in the following way: \f[ w_K =
172      * \sum_i a_i \, n_K^{b_i} \f]
173      *
174      * The right hand side will be rounded to the nearest integer since cell
175      * weights are required to be integers.
176      */
177     static WeightingFunction
178     ndofs_weighting(const std::vector<std::pair<float, float>> &coefficients);
179 
180     /**
181      * @}
182      */
183 
184     /**
185      * @name Deprecated functions
186      * @{
187      */
188 
189     /**
190      * Constructor.
191      *
192      * @param[in] dof_handler The hp::DoFHandler which will be used to
193      *    determine each cell's finite element.
194      */
195     DEAL_II_DEPRECATED
196     CellWeights(const dealii::DoFHandler<dim, spacedim> &dof_handler);
197 
198     /**
199      * Choose a constant weight @p factor on each cell.
200      *
201      * @deprecated Use CellWeights::constant_weighting() instead.
202      */
203     DEAL_II_DEPRECATED void
204     register_constant_weighting(const unsigned int factor = 1000);
205 
206     /**
207      * Choose a weight for each cell that is proportional to its number of
208      * degrees of freedom with a factor @p factor.
209      *
210      * @deprecated Use CellWeights::ndofs_weighting() instead.
211      */
212     DEAL_II_DEPRECATED void
213     register_ndofs_weighting(const unsigned int factor = 1000);
214 
215     /**
216      * Choose a weight for each cell that is proportional to its number of
217      * degrees of freedom <i>squared</i> with a factor @p factor.
218      *
219      * @deprecated Use CellWeights::ndofs_weighting() instead.
220      */
221     DEAL_II_DEPRECATED void
222     register_ndofs_squared_weighting(const unsigned int factor = 1000);
223 
224     /**
225      * Register a custom weight for each cell by providing a function as a
226      * parameter.
227      *
228      * @param[in] custom_function The custom weighting function returning
229      *    the weight of each cell as an unsigned integer. It is required
230      *    to have two arguments, namely the FiniteElement that will be
231      *    active on the particular cell, and the cell itself of type
232      *    hp::DoFHandler::cell_iterator. We require both to make sure to
233      *    get the right active FiniteElement on each cell in case that we
234      *    coarsen the Triangulation.
235      */
236     DEAL_II_DEPRECATED void
237     register_custom_weighting(
238       const std::function<
239         unsigned int(const FiniteElement<dim, spacedim> &,
240                      const typename DoFHandler<dim, spacedim>::cell_iterator &)>
241         custom_function);
242 
243     /**
244      * @}
245      */
246 
247   private:
248     /**
249      * @name Deprecated members
250      * @{
251      */
252 
253     /**
254      * Pointer to the degree of freedom handler.
255      */
256     DEAL_II_DEPRECATED
257     SmartPointer<const dealii::DoFHandler<dim, spacedim>, CellWeights>
258       dof_handler;
259 
260     /**
261      * Pointer to the Triangulation associated with the degree of freedom
262      * handler.
263      *
264      * We store both to make sure to always work on the correct combination of
265      * both.
266      */
267     DEAL_II_DEPRECATED
268     SmartPointer<const parallel::TriangulationBase<dim, spacedim>, CellWeights>
269       triangulation;
270 
271     /**
272      * @}
273      */
274 
275     /**
276      * A connection to the corresponding cell_weight signal of the Triangulation
277      * which is attached to the DoFHandler.
278      */
279     boost::signals2::connection connection;
280 
281     /**
282      * A callback function that will be connected to the cell_weight signal of
283      * the @p triangulation, to which the @p dof_handler is attached. Ultimately
284      * returns the weight for each cell, determined by the @p weighting_function
285      * provided as a parameter.
286      */
287     static unsigned int
288     weighting_callback(
289       const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
290       const typename dealii::Triangulation<dim, spacedim>::CellStatus status,
291       const DoFHandler<dim, spacedim> &                 dof_handler,
292       const parallel::TriangulationBase<dim, spacedim> &triangulation,
293       const WeightingFunction &                         weighting_function);
294   };
295 } // namespace parallel
296 
297 
298 DEAL_II_NAMESPACE_CLOSE
299 
300 #endif
301