ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
discretization.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
10#include "convenience_macros.h"
11#include "geometries/geometry.h"
12#include "mpi_ensemble.h"
13#include "patterns_conversion.h"
14
15#include <deal.II/base/parameter_acceptor.h>
16#include <deal.II/distributed/fully_distributed_tria.h>
17#include <deal.II/distributed/shared_tria.h>
18#include <deal.II/distributed/tria.h>
19#include <deal.II/hp/fe_collection.h>
20#include <deal.II/hp/mapping_collection.h>
21#include <deal.II/hp/q_collection.h>
22
23#include <memory>
24#include <set>
25
26namespace ryujin
27{
53 enum Boundary : dealii::types::boundary_id {
63
70
80 slip = 2,
81
91
103
119
128 };
129
130
136 enum class Ansatz {
138 cg_q1,
139
141 cg_q2,
142
144 cg_q3,
145
147 dg_q1,
148
150 dg_q2,
151
153 dg_q3
154 };
155
162 enum class MeshType {
164 serial,
171 };
172} // namespace ryujin
173
174#ifndef DOXYGEN
175DECLARE_ENUM(ryujin::Boundary,
176 LIST({ryujin::Boundary::do_nothing, "do nothing"},
177 {ryujin::Boundary::periodic, "periodic"},
178 {ryujin::Boundary::slip, "slip"},
179 {ryujin::Boundary::no_slip, "no slip"},
180 {ryujin::Boundary::dirichlet, "dirichlet"},
181 {ryujin::Boundary::dynamic, "dynamic"},
183 "dirichlet momentum"}));
184
185DECLARE_ENUM(ryujin::Ansatz,
186 LIST({ryujin::Ansatz::cg_q1, "cG Q1"},
187 {ryujin::Ansatz::cg_q2, "cG Q2"},
188 {ryujin::Ansatz::cg_q3, "cG Q3"},
189 {ryujin::Ansatz::dg_q1, "dG Q1"},
190 {ryujin::Ansatz::dg_q2, "dG Q2"},
191 {ryujin::Ansatz::dg_q3, "dG Q3"}));
192
193DECLARE_ENUM(ryujin::MeshType,
194 LIST({ryujin::MeshType::serial, "serial"},
195 {ryujin::MeshType::parallel_shared, "parallel shared"},
197 "parallel distributed"},
199 "parallel fullydistributed"}));
200#endif
201
202namespace ryujin
203{
219 template <int dim>
220 class Discretization : public dealii::ParameterAcceptor
221 {
222 public:
226 Discretization(const MPIEnsemble &mpi_ensemble,
227 const std::string &subsection = "/Discretization");
228
233 void prepare(const std::string &base_name);
234
241 struct Collection {
242 std::unique_ptr<const dealii::hp::MappingCollection<dim>> mapping;
243 std::unique_ptr<const dealii::hp::FECollection<dim>> finite_element;
244 std::unique_ptr<const dealii::hp::FECollection<dim>> finite_element_cg;
245 std::unique_ptr<const dealii::hp::QCollection<dim>> quadrature;
246 std::unique_ptr<const dealii::hp::QCollection<dim>> quadrature_high_order;
247 std::unique_ptr<const dealii::hp::QCollection<dim>> nodal_quadrature;
248 std::unique_ptr<const dealii::hp::QCollection<1>> quadrature_1d;
249 std::unique_ptr<const dealii::hp::QCollection<dim - 1>> face_quadrature;
250 std::unique_ptr<const dealii::hp::QCollection<dim - 1>>
252 };
253
258
263
264
268
274 {
275 switch (ansatz_) {
276 /* Continuous Ansatz: */
277 case Ansatz::cg_q1:
278 [[fallthrough]];
279 case Ansatz::cg_q2:
280 [[fallthrough]];
281 case Ansatz::cg_q3:
282 return false;
283
284 /* Discontinuous Ansatz: */
285 case Ansatz::dg_q1:
286 [[fallthrough]];
287 case Ansatz::dg_q2:
288 [[fallthrough]];
289 case Ansatz::dg_q3:
290 return true;
291 }
292
293 AssertThrow(false, dealii::ExcInternalError());
294 __builtin_trap();
295 }
296
300 unsigned int polynomial_degree() const
301 {
302 switch (ansatz_) {
303 case Ansatz::cg_q1:
304 [[fallthrough]];
305 case Ansatz::dg_q1:
306 return 1;
307 case Ansatz::cg_q2:
308 [[fallthrough]];
309 case Ansatz::dg_q2:
310 return 2;
311 case Ansatz::cg_q3:
312 [[fallthrough]];
313 case Ansatz::dg_q3:
314 return 3;
315 }
316
317 AssertThrow(false, dealii::ExcInternalError());
318 __builtin_trap();
319 }
320
324 ACCESSOR(refinement)
325
326
330
331
340 ACCESSOR_CONTAINER_READ_ONLY(collection_, mapping)
341
347 ACCESSOR_CONTAINER_READ_ONLY(collection_, finite_element)
348
355 ACCESSOR_CONTAINER_READ_ONLY(collection_, finite_element_cg)
356
362 ACCESSOR_CONTAINER_READ_ONLY(collection_, quadrature)
363
370 ACCESSOR_CONTAINER_READ_ONLY(collection_, quadrature_high_order)
371
378 ACCESSOR_CONTAINER_READ_ONLY(collection_, nodal_quadrature)
379
385 ACCESSOR_CONTAINER_READ_ONLY(collection_, quadrature_1d)
386
392 ACCESSOR_CONTAINER_READ_ONLY(collection_, face_quadrature)
393
400 ACCESSOR_CONTAINER_READ_ONLY(collection_, face_nodal_quadrature)
401
402 private:
404
408
409 Ansatz ansatz_;
410 MeshType mesh_type_;
411
412 std::string geometry_;
413
414 unsigned int refinement_;
415
416 bool mesh_writeout_;
417 double mesh_distortion_;
418
420
424 //
425 const MPIEnsemble &mpi_ensemble_;
426
427 std::unique_ptr<dealii::Triangulation<dim>> triangulation_;
428
429 Collection collection_;
430
431 std::set<std::shared_ptr<Geometry<dim>>> geometry_list_;
432 std::shared_ptr<Geometry<dim>> selected_geometry_;
433
435
443 template <typename Discretization, int dim_, typename Number_>
444 friend class SolutionTransfer;
445
453 template <int dim_>
454 friend class Geometry;
455 };
456} /* namespace ryujin */
const auto & triangulation() const
const auto & selected_geometry() const
const auto & ansatz() const
Discretization(const MPIEnsemble &mpi_ensemble, const std::string &subsection="/Discretization")
void prepare(const std::string &base_name)
unsigned int polynomial_degree() const
bool have_discontinuous_ansatz() const
@ dirichlet_momentum
#define ACCESSOR(member)
#define ACCESSOR_CONTAINER_READ_ONLY(container, member)
#define ACCESSOR_READ_ONLY(member)
std::unique_ptr< const dealii::hp::FECollection< dim > > finite_element
std::unique_ptr< const dealii::hp::QCollection< 1 > > quadrature_1d
std::unique_ptr< const dealii::hp::FECollection< dim > > finite_element_cg
std::unique_ptr< const dealii::hp::MappingCollection< dim > > mapping
std::unique_ptr< const dealii::hp::QCollection< dim > > quadrature_high_order
std::unique_ptr< const dealii::hp::QCollection< dim > > quadrature
std::unique_ptr< const dealii::hp::QCollection< dim - 1 > > face_quadrature
std::unique_ptr< const dealii::hp::QCollection< dim - 1 > > face_nodal_quadrature
std::unique_ptr< const dealii::hp::QCollection< dim > > nodal_quadrature