ExaDG
Loading...
Searching...
No Matches
balanced_granularity_partition_policy.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_GRID_BALANCED_GRANULARITY_PARTITION_POLICY_H_
23#define INCLUDE_EXADG_GRID_BALANCED_GRANULARITY_PARTITION_POLICY_H_
24
25// deal.II
26#include <deal.II/distributed/fully_distributed_tria.h>
27
28namespace ExaDG
29{
35template<int dim, int spacedim = dim>
37 : public dealii::RepartitioningPolicyTools::Base<dim, spacedim>
38{
39public:
40 BalancedGranularityPartitionPolicy(unsigned int const n_mpi_processes)
41 : n_mpi_processes_per_level{n_mpi_processes}
42 {
43 }
44
46
47 virtual dealii::LinearAlgebra::distributed::Vector<double>
48 partition(dealii::Triangulation<dim, spacedim> const & tria_coarse_in) const override
49 {
50 dealii::types::global_cell_index const n_cells = tria_coarse_in.n_global_active_cells();
51
52 // TODO: We hard-code a grain-size limit of 200 cells per processor
53 // (assuming linear finite elements and typical behavior of
54 // supercomputers). In case we have fewer cells on the fine level, we do
55 // not immediately go to 200 cells per rank, but limit the growth by a
56 // factor of 8, which limits makes sure that we do not create too many
57 // messages for individual MPI processes.
58 unsigned int const grain_size_limit =
59 std::min<unsigned int>(200, 8 * n_cells / n_mpi_processes_per_level.back() + 1);
60
61 dealii::RepartitioningPolicyTools::MinimalGranularityPolicy<dim, spacedim> partitioning_policy(
62 grain_size_limit);
63 dealii::LinearAlgebra::distributed::Vector<double> const partitions =
64 partitioning_policy.partition(tria_coarse_in);
65
66 // The vector 'partitions' contains the partition numbers. To get the
67 // number of partitions, we take the infinity norm.
68 n_mpi_processes_per_level.push_back(static_cast<unsigned int>(partitions.linfty_norm()) + 1);
69 return partitions;
70 }
71
72private:
73 mutable std::vector<unsigned int> n_mpi_processes_per_level;
74};
75} // namespace ExaDG
76
77
78#endif /* INCLUDE_EXADG_GRID_BALANCED_GRANULARITY_PARTITION_POLICY_H_ */
Definition balanced_granularity_partition_policy.h:38
Definition driver.cpp:33