ExaDG
Loading...
Searching...
No Matches
grid_utilities.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_GRID_UTILITIES_H_
23#define INCLUDE_EXADG_GRID_GRID_UTILITIES_H_
24
25// deal.II
26#include <deal.II/fe/fe_simplex_p.h>
27#include <deal.II/fe/mapping_fe.h>
28#include <deal.II/fe/mapping_q.h>
29#include <deal.II/grid/grid_in.h>
30#include <deal.II/grid/grid_tools.h>
31#include <deal.II/multigrid/mg_transfer_global_coarsening.h>
32
33// ExaDG
34#include <exadg/grid/balanced_granularity_partition_policy.h>
35#include <exadg/grid/grid.h>
36#include <exadg/grid/grid_data.h>
37#include <exadg/grid/perform_local_refinements.h>
38
39namespace ExaDG
40{
41namespace GridUtilities
42{
43template<int dim>
44using PeriodicFacePairs = std::vector<
45 dealii::GridTools::PeriodicFacePair<typename dealii::Triangulation<dim>::cell_iterator>>;
46
50template<int dim>
51void
52create_mapping(std::shared_ptr<dealii::Mapping<dim>> & mapping,
53 ElementType const & element_type,
54 unsigned int const & mapping_degree)
55{
56 if(element_type == ElementType::Hypercube)
57 {
58 mapping = std::make_shared<dealii::MappingQ<dim>>(mapping_degree);
59 }
60 else if(element_type == ElementType::Simplex)
61 {
62 mapping = std::make_shared<dealii::MappingFE<dim>>(dealii::FE_SimplexP<dim>(mapping_degree));
63 }
64 else
65 {
66 AssertThrow(false, dealii::ExcMessage("Invalid parameter element_type."));
67 }
68}
69
76template<int dim, typename Number>
77void
78create_mapping_with_multigrid(std::shared_ptr<dealii::Mapping<dim>> & mapping,
79 std::shared_ptr<MultigridMappings<dim, Number>> & multigrid_mappings,
80 ElementType const & element_type,
81 unsigned int const & mapping_degree_fine,
82 unsigned int const & mapping_degree_coarse,
83 bool const involves_h_multigrid)
84{
85 // create fine mapping
86 create_mapping(mapping, element_type, mapping_degree_fine);
87
88 // create coarse mappings if needed
89 std::shared_ptr<dealii::Mapping<dim>> coarse_mapping;
90 if(involves_h_multigrid and (mapping_degree_coarse != mapping_degree_fine))
91 {
92 create_mapping(coarse_mapping, element_type, mapping_degree_coarse);
93 }
94 multigrid_mappings = std::make_shared<MultigridMappings<dim, Number>>(mapping, coarse_mapping);
95}
96
102template<int dim>
103std::vector<dealii::GridTools::PeriodicFacePair<typename dealii::DoFHandler<dim>::cell_iterator>>
104transform_periodic_face_pairs_to_dof_cell_iterator(
105 std::vector<dealii::GridTools::PeriodicFacePair<
106 typename dealii::Triangulation<dim>::cell_iterator>> const & periodic_faces,
107 dealii::DoFHandler<dim> const & dof_handler)
108{
109 typedef typename std::vector<
110 dealii::GridTools::PeriodicFacePair<typename dealii::DoFHandler<dim>::cell_iterator>>
111 PeriodicFacesDoF;
112
113 PeriodicFacesDoF periodic_faces_dof;
114
115 for(auto it : periodic_faces)
116 {
117 dealii::GridTools::PeriodicFacePair<typename dealii::DoFHandler<dim>::cell_iterator>
118 face_pair_dof_hander;
119
120 face_pair_dof_hander.cell[0] = it.cell[0]->as_dof_handler_iterator(dof_handler);
121 face_pair_dof_hander.cell[1] = it.cell[1]->as_dof_handler_iterator(dof_handler);
122
123 face_pair_dof_hander.face_idx[0] = it.face_idx[0];
124 face_pair_dof_hander.face_idx[1] = it.face_idx[1];
125
126 face_pair_dof_hander.orientation = it.orientation;
127 face_pair_dof_hander.matrix = it.matrix;
128
129 periodic_faces_dof.push_back(face_pair_dof_hander);
130 }
131
132 return periodic_faces_dof;
133}
134
143template<int dim>
144inline void
145create_triangulation(
146 std::shared_ptr<dealii::Triangulation<dim>> & triangulation,
147 PeriodicFacePairs<dim> & periodic_face_pairs,
148 MPI_Comm const & mpi_comm,
149 GridData const & data,
150 bool const construct_multigrid_hierarchy,
151 std::function<void(dealii::Triangulation<dim> &,
152 PeriodicFacePairs<dim> &,
153 unsigned int const,
154 std::vector<unsigned int> const &)> const & lambda_create_triangulation,
155 unsigned int const global_refinements,
156 std::vector<unsigned int> const & vector_local_refinements)
157{
158 if(vector_local_refinements.size() != 0)
159 {
160 AssertThrow(
161 data.element_type == ElementType::Hypercube,
162 dealii::ExcMessage(
163 "Local refinements are currently only supported for meshes composed of hypercube elements."));
164 }
165
166 // mesh smoothing
167 auto mesh_smoothing = dealii::Triangulation<dim>::none;
168
169 // the mesh used for a simulation should not depend on the preconditioner used to solve the
170 // problem (like whether multigrid is used or not). Hence, we always set this parameter for
171 // ElementType::Hypercube. It is not posible to set this parameter for ElementType::Simplex due to
172 // the implementation in deal.II.
173 if(data.element_type == ElementType::Hypercube)
174 {
175 mesh_smoothing = dealii::Triangulation<dim>::limit_level_difference_at_vertices;
176 }
177
178 if(data.triangulation_type == TriangulationType::Serial)
179 {
180 AssertDimension(dealii::Utilities::MPI::n_mpi_processes(mpi_comm), 1);
181 triangulation = std::make_shared<dealii::Triangulation<dim>>(mesh_smoothing);
182
183 lambda_create_triangulation(*triangulation,
184 periodic_face_pairs,
185 global_refinements,
186 vector_local_refinements);
187 }
188 else if(data.triangulation_type == TriangulationType::Distributed)
189 {
190 typename dealii::parallel::distributed::Triangulation<dim>::Settings distributed_settings;
191
192 if(construct_multigrid_hierarchy)
193 {
194 distributed_settings =
195 dealii::parallel::distributed::Triangulation<dim>::construct_multigrid_hierarchy;
196 }
197
198 triangulation =
199 std::make_shared<dealii::parallel::distributed::Triangulation<dim>>(mpi_comm,
200 mesh_smoothing,
201 distributed_settings);
202
203 lambda_create_triangulation(*triangulation,
204 periodic_face_pairs,
205 global_refinements,
206 vector_local_refinements);
207 }
208 else if(data.triangulation_type == TriangulationType::FullyDistributed)
209 {
210 auto const serial_grid_generator = [&](dealii::Triangulation<dim, dim> & tria_serial) {
211 lambda_create_triangulation(tria_serial,
212 periodic_face_pairs,
213 global_refinements,
214 vector_local_refinements);
215 };
216
217 auto const serial_grid_partitioner = [&](dealii::Triangulation<dim, dim> & tria_serial,
218 MPI_Comm const comm,
219 unsigned int const group_size) {
220 (void)group_size;
221 if(data.partitioning_type == PartitioningType::Metis)
222 {
223 dealii::GridTools::partition_triangulation(dealii::Utilities::MPI::n_mpi_processes(comm),
224 tria_serial);
225 }
226 else if(data.partitioning_type == PartitioningType::z_order)
227 {
228 dealii::GridTools::partition_triangulation_zorder(
229 dealii::Utilities::MPI::n_mpi_processes(comm), tria_serial);
230 }
231 else
232 {
233 AssertThrow(false, dealii::ExcNotImplemented());
234 }
235 };
236
237 unsigned int const group_size = 1;
238
239 typename dealii::TriangulationDescription::Settings triangulation_description_setting =
240 dealii::TriangulationDescription::default_setting;
241
242 if(construct_multigrid_hierarchy)
243 {
244 triangulation_description_setting =
245 dealii::TriangulationDescription::construct_multigrid_hierarchy;
246 }
247
248 triangulation =
249 std::make_shared<dealii::parallel::fullydistributed::Triangulation<dim>>(mpi_comm);
250
251 auto const description = dealii::TriangulationDescription::Utilities::
252 create_description_from_triangulation_in_groups<dim, dim>(serial_grid_generator,
253 serial_grid_partitioner,
254 triangulation->get_communicator(),
255 group_size,
256 mesh_smoothing,
257 triangulation_description_setting);
258
259 triangulation->create_triangulation(description);
260 }
261 else
262 {
263 AssertThrow(false, dealii::ExcMessage("Invalid parameter triangulation_type."));
264 }
265}
266
278template<int dim>
279inline void
280create_coarse_triangulations_automatically_from_fine_triangulation(
281 dealii::Triangulation<dim> const & fine_triangulation,
282 PeriodicFacePairs<dim> const & fine_periodic_face_pairs,
283 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> & coarse_triangulations_const,
284 std::vector<PeriodicFacePairs<dim>> & coarse_periodic_face_pairs,
285 GridData const & data)
286{
287 // In case of a serial or distributed triangulation, deal.II can automatically generate the
288 // coarse triangulations.
289 AssertThrow(data.triangulation_type == TriangulationType::Serial or
290 data.triangulation_type == TriangulationType::Distributed,
291 dealii::ExcMessage("Invalid parameter triangulation_type."));
292
293 if(data.triangulation_type == TriangulationType::Serial)
294 {
295 AssertThrow(
296 fine_triangulation.all_reference_cells_are_hyper_cube(),
297 dealii::ExcMessage(
298 "The create_geometric_coarsening_sequence function of dealii does currently not support "
299 "simplicial elements."));
300
301 coarse_triangulations_const =
302 dealii::MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence(
303 fine_triangulation);
304 }
305 else if(data.triangulation_type == TriangulationType::Distributed)
306 {
307 AssertThrow(
308 fine_triangulation.all_reference_cells_are_hyper_cube(),
309 dealii::ExcMessage(
310 "dealii::parallel::distributed::Triangulation does not support simplicial elements."));
311
312 coarse_triangulations_const =
313 dealii::MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence(
314 fine_triangulation,
315 BalancedGranularityPartitionPolicy<dim>(
316 dealii::Utilities::MPI::n_mpi_processes(fine_triangulation.get_communicator())));
317 }
318 else
319 {
320 AssertThrow(false, dealii::ExcMessage("Invalid parameter triangulation_type."));
321 }
322
323 // deal.II adds the fine triangulation as the last entry of the vector. According to our
324 // convention in ExaDG, we only include the triangulations coarser than the fine level. Hence,
325 // we remove the last entry of the vector.
326 coarse_triangulations_const.pop_back();
327
328 coarse_periodic_face_pairs.resize(coarse_triangulations_const.size());
329 for(unsigned int level = 0; level < coarse_periodic_face_pairs.size(); ++level)
330 {
331 coarse_periodic_face_pairs[level] = fine_periodic_face_pairs;
332 }
333}
334
335template<int dim>
336inline void
337create_coarse_triangulations_for_fully_distributed_triangulation(
338 dealii::Triangulation<dim> const & fine_triangulation,
339 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> & coarse_triangulations_const,
340 std::vector<PeriodicFacePairs<dim>> & coarse_periodic_face_pairs,
341 GridData const & data,
342 std::function<void(dealii::Triangulation<dim> &,
343 PeriodicFacePairs<dim> &,
344 unsigned int const,
345 std::vector<unsigned int> const &)> const & lambda_create_triangulation,
346 std::vector<unsigned int> const vector_local_refinements)
347{
348 // In case of a fully distributed triangulation, deal.II cannot automatically generate the
349 // coarse triangulations. Create the coarse traingulations using the lambda function and
350 // the vector of local refinements.
351 AssertThrow(data.triangulation_type == TriangulationType::FullyDistributed,
352 dealii::ExcMessage("Invalid parameter triangulation_type."));
353
354 if(fine_triangulation.n_global_levels() >= 2)
355 {
356 // Resize the empty coarse_triangulations and coarse_periodic_face_pairs vectors.
357 std::vector<std::shared_ptr<dealii::Triangulation<dim>>> coarse_triangulations =
358 std::vector<std::shared_ptr<dealii::Triangulation<dim>>>(
359 fine_triangulation.n_global_levels() - 1);
360
361 coarse_periodic_face_pairs = std::vector<PeriodicFacePairs<dim>>(coarse_triangulations.size());
362
363 // Start one level below the fine triangulation.
364 unsigned int level = fine_triangulation.n_global_levels() - 2;
365 std::vector<unsigned int> refine_local = vector_local_refinements;
366
367 // Undo global refinements.
368 if(data.n_refine_global >= 1)
369 {
370 unsigned int const n_refine_global_start = (unsigned int)(data.n_refine_global - 1);
371 for(int refine_global = n_refine_global_start; refine_global >= 0; --refine_global)
372 {
373 GridUtilities::create_triangulation<dim>(coarse_triangulations[level],
374 coarse_periodic_face_pairs[level],
375 fine_triangulation.get_communicator(),
376 data,
377 false /*construct_multigrid_hierarchy */,
378 lambda_create_triangulation,
379 refine_global,
381
382 if(level > 0)
383 {
384 level--;
385 }
386 }
387 }
388
389 // Undo local refinements.
390 if(refine_local.size() > 0)
391 {
392 while(*std::max_element(refine_local.begin(), refine_local.end()) != 0)
393 {
394 for(size_t material_id = 0; material_id < refine_local.size(); material_id++)
395 {
396 if(refine_local[material_id] > 0)
397 {
398 refine_local[material_id]--;
399 }
400 }
401
402 GridUtilities::create_triangulation<dim>(coarse_triangulations[level],
403 coarse_periodic_face_pairs[level],
404 fine_triangulation.get_communicator(),
405 data,
406 false /*construct_multigrid_hierarchy */,
407 lambda_create_triangulation,
408 0 /*refine_global*/,
410
411 if(level > 0)
412 {
413 level--;
414 }
415 }
416 }
417
418 AssertThrow(
419 level == 0,
420 dealii::ExcMessage(
421 "There occurred a logical error when creating the geometric coarsening sequence."));
422
423 // Make all entries in the vector of shared pointers const.
424 for(auto const & it : coarse_triangulations)
425 {
426 coarse_triangulations_const.push_back(it);
427 }
428 }
429}
430
431template<int dim>
432inline void
433create_coarse_triangulations(
434 dealii::Triangulation<dim> const & fine_triangulation,
435 PeriodicFacePairs<dim> const & fine_periodic_face_pairs,
436 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> & coarse_triangulations_const,
437 std::vector<PeriodicFacePairs<dim>> & coarse_periodic_face_pairs,
438 GridData const & data,
439 std::function<void(dealii::Triangulation<dim> &,
440 PeriodicFacePairs<dim> &,
441 unsigned int const,
442 std::vector<unsigned int> const &)> const & lambda_create_triangulation,
443 std::vector<unsigned int> const vector_local_refinements)
444{
445 // In case of a serial or distributed triangulation, deal.II can automatically generate the
446 // coarse triangulations, otherwise, the coarse triangulations have to be explicitily created
447 // using the provided lambda function and the vector of local refinements.
448 if(data.triangulation_type == TriangulationType::Serial or
449 data.triangulation_type == TriangulationType::Distributed)
450 {
451 create_coarse_triangulations_automatically_from_fine_triangulation(fine_triangulation,
452 fine_periodic_face_pairs,
453 coarse_triangulations_const,
454 coarse_periodic_face_pairs,
455 data);
456 }
457 else if(data.triangulation_type == TriangulationType::FullyDistributed)
458 {
459 create_coarse_triangulations_for_fully_distributed_triangulation(fine_triangulation,
460 coarse_triangulations_const,
461 coarse_periodic_face_pairs,
462 data,
463 lambda_create_triangulation,
464 vector_local_refinements);
465 }
466 else
467 {
468 AssertThrow(false, dealii::ExcMessage("Invalid parameter triangulation_type."));
469 }
470}
471
478template<int dim>
479inline void
480create_triangulation(
481 Grid<dim> & grid,
482 MPI_Comm const & mpi_comm,
483 GridData const & data,
484 std::function<void(dealii::Triangulation<dim> &,
485 PeriodicFacePairs<dim> &,
486 unsigned int const,
487 std::vector<unsigned int> const &)> const & lambda_create_triangulation,
488 std::vector<unsigned int> const vector_local_refinements)
489{
490 GridUtilities::create_triangulation(grid.triangulation,
491 grid.periodic_face_pairs,
492 mpi_comm,
493 data,
494 false /*construct_multigrid_hierarchy */,
495 lambda_create_triangulation,
496 data.n_refine_global,
497 vector_local_refinements);
498}
499
505template<int dim>
506inline void
507create_triangulation_with_multigrid(
508 Grid<dim> & grid,
509 MPI_Comm const & mpi_comm,
510 GridData const & data,
511 bool const involves_h_multigrid,
512 std::function<void(dealii::Triangulation<dim> &,
513 PeriodicFacePairs<dim> &,
514 unsigned int const,
515 std::vector<unsigned int> const &)> const & lambda_create_triangulation,
516 std::vector<unsigned int> const vector_local_refinements)
517{
518 if(involves_h_multigrid)
519 {
520 // Make sure that we create coarse triangulations in case of simplex meshes.
521 if(data.element_type == ElementType::Simplex)
522 {
523 AssertThrow(data.create_coarse_triangulations == true,
524 dealii::ExcMessage(
525 "You need to set GridData::create_coarse_triangulations = true "
526 "in order to use h-multigrid for simplex meshes."));
527 }
528
529 // create fine triangulation
530 GridUtilities::create_triangulation(grid.triangulation,
531 grid.periodic_face_pairs,
532 mpi_comm,
533 data,
534 not data.create_coarse_triangulations,
535 lambda_create_triangulation,
536 data.n_refine_global,
537 vector_local_refinements);
538
539 // Make sure that we create coarse triangulations in case of meshes with hanging nodes.
540 if(grid.triangulation->has_hanging_nodes())
541 {
542 AssertThrow(data.create_coarse_triangulations == true,
543 dealii::ExcMessage(
544 "You need to set GridData::create_coarse_triangulations = true "
545 "in order to use h-multigrid for meshes with hanging nodes."));
546 }
547
548 // create coarse triangulations
549 if(data.create_coarse_triangulations)
550 {
551 GridUtilities::create_coarse_triangulations(*grid.triangulation,
552 grid.periodic_face_pairs,
553 grid.coarse_triangulations,
554 grid.coarse_periodic_face_pairs,
555 data,
556 lambda_create_triangulation,
557 vector_local_refinements);
558 }
559 }
560 else
561 {
562 // If no h-multigrid is involved, simply re-direct to the other function that creates the fine
563 // triangulation only.
564 GridUtilities::create_triangulation<dim>(
565 grid, mpi_comm, data, lambda_create_triangulation, vector_local_refinements);
566 }
567}
568
573template<int dim>
574inline void
575create_coarse_triangulations_after_coarsening_and_refinement(
576 dealii::Triangulation<dim> const & fine_triangulation,
577 PeriodicFacePairs<dim> const & fine_periodic_face_pairs,
578 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> & coarse_triangulations_const,
579 std::vector<PeriodicFacePairs<dim>> & coarse_periodic_face_pairs,
580 GridData const & data,
581 bool const amr_preserves_boundary_cells)
582{
583 if(data.triangulation_type == TriangulationType::Serial or
584 data.triangulation_type == TriangulationType::Distributed)
585 {
586 // Update periodic face pairs if existent.
587 if(fine_periodic_face_pairs.size() > 0)
588 {
589 AssertThrow(amr_preserves_boundary_cells,
590 dealii::ExcMessage(
591 "Combination of adaptive mesh refinement and periodic face pairs"
592 " requires boundary cells to be preserved."));
593 }
594
595 create_coarse_triangulations_automatically_from_fine_triangulation(fine_triangulation,
596 fine_periodic_face_pairs,
597 coarse_triangulations_const,
598 coarse_periodic_face_pairs,
599 data);
600 }
601 else if(data.triangulation_type == TriangulationType::FullyDistributed)
602 {
603 AssertThrow(false,
604 dealii::ExcMessage("Combination of adaptive mesh refinement and "
605 "TriangulationType::FullyDistributed not implemented."));
606 }
607 else
608 {
609 AssertThrow(false, dealii::ExcMessage("Invalid parameter triangulation_type."));
610 }
611}
612
617template<int dim>
618inline void
619read_external_triangulation(dealii::Triangulation<dim, dim> & tria, GridData const & data)
620{
621 AssertThrow(not data.file_name.empty(),
622 dealii::ExcMessage(
623 "You are trying to read a grid file, but the string, which is supposed to contain"
624 " the file name, is empty."));
625
626 dealii::GridIn<dim> grid_in;
627
628 grid_in.attach_triangulation(tria);
629
630 // find the file extension from the file name
631 std::string extension = data.file_name.substr(data.file_name.find_last_of('.') + 1);
632
633 AssertThrow(not extension.empty(),
634 dealii::ExcMessage("You are trying to read a grid file, but the file extension is"
635 " empty."));
636
637 // decide the file format
638 typename dealii::GridIn<dim>::Format format;
639 if(extension == "e" || extension == "exo")
640 format = dealii::GridIn<dim>::Format::exodusii;
641 else
642 format = grid_in.parse_format(extension);
643
644 // TODO: check if the exodusIIData is needed
645 // typename dealii::GridIn<dim>::ExodusIIData exodusIIData;
646
647 grid_in.read(data.file_name, format);
648
649 AssertThrow(get_element_type(tria) == data.element_type,
650 dealii::ExcMessage("You are trying to read a grid file, but the element type of the"
651 " external grid file and the element type specified in GridData"
652 " don't match. Most likely, you forgot to change the element_type"
653 " parameter of GridData to the desired element type."));
654}
655
656} // namespace GridUtilities
657} // namespace ExaDG
658
659
660#endif /* INCLUDE_EXADG_GRID_GRID_UTILITIES_H_ */
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
void refine_local(dealii::Triangulation< dim > &tria, std::vector< unsigned int > const &vector_local_refinements)
Definition perform_local_refinements.h:43