ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
discretization.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include <boost/random/detail/polynomial.hpp>
9#include <compile_time_options.h>
10
11#include "discretization.h"
13
14#include <deal.II/base/quadrature_lib.h>
15#include <deal.II/fe/fe_dgq.h>
16#include <deal.II/fe/fe_nothing.h>
17#include <deal.II/fe/fe_q.h>
18#include <deal.II/fe/fe_simplex_p.h>
19#include <deal.II/fe/fe_tools.h>
20#include <deal.II/fe/mapping_fe.h>
21#include <deal.II/fe/mapping_q.h>
22#include <deal.II/grid/grid_out.h>
23
24#include <random>
25
26namespace ryujin
27{
28 using namespace dealii;
29
30 template <int dim>
32 const std::string &subsection)
33 : ParameterAcceptor(subsection)
34 , mpi_ensemble_(mpi_ensemble)
35 {
36 /* Options: */
37
38 ansatz_ = Ansatz::cg_q1;
39 add_parameter("finite element ansatz",
40 ansatz_,
41 "The finite element ansatz used for discretization. Valid "
42 "choices are cG Q1, cG Q2, cG Q3.");
43
44 mesh_type_ =
46 add_parameter("mesh type",
47 mesh_type_,
48 "The triangulation class used. Valid choices are \"serial\", "
49 "\"parallel shared\", \"parallel distributed\", \"parallel "
50 "fullydistributed\".");
51
52 if constexpr (dim == 1) {
53 geometry_ = "rectangular domain";
54 } else {
55 geometry_ = "cylinder";
56 }
57 add_parameter("geometry",
58 geometry_,
59 "Name of the geometry used to create the mesh. Valid names "
60 "are given by any of the subsections defined below.");
61
62 refinement_ = 5;
63 add_parameter("mesh refinement",
64 refinement_,
65 "number of refinement of global refinement steps");
66
67 mesh_writeout_ = true;
68 add_parameter("mesh writeout",
69 mesh_writeout_,
70 "Write out shared coarse mesh to a GMSH *.msh file.");
71
72 mesh_distortion_ = 0.;
73 add_parameter(
74 "mesh distortion", mesh_distortion_, "Strength of mesh distortion");
75
76 Geometries::populate_geometry_list<dim>(geometry_list_, subsection);
77 }
78
79
80 template <int dim>
81 void Discretization<dim>::prepare(const std::string &base_name)
82 {
83#ifdef DEBUG_OUTPUT
84 std::cout << "Discretization<dim>::prepare()" << std::endl;
85#endif
86
87 /* Select geometry: */
88
89 {
90 bool initialized = false;
91 for (auto &it : geometry_list_)
92 if (it->name() == geometry_) {
93 selected_geometry_ = it;
94 initialized = true;
95 break;
96 }
97
98 AssertThrow(
99 initialized,
100 ExcMessage("Could not find a geometry description with name \"" +
101 geometry_ + "\""));
102 }
103
104 /* Set up Triangulation object: */
105
106 const auto smoothing =
107 dealii::Triangulation<dim>::limit_level_difference_at_vertices;
108
109 switch (mesh_type_) {
111 triangulation_ = std::make_unique<
112 dealii::parallel::fullydistributed::Triangulation<dim>>(
113 mpi_ensemble_.ensemble_communicator());
114 triangulation_->set_mesh_smoothing(smoothing);
115 } break;
116
118 const auto settings = dealii::parallel::distributed::Triangulation<
119 dim>::Settings::construct_multigrid_hierarchy;
120 triangulation_ =
121 std::make_unique<dealii::parallel::distributed::Triangulation<dim>>(
122 mpi_ensemble_.ensemble_communicator(), smoothing, settings);
123 } break;
124
126 const auto settings = static_cast<
127 typename dealii::parallel::shared::Triangulation<dim>::Settings>(
128 dealii::parallel::shared::Triangulation<dim>::partition_auto |
129 dealii::parallel::shared::Triangulation<
130 dim>::construct_multigrid_hierarchy);
131 /* Beware of the boolean: */
132 triangulation_ =
133 std::make_unique<dealii::parallel::shared::Triangulation<dim>>(
134 mpi_ensemble_.ensemble_communicator(),
135 smoothing,
136 /*artificial cells*/ true,
137 settings);
138 } break;
139
140 case MeshType::serial: {
141 AssertThrow(
142 mpi_ensemble_.n_ensemble_ranks() == 1,
143 ExcMessage(
144 "The serial triangulation can only be used for serial "
145 "computations. If you want to run simulations with more than one "
146 "rank per ensemble, then please set \"mesh type\" to one of the "
147 "parallel triangulations supported by deal.II"));
148
149 triangulation_ = std::make_unique<dealii::Triangulation<dim>>(smoothing);
150
151 } break;
152
153 default:
154 __builtin_trap();
155 }
156
157 /* Create and distribute mesh: */
158
159 auto &triangulation = *triangulation_;
160 selected_geometry_->create_coarse_triangulation(triangulation);
161
162 if (mesh_writeout_ && dealii::Utilities::MPI::this_mpi_process(
163 mpi_ensemble_.ensemble_communicator()) == 0) {
164#ifdef DEAL_II_GMSH_WITH_API
165 GridOut grid_out;
166 grid_out.write_msh(triangulation, base_name + "-coarse_grid.msh");
167#else
168 GridOut grid_out;
169 GridOutFlags::Msh flags(/* write faces */ true, /* write lines */ true);
170 grid_out.set_flags(flags);
171 std::ofstream file(base_name + "-coarse_grid.msh");
172 grid_out.write_msh(triangulation, file);
173#endif
174 }
175
176 triangulation.refine_global(refinement_);
177
178 if (std::abs(mesh_distortion_) > 1.0e-10)
179 GridTools::distort_random(
180 mesh_distortion_, triangulation, false, std::random_device()());
181
182 const auto fe_degree = polynomial_degree();
183 const auto mapping_degree = fe_degree;
184 const auto quadrature_degree = fe_degree + 1;
185
186 /*
187 * First, let the selected geometry populate our hp::*Collection
188 * objects. If the method returns standard_quarilaterls, or
189 * standard_simplices, however, we need to do the setup ourselves:
190 */
191
192 const auto collection_type = selected_geometry_->populate_hp_collections(
193 fe_degree, have_discontinuous_ansatz(), collection_);
194
195 switch (collection_type) {
197 /*
198 * The geometry already populated the hp::*Collections
199 */
200
201 Assert(collection_.mapping, dealii::ExcInternalError());
202 Assert(collection_.finite_element, dealii::ExcInternalError());
203 Assert(collection_.finite_element_cg, dealii::ExcInternalError());
204 Assert(collection_.quadrature, dealii::ExcInternalError());
205 Assert(collection_.quadrature_high_order, dealii::ExcInternalError());
206 Assert(collection_.nodal_quadrature, dealii::ExcInternalError());
207 Assert(collection_.quadrature_1d, dealii::ExcInternalError());
208 Assert(collection_.face_quadrature, dealii::ExcInternalError());
209 Assert(collection_.face_nodal_quadrature, dealii::ExcInternalError());
210 } break;
211
213 /*
214 * Populate all collections with appropriate objects for the cG Qk, dG
215 * Qk finite element on purely quadrilateral, or hexahedral meshes:
216 */
217
218 if (have_discontinuous_ansatz()) {
219 collection_.finite_element =
220 std::make_unique<hp::FECollection<dim>>(FE_DGQ<dim>(fe_degree));
221 collection_.finite_element_cg =
222 std::make_unique<hp::FECollection<dim>>(FE_Q<dim>(fe_degree));
223 } else {
224 collection_.finite_element =
225 std::make_unique<hp::FECollection<dim>>(FE_Q<dim>(fe_degree));
226 collection_.finite_element_cg =
227 std::make_unique<hp::FECollection<dim>>(FE_Q<dim>(fe_degree));
228 }
229
230 collection_.mapping =
231 std::make_unique<dealii::hp::MappingCollection<dim>>(
232 MappingQ<dim>(mapping_degree));
233
234 collection_.quadrature = std::make_unique<hp::QCollection<dim>>(
235 QGauss<dim>(quadrature_degree));
236 collection_.quadrature_high_order =
237 std::make_unique<hp::QCollection<dim>>(
238 QGauss<dim>(quadrature_degree + 1));
239 collection_.nodal_quadrature = std::make_unique<hp::QCollection<dim>>(
240 QGaussLobatto<dim>(quadrature_degree));
241 collection_.quadrature_1d =
242 std::make_unique<hp::QCollection<1>>(QGauss<1>(quadrature_degree));
243 collection_.face_quadrature = std::make_unique<hp::QCollection<dim - 1>>(
244 QGauss<dim - 1>(quadrature_degree));
245 collection_.face_nodal_quadrature =
246 std::make_unique<hp::QCollection<dim - 1>>(
247 QGaussLobatto<dim - 1>(quadrature_degree));
248 } break;
249
251 /*
252 * Populate all collections with appropriate objects for the cG Pk, dG
253 * Pk finite element on purely quadrilateral, or hexahedral meshes:
254 */
255
256 if (have_discontinuous_ansatz()) {
257 collection_.finite_element = std::make_unique<hp::FECollection<dim>>(
258 FE_SimplexDGP<dim>(fe_degree));
259 collection_.finite_element_cg = std::make_unique<hp::FECollection<dim>>(
260 FE_SimplexP<dim>(fe_degree));
261 } else {
262 collection_.finite_element = std::make_unique<hp::FECollection<dim>>(
263 FE_SimplexP<dim>(fe_degree));
264 collection_.finite_element_cg = std::make_unique<hp::FECollection<dim>>(
265 FE_SimplexP<dim>(fe_degree));
266 }
267
268 collection_.mapping =
269 std::make_unique<dealii::hp::MappingCollection<dim>>(
270 MappingFE<dim>(FE_SimplexP<dim>(mapping_degree)));
271
272 collection_.quadrature = std::make_unique<hp::QCollection<dim>>(
273 QGaussSimplex<dim>(quadrature_degree));
274 collection_.quadrature_high_order =
275 std::make_unique<hp::QCollection<dim>>(
276 QGaussSimplex<dim>(quadrature_degree + 1));
277#if DEAL_II_VERSION_GTE(9, 7, 0)
278 collection_.nodal_quadrature = std::make_unique<hp::QCollection<dim>>(
279 FETools::compute_nodal_quadrature(FE_SimplexP<dim>(mapping_degree)));
280#else
281 AssertThrow(false,
282 dealii::ExcMessage("Discretization: Simplex support requires "
283 "deal.II version 9.7.0 or newer"));
284
285#endif
286 collection_.quadrature_1d = std::make_unique<hp::QCollection<1>>(
287 QGaussSimplex<1>(quadrature_degree));
288 collection_.face_quadrature = std::make_unique<hp::QCollection<dim - 1>>(
289 QGaussSimplex<dim - 1>(quadrature_degree));
290 if constexpr (dim == 1) {
291 collection_.face_nodal_quadrature =
292 std::make_unique<hp::QCollection<dim - 1>>(
293 QGaussLobatto<dim - 1>(quadrature_degree));
294 } else {
295#if DEAL_II_VERSION_GTE(9, 7, 0)
296 collection_.face_nodal_quadrature =
297 std::make_unique<hp::QCollection<dim - 1>>(
298 FETools::compute_nodal_quadrature(
299 FE_SimplexP<dim - 1>(mapping_degree)));
300#else
301 AssertThrow(
302 false,
303 dealii::ExcMessage("Discretization: Simplex support requires "
304 "deal.II version 9.7.0 or newer"));
305
306#endif
307 }
308
309 return;
310 } break;
311 default:
312 __builtin_trap();
313 }
314 }
315
316} /* namespace ryujin */
Discretization(const MPIEnsemble &mpi_ensemble, const std::string &subsection="/Discretization")
void prepare(const std::string &base_name)