Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
repartitioning_policy_tools.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2021 - 2022 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
19
22
25
27
28
30{
31 namespace
32 {
33 template <int dim, int spacedim>
34 void
35 add_indices_recursively_for_first_child_policy(
37 const internal::CellIDTranslator<dim> & cell_id_translator,
38 IndexSet & is_fine)
39 {
40 is_fine.add_index(cell_id_translator.translate(cell));
41
42 if (cell->level() > 0 && cell->parent()->child(0) == cell)
43 add_indices_recursively_for_first_child_policy(cell->parent(),
44 cell_id_translator,
45 is_fine);
46 }
47 } // namespace
48
49
50 template <int dim, int spacedim>
52 : tighten(tighten)
53 {}
54
55 template <int dim, int spacedim>
58 const Triangulation<dim, spacedim> &tria_in) const
59 {
60 if (tighten == false)
61 return {}; // nothing to do
62
63 const auto tria =
65 &tria_in);
66
67 if (tria == nullptr)
68 return {}; // nothing to do, since serial triangulation
69
70#ifndef DEAL_II_WITH_MPI
71 Assert(false, ExcNeedsMPI());
72 return {};
73#else
74
75 const auto comm = tria->get_communicator();
76
77 const unsigned int process_has_active_locally_owned_cells =
78 tria->n_locally_owned_active_cells() > 0;
79 const unsigned int n_processes_with_active_locally_owned_cells =
80 Utilities::MPI::sum(process_has_active_locally_owned_cells, comm);
81
82 if (n_processes_with_active_locally_owned_cells ==
84 return {}; // nothing to do, since all processes have cells
85
86 unsigned int offset = 0;
87
88 const int ierr = MPI_Exscan(&process_has_active_locally_owned_cells,
89 &offset,
90 1,
92 process_has_active_locally_owned_cells)>,
93 MPI_SUM,
94 comm);
95 AssertThrowMPI(ierr);
96
99
100 partition = offset;
101
102 return partition;
103#endif
104 }
105
106
107
108 template <int dim, int spacedim>
110 const Triangulation<dim, spacedim> &tria_fine)
111 : n_coarse_cells(tria_fine.n_global_coarse_cells())
112 , n_global_levels(tria_fine.n_global_levels())
113 {
114 Assert(
117 "FirstChildPolicy is only working for pure hex meshes at the moment."))
118
120 cell_id_translator(n_coarse_cells, n_global_levels);
121 is_level_partitions.set_size(cell_id_translator.size());
122
123 for (const auto &cell : tria_fine.active_cell_iterators())
124 if (cell->is_locally_owned())
125 add_indices_recursively_for_first_child_policy(cell,
126 cell_id_translator,
128 }
129
130
131
132 template <int dim, int spacedim>
135 const Triangulation<dim, spacedim> &tria_coarse_in) const
136 {
137 const auto communicator = tria_coarse_in.get_communicator();
138
139 const internal::CellIDTranslator<dim> cell_id_translator(n_coarse_cells,
140 n_global_levels);
141
142 IndexSet is_coarse(cell_id_translator.size());
143
144 for (const auto &cell : tria_coarse_in.active_cell_iterators())
145 if (cell->is_locally_owned())
146 is_coarse.add_index(cell_id_translator.translate(cell));
147
148 std::vector<unsigned int> owning_ranks_of_coarse_cells(
149 is_coarse.n_elements());
150 {
152 process(is_level_partitions,
153 is_coarse,
154 communicator,
155 owning_ranks_of_coarse_cells,
156 false);
157
159 std::vector<
160 std::pair<types::global_cell_index, types::global_cell_index>>,
161 std::vector<unsigned int>>
162 consensus_algorithm;
163 consensus_algorithm.run(process, communicator);
164 }
165
166 const auto tria =
167 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
168 &tria_coarse_in);
169
171
174
175 for (const auto &cell : tria_coarse_in.active_cell_iterators() |
177 partition[cell->global_active_cell_index()] =
178 owning_ranks_of_coarse_cells[is_coarse.index_within_set(
179 cell_id_translator.translate(cell))];
180
181 return partition;
182 }
183
184
185 template <int dim, int spacedim>
187 const unsigned int n_min_cells)
188 : n_min_cells(n_min_cells)
189 {}
190
191
192
193 template <int dim, int spacedim>
196 const Triangulation<dim, spacedim> &tria_in) const
197 {
198 const auto tria =
199 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
200 &tria_in);
201
203
206
207 // step 1) check if all processes have enough cells
208
209 const unsigned int n_locally_owned_active_cells =
210 std::count_if(tria_in.begin_active(),
212 tria_in.end()),
213 [](const auto &cell) { return cell.is_locally_owned(); });
214
215 const auto comm = tria_in.get_communicator();
216
217 if (Utilities::MPI::min(n_locally_owned_active_cells, comm) >= n_min_cells)
218 return {}; // all processes have enough cells
219
220 // step 2) there are processes which do not have enough cells so that
221 // a repartitioning kicks in with the aim that all processes that own
222 // cells have at least the specified number of cells
223
224 const types::global_cell_index n_global_active_cells =
225 tria_in.n_global_active_cells();
226
227 const unsigned int n_partitions =
228 std::max<unsigned int>(1,
229 std::min<types::global_cell_index>(
230 n_global_active_cells / n_min_cells,
232
233 const unsigned int min_cells = n_global_active_cells / n_partitions;
234
235 const auto convert = [&](const unsigned int i) {
236 // determine the owner of a given index: we try to assign each process
237 // the same number of cells; if there is a remainder, we assign the
238 // first processes an additional cell (so that the difference of number of
239 // locally owned cells is never larger than one between processes).
240 const unsigned int n_partitions_with_additional_cell =
241 n_global_active_cells - min_cells * n_partitions;
242
243 const unsigned int rank =
244 (i < (min_cells + 1) * n_partitions_with_additional_cell) ?
245 (i / (min_cells + 1)) :
246 ((i - n_partitions_with_additional_cell) / min_cells);
247
248 AssertIndexRange(rank, n_partitions);
249
250 return rank;
251 };
252
253 for (const auto i : partition.locally_owned_elements())
254 partition[i] = convert(i);
255
256 return partition;
257 }
258
259
260
261 template <int dim, int spacedim>
263 const std::function<
264 unsigned int(const typename Triangulation<dim, spacedim>::cell_iterator &,
266 &weighting_function)
267 : weighting_function(weighting_function)
268 {}
269
270
271
272 template <int dim, int spacedim>
275 const Triangulation<dim, spacedim> &tria_in) const
276 {
277#ifndef DEAL_II_WITH_MPI
278 (void)tria_in;
279 return {};
280#else
281
282 const auto tria =
283 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
284 &tria_in);
285
287
288 const auto partitioner =
290
291 std::vector<unsigned int> weights(partitioner->locally_owned_size());
292
293 const auto mpi_communicator = tria_in.get_communicator();
294 const auto n_subdomains = Utilities::MPI::n_mpi_processes(mpi_communicator);
295
296 // determine weight of each cell
297 for (const auto &cell :
299 weights[partitioner->global_to_local(cell->global_active_cell_index())] =
300 weighting_function(
302
303 // determine weight of all the cells locally owned by this process
304 std::uint64_t process_local_weight = 0;
305 for (const auto &weight : weights)
306 process_local_weight += weight;
307
308 // determine partial sum of weights of this process
309 std::uint64_t process_local_weight_offset = 0;
310
311 int ierr = MPI_Exscan(
312 &process_local_weight,
313 &process_local_weight_offset,
314 1,
315 Utilities::MPI::mpi_type_id_for_type<decltype(process_local_weight)>,
316 MPI_SUM,
318 AssertThrowMPI(ierr);
319
320 // total weight of all processes
321 std::uint64_t total_weight =
322 process_local_weight_offset + process_local_weight;
323
324 ierr =
325 MPI_Bcast(&total_weight,
326 1,
327 Utilities::MPI::mpi_type_id_for_type<decltype(total_weight)>,
328 n_subdomains - 1,
329 mpi_communicator);
330 AssertThrowMPI(ierr);
331
332 // set up partition
333 LinearAlgebra::distributed::Vector<double> partition(partitioner);
334
335 for (std::uint64_t i = 0, weight = process_local_weight_offset;
336 i < partition.locally_owned_size();
337 weight += weights[i], ++i)
338 partition.local_element(i) =
339 static_cast<double>(weight * n_subdomains / total_weight);
340
341 return partition;
342#endif
343 }
344
345
346} // namespace RepartitioningPolicyTools
347
348
349
350/*-------------- Explicit Instantiations -------------------------------*/
351#include "repartitioning_policy_tools.inst"
352
size_type index_within_set(const size_type global_index) const
Definition index_set.h:1883
size_type n_elements() const
Definition index_set.h:1816
void set_size(const size_type size)
Definition index_set.h:1649
void add_index(const size_type index)
Definition index_set.h:1680
CellWeightPolicy(const std::function< unsigned int(const typename Triangulation< dim, spacedim >::cell_iterator &, const typename Triangulation< dim, spacedim >::CellStatus)> &weighting_function)
virtual LinearAlgebra::distributed::Vector< double > partition(const Triangulation< dim, spacedim > &tria_in) const override
virtual LinearAlgebra::distributed::Vector< double > partition(const Triangulation< dim, spacedim > &tria_coarse_in) const override
virtual LinearAlgebra::distributed::Vector< double > partition(const Triangulation< dim, spacedim > &tria_coarse_in) const override
FirstChildPolicy(const Triangulation< dim, spacedim > &tria_fine)
virtual LinearAlgebra::distributed::Vector< double > partition(const Triangulation< dim, spacedim > &tria_in) const override
virtual types::global_cell_index n_global_active_cells() const
virtual MPI_Comm get_communicator() const
bool all_reference_cells_are_hyper_cube() const
cell_iterator end() const
virtual const std::weak_ptr< const Utilities::MPI::Partitioner > global_active_cell_index_partitioner() const
active_cell_iterator begin_active(const unsigned int level=0) const
virtual std::vector< unsigned int > run(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm) override
types::global_cell_index translate(const TriaIterator< Accessor > &cell) const
types::global_cell_index size() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:150
T min(const T &t, const MPI_Comm mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
Definition mpi.h:1732
const MPI_Comm comm
const ::Triangulation< dim, spacedim > & tria