ExaDG
Loading...
Searching...
No Matches
grid_data.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_DATA_H_
23#define INCLUDE_EXADG_GRID_GRID_DATA_H_
24
25// C/C++
26#include <string>
27
28// deal.II
29#include <deal.II/grid/tria.h>
30
31// ExaDG
32#include <exadg/utilities/print_functions.h>
33
34namespace ExaDG
35{
36/*
37 * Triangulation type
38 */
39enum class TriangulationType
40{
41 Serial,
42 Distributed,
43 FullyDistributed
44};
45
46/*
47 * Element type
48 */
49enum class ElementType
50{
51 Hypercube,
52 Simplex
53};
54
55/*
56 * Partitioning type (relevant for fully-distributed triangulation)
57 */
58enum class PartitioningType
59{
60 Metis,
61 z_order
62};
63
68template<int dim>
69ElementType
70get_element_type(dealii::Triangulation<dim> const & tria)
71{
72 if(tria.all_reference_cells_are_simplex())
73 {
74 return ElementType::Simplex;
75 }
76 else if(tria.all_reference_cells_are_hyper_cube())
77 {
78 return ElementType::Hypercube;
79 }
80 else
81 {
82 AssertThrow(false, dealii::ExcMessage("Invalid parameter element_type."));
83 return ElementType::Hypercube;
84 }
85}
86
88{
89 GridData()
90 : triangulation_type(TriangulationType::Distributed),
91 element_type(ElementType::Hypercube),
92 partitioning_type(PartitioningType::Metis),
93 n_refine_global(0),
94 file_name(),
95 create_coarse_triangulations(false)
96 {
97 }
98
99 void
100 check() const
101 {
102 }
103
104 void
105 print(dealii::ConditionalOStream const & pcout) const
106 {
107 print_parameter(pcout, "Triangulation type", triangulation_type);
108
109 print_parameter(pcout, "Element type", element_type);
110
111 if(triangulation_type == TriangulationType::FullyDistributed)
112 print_parameter(pcout, "Partitioning type (fully-distributed)", partitioning_type);
113
114 print_parameter(pcout, "Number of global refinements", n_refine_global);
115
116 if(not file_name.empty())
117 print_parameter(pcout, "Grid file name", file_name);
118
119 print_parameter(pcout, "Create coarse triangulations", create_coarse_triangulations);
120 }
121
122 TriangulationType triangulation_type;
123
124 ElementType element_type;
125
126 // only relevant for TriangulationType::FullyDistributed
127 PartitioningType partitioning_type;
128
129 unsigned int n_refine_global;
130
131 // path to a grid file
132 // the filename needs to include a proper filename ending/extension so that we can internally
133 // deduce the correct type of the file format
134 std::string file_name;
135
136 // In case of a hypercube mesh that is globally refined, i.e. without hanging nodes, the fine
137 // triangulation can be used for all multigrid h-levels without the need to create coarse
138 // triangulations explicitly. Hence, this parameter is typically set to false for globally-refined
139 // hypercube meshes.
140 // Nevertheless, it is possible to set this parameter to true for globally-refined hypercube
141 // meshes. In that case, the coarse triangulations are created explicitly for use in
142 // h-multigrid methods.
143 // This parameter needs to be set to true if one wants to use h-multigrid methods for
144 // locally-refined hypercube meshes or non-hypercube meshes.
145 bool create_coarse_triangulations;
146};
147
148} // namespace ExaDG
149
150
151
152#endif /* INCLUDE_EXADG_GRID_GRID_DATA_H_ */
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
Definition grid_data.h:88