Reference documentation for deal.II version 9.2.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\}}\)
fully_distributed_tria.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2019 - 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 
18 #include <deal.II/base/mpi.h>
20 
22 
24 
26 
27 // Forward declarations
28 namespace GridGenerator
29 {
30  template <int dim, int spacedim>
31  void
33  const double left,
34  const double right,
35  const bool colorize);
36 } // namespace GridGenerator
37 
38 namespace parallel
39 {
40  namespace fullydistributed
41  {
42  template <int dim, int spacedim>
44  : parallel::DistributedTriangulationBase<dim, spacedim>(mpi_communicator)
46  , partitioner([](::Triangulation<dim, spacedim> &tria,
47  const unsigned int n_partitions) {
49  })
50  , currently_processing_create_triangulation_for_internal_usage(false)
51  , currently_processing_prepare_coarsening_and_refinement_for_internal_usage(
52  false)
53  {}
54 
55 
56 
57  template <int dim, int spacedim>
58  void
61  &construction_data)
62  {
63  // check if the communicator of this parallel triangulation has been used
64  // to construct the TriangulationDescription::Description
65  Assert(construction_data.comm == this->mpi_communicator,
66  ExcMessage("MPI communicators do not match!"));
67 
68  // store internally the settings
69  settings = construction_data.settings;
70 
71  // set the smoothing properties
72  if (settings &
74  this->set_mesh_smoothing(
75  static_cast<
76  typename ::Triangulation<dim, spacedim>::MeshSmoothing>(
79  else
80  this->set_mesh_smoothing(
81  static_cast<
82  typename ::Triangulation<dim, spacedim>::MeshSmoothing>(
84 
85  this->set_mesh_smoothing(construction_data.smoothing);
86 
87  // clear internal data structures
88  this->coarse_cell_id_to_coarse_cell_index_vector.clear();
89  this->coarse_cell_index_to_coarse_cell_id_vector.clear();
90 
91  // check if no locally relevant coarse-grid cells have been provided
92  if (construction_data.coarse_cell_vertices.empty())
93  {
94  // 1) create a dummy hypercube
95  currently_processing_create_triangulation_for_internal_usage = true;
96  GridGenerator::hyper_cube(*this, 0, 1, false);
97  currently_processing_create_triangulation_for_internal_usage = false;
98 
99  // 2) mark cell as artificial
100  auto cell = this->begin();
101  cell->set_subdomain_id(::numbers::artificial_subdomain_id);
102  cell->set_level_subdomain_id(
104 
105  // 3) set up dummy mapping between locally relevant coarse-grid cells
106  // and global cells
107  this->coarse_cell_id_to_coarse_cell_index_vector.emplace_back(
109  this->coarse_cell_index_to_coarse_cell_id_vector.emplace_back(
111  }
112  else
113  {
114  // 1) store `coarse-cell index to coarse-cell id`-mapping
115  this->coarse_cell_index_to_coarse_cell_id_vector =
116  construction_data.coarse_cell_index_to_coarse_cell_id;
117 
118  // 2) set up `coarse-cell id to coarse-cell index`-mapping
119  std::map<types::coarse_cell_id, unsigned int>
120  coarse_cell_id_to_coarse_cell_index_vector;
121  for (unsigned int i = 0;
122  i < construction_data.coarse_cell_index_to_coarse_cell_id.size();
123  ++i)
124  coarse_cell_id_to_coarse_cell_index_vector
125  [construction_data.coarse_cell_index_to_coarse_cell_id[i]] = i;
126 
127  for (auto i : coarse_cell_id_to_coarse_cell_index_vector)
128  this->coarse_cell_id_to_coarse_cell_index_vector.emplace_back(i);
129 
130  // create locally-relevant
131  currently_processing_prepare_coarsening_and_refinement_for_internal_usage =
132  true;
133  currently_processing_create_triangulation_for_internal_usage = true;
135  construction_data);
136  currently_processing_prepare_coarsening_and_refinement_for_internal_usage =
137  false;
138  currently_processing_create_triangulation_for_internal_usage = false;
139 
140  // create a copy of cell_infos such that we can sort them
141  auto cell_infos = construction_data.cell_infos;
142 
143  // sort cell_infos on each level separately (as done in
144  // ::Triangulation::create_triangulation())
145  for (auto &cell_info : cell_infos)
146  std::sort(cell_info.begin(),
147  cell_info.end(),
150  const CellId a_id(a.id);
151  const CellId b_id(b.id);
152 
153  const auto a_coarse_cell_index =
154  this->coarse_cell_id_to_coarse_cell_index(
155  a_id.get_coarse_cell_id());
156  const auto b_coarse_cell_index =
157  this->coarse_cell_id_to_coarse_cell_index(
158  b_id.get_coarse_cell_id());
159 
160  // according to their coarse-cell index and if that is
161  // same according to their cell id (the result is that
162  // cells on each level are sorted according to their
163  // index on that level - what we need in the following
164  // operations)
165  if (a_coarse_cell_index != b_coarse_cell_index)
166  return a_coarse_cell_index < b_coarse_cell_index;
167  else
168  return a_id < b_id;
169  });
170 
171  // 4a) set all cells artificial (and set the actual
172  // (level_)subdomain_ids in the next step)
173  for (auto cell = this->begin(); cell != this->end(); ++cell)
174  {
175  if (cell->is_active())
176  cell->set_subdomain_id(
178 
179  cell->set_level_subdomain_id(
181  }
182 
183  // 4b) set actual (level_)subdomain_ids
184  for (unsigned int level = 0; level < cell_infos.size(); ++level)
185  {
186  auto cell = this->begin(level);
187  auto cell_info = cell_infos[level].begin();
188  for (; cell_info != cell_infos[level].end(); ++cell_info)
189  {
190  // find cell that has the correct cell
191  while (cell_info->id != cell->id().template to_binary<dim>())
192  ++cell;
193 
194  // subdomain id
195  if (cell->is_active())
196  cell->set_subdomain_id(cell_info->subdomain_id);
197 
198  // level subdomain id
199  if (settings & TriangulationDescription::Settings::
201  cell->set_level_subdomain_id(cell_info->level_subdomain_id);
202  }
203  }
204  }
205 
206  this->update_number_cache();
207  }
208 
209 
210 
211  template <int dim, int spacedim>
212  void
214  const std::vector<Point<spacedim>> & vertices,
215  const std::vector<::CellData<dim>> &cells,
216  const SubCellData & subcelldata)
217  {
218  AssertThrow(
219  currently_processing_create_triangulation_for_internal_usage,
220  ExcMessage(
221  "Use the other create_triangulation() function to create triangulations of type parallel::fullydistributed::Triangulation.!"));
222 
224  cells,
225  subcelldata);
226  }
227 
228 
229 
230  template <int dim, int spacedim>
231  void
233  const ::Triangulation<dim, spacedim> &other_tria)
234  {
235  // pointer to the triangulation for which the construction data
236  // should be created (normally it is the input triangulation but
237  // in the case of a serial triangulation we create a copy which should
238  // be used)
239  const ::Triangulation<dim, spacedim> *other_tria_ptr = &other_tria;
240 
241  // temporary serial triangulation (since the input triangulation is const
242  // and we might modify its subdomain_ids
243  // and level_subdomain_ids during partitioning); this pointer only points
244  // to anything if the source triangulation is serial, and ensures that our
245  // copy is eventually deleted
246  std::unique_ptr<::Triangulation<dim, spacedim>> serial_tria;
247 
248  // check if other triangulation is not a parallel one, which needs to be
249  // partitioned
250  if (dynamic_cast<const ::parallel::TriangulationBase<dim, spacedim>
251  *>(&other_tria) == nullptr)
252  {
253  serial_tria =
254  std_cxx14::make_unique<::Triangulation<dim, spacedim>>();
255 
256  // actually copy the serial triangulation
257  serial_tria->copy_triangulation(other_tria);
258 
259  // partition triangulation
260  this->partitioner(*serial_tria,
262  this->mpi_communicator));
263 
264  // partition multigrid levels
265  if (this->is_multilevel_hierarchy_constructed())
267 
268  // use the new serial triangulation to create the construction data
269  other_tria_ptr = serial_tria.get();
270  }
271 
272  // create construction data
273  const auto construction_data = TriangulationDescription::Utilities::
275  this->mpi_communicator,
276  this->settings);
277 
278  // finally create triangulation
279  this->create_triangulation(construction_data);
280  }
281 
282 
283 
284  template <int dim, int spacedim>
285  void
287  const std::function<void(::Triangulation<dim, spacedim> &,
288  const unsigned int)> &partitioner,
289  const TriangulationDescription::Settings & settings)
290  {
291  this->partitioner = partitioner;
292  this->settings = settings;
293  }
294 
295 
296 
297  template <int dim, int spacedim>
298  void
300  {
301  Assert(false, ExcNotImplemented());
302  }
303 
304 
305 
306  template <int dim, int spacedim>
307  bool
309  {
310  Assert(
311  currently_processing_prepare_coarsening_and_refinement_for_internal_usage,
312  ExcMessage("No coarsening and refinement is supported!"));
313 
314  return ::Triangulation<dim, spacedim>::
315  prepare_coarsening_and_refinement();
316  }
317 
318 
319 
320  template <int dim, int spacedim>
321  bool
323  {
324  Assert(false, ExcNotImplemented());
325  return false;
326  }
327 
328 
329 
330  template <int dim, int spacedim>
331  std::size_t
333  {
334  const std::size_t mem =
338  coarse_cell_id_to_coarse_cell_index_vector) +
340  coarse_cell_index_to_coarse_cell_id_vector);
341  return mem;
342  }
343 
344 
345 
346  template <int dim, int spacedim>
347  bool
349  {
350  return (
351  settings &
353  }
354 
355 
356 
357  template <int dim, int spacedim>
358  unsigned int
361  {
362  const auto coarse_cell_index = std::lower_bound(
363  coarse_cell_id_to_coarse_cell_index_vector.begin(),
364  coarse_cell_id_to_coarse_cell_index_vector.end(),
366  [](const std::pair<types::coarse_cell_id, unsigned int> &pair,
367  const types::coarse_cell_id &val) { return pair.first < val; });
368  Assert(coarse_cell_index !=
369  coarse_cell_id_to_coarse_cell_index_vector.cend(),
370  ExcMessage("Coarse cell index not found!"));
371  return coarse_cell_index->second;
372  }
373 
374 
375 
376  template <int dim, int spacedim>
379  const unsigned int coarse_cell_index) const
380  {
381  AssertIndexRange(coarse_cell_index,
382  coarse_cell_index_to_coarse_cell_id_vector.size());
383 
384  const auto coarse_cell_id =
385  coarse_cell_index_to_coarse_cell_id_vector[coarse_cell_index];
387  ExcMessage("You are trying to access a dummy cell!"));
388  return coarse_cell_id;
389  }
390 
391 
392 
393  } // namespace fullydistributed
394 } // namespace parallel
395 
396 
397 
398 /*-------------- Explicit Instantiations -------------------------------*/
399 #include "fully_distributed_tria.inst"
400 
401 
TriangulationDescription::default_setting
@ default_setting
Definition: tria_description.h:258
TriangulationDescription::Description::cell_infos
std::vector< std::vector< CellData< dim > > > cell_infos
Definition: tria_description.h:384
TriangulationDescription::Description::comm
MPI_Comm comm
Definition: tria_description.h:396
GridTools::partition_triangulation_zorder
void partition_triangulation_zorder(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const bool group_siblings=true)
Definition: grid_tools.cc:2816
parallel::DistributedTriangulationBase
Definition: tria_base.h:352
parallel::fullydistributed::Triangulation::set_partitioner
void set_partitioner(const std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> &partitioner, const TriangulationDescription::Settings &settings)
Definition: fully_distributed_tria.cc:286
parallel::fullydistributed::Triangulation
Definition: fully_distributed_tria.h:115
CellData
Definition: tria_description.h:67
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
parallel::fullydistributed::Triangulation::execute_coarsening_and_refinement
virtual void execute_coarsening_and_refinement() override
Definition: fully_distributed_tria.cc:299
Triangulation
Definition: tria.h:1109
memory_consumption.h
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
MPI_Comm
parallel::TriangulationBase::memory_consumption
virtual std::size_t memory_consumption() const override
Definition: tria_base.cc:88
parallel::fullydistributed::Triangulation::Triangulation
Triangulation(MPI_Comm mpi_communicator)
Definition: fully_distributed_tria.cc:43
GridGenerator::Airfoil::create_triangulation
void create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
TriangulationDescription::Description::coarse_cell_index_to_coarse_cell_id
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id
Definition: tria_description.h:377
level
unsigned int level
Definition: grid_out.cc:4355
parallel::fullydistributed::Triangulation::memory_consumption
virtual std::size_t memory_consumption() const override
Definition: fully_distributed_tria.cc:332
Triangulation::create_triangulation
virtual void create_triangulation(const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
Definition: tria.cc:10521
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
mpi.h
TriangulationDescription::Description::coarse_cell_vertices
std::vector< Point< spacedim > > coarse_cell_vertices
Definition: tria_description.h:370
grid_tools.h
TriangulationDescription::Utilities::create_description_from_triangulation
Description< dim, spacedim > create_description_from_triangulation(const ::Triangulation< dim, spacedim > &tria, const MPI_Comm comm, const TriangulationDescription::Settings settings=TriangulationDescription::Settings::default_setting, const unsigned int my_rank_in=numbers::invalid_unsigned_int)
Definition: tria_description.cc:104
Utilities::lower_bound
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1102
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
GridTools::partition_multigrid_levels
void partition_multigrid_levels(Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2921
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
parallel::fullydistributed::Triangulation::prepare_coarsening_and_refinement
virtual bool prepare_coarsening_and_refinement() override
Definition: fully_distributed_tria.cc:308
TriangulationDescription::CellData
Definition: tria_description.h:282
TriangulationDescription::Description::smoothing
Triangulation< dim, spacedim >::MeshSmoothing smoothing
Definition: tria_description.h:406
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
TriangulationDescription
Definition: tria_description.h:245
parallel::fullydistributed::Triangulation::create_triangulation
void create_triangulation(const TriangulationDescription::Description< dim, spacedim > &construction_data) override
Definition: fully_distributed_tria.cc:59
GridGenerator
Definition: grid_generator.h:57
parallel::fullydistributed::Triangulation::copy_triangulation
void copy_triangulation(const ::Triangulation< dim, spacedim > &other_tria) override
Definition: fully_distributed_tria.cc:232
unsigned int
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
types::coarse_cell_id
global_cell_index coarse_cell_id
Definition: types.h:114
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
TriangulationDescription::Description::settings
Settings settings
Definition: tria_description.h:401
GridGenerator::hyper_cube
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
parallel::fullydistributed::Triangulation::coarse_cell_index_to_coarse_cell_id
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const override
Definition: fully_distributed_tria.cc:378
Point< spacedim >
parallel::fullydistributed::Triangulation::has_hanging_nodes
virtual bool has_hanging_nodes() const override
Definition: fully_distributed_tria.cc:322
parallel::fullydistributed::Triangulation::is_multilevel_hierarchy_constructed
virtual bool is_multilevel_hierarchy_constructed() const override
Definition: fully_distributed_tria.cc:348
fully_distributed_tria.h
parallel::fullydistributed::Triangulation::coarse_cell_id_to_coarse_cell_index
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const override
Definition: fully_distributed_tria.cc:359
colorize
bool colorize
Definition: grid_out.cc:4354
numbers::artificial_subdomain_id
const types::subdomain_id artificial_subdomain_id
Definition: types.h:302
memory.h
TriangulationDescription::Description
Definition: tria_description.h:347
SubCellData
Definition: tria_description.h:199
TriangulationDescription::Settings
Settings
Definition: tria_description.h:253
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
numbers::invalid_coarse_cell_id
const types::coarse_cell_id invalid_coarse_cell_id
Definition: types.h:215
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
parallel
Definition: distributed.h:416
TriangulationDescription::construct_multigrid_hierarchy
@ construct_multigrid_hierarchy
Definition: tria_description.h:264
int