ryujin 2.1.1 revision dbf0e3ba7acdb60b6d558e4257815df4a8f8daf9
geometry_rectangular_domain.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2022 - 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
11
12namespace ryujin
13{
14 namespace Geometries
15 {
36 template <int dim>
37 class RectangularDomain : public Geometry<dim>
38 {
39 public:
40 RectangularDomain(const std::string &name, const std::string &subsection)
41 : Geometry<dim>(name, subsection)
42 {
43 this->add_parameter("position bottom left",
44 point_left_,
45 "Position of bottom left corner");
46
47 for (unsigned int d = 0; d < dim; ++d)
48 point_right_[d] = 20.0;
49 this->add_parameter(
50 "position top right", point_right_, "Position of top right corner");
51
52 grading_push_forward_ = dim == 1 ? "x" : (dim == 2 ? "x;y" : "x;y;z;");
53 this->add_parameter("grading push forward",
54 grading_push_forward_,
55 "push forward of grading manifold");
56
57 grading_pull_back_ = grading_push_forward_;
58 this->add_parameter("grading pull back",
59 grading_pull_back_,
60 "pull back of grading manifold");
61
62 subdivisions_x_ = 1;
63 subdivisions_y_ = 1;
64 subdivisions_z_ = 1;
65 boundary_back_ = Boundary::dirichlet;
66 boundary_bottom_ = Boundary::dirichlet;
67 boundary_front_ = Boundary::dirichlet;
68 boundary_left_ = Boundary::dirichlet;
69 boundary_right_ = Boundary::dirichlet;
70 boundary_top_ = Boundary::dirichlet;
71
72 this->add_parameter("subdivisions x",
73 subdivisions_x_,
74 "number of subdivisions in x direction");
75 this->add_parameter(
76 "boundary condition left",
77 boundary_left_,
78 "Type of boundary condition enforced on the left side of the "
79 "domain (faces with normal in negative x direction)");
80 this->add_parameter(
81 "boundary condition right",
82 boundary_right_,
83 "Type of boundary condition enforced on the right side of the "
84 "domain (faces with normal in positive x direction)");
85
86 if constexpr (dim >= 2) {
87 this->add_parameter("subdivisions y",
88 subdivisions_y_,
89 "number of subdivisions in y direction");
90 this->add_parameter(
91 "boundary condition bottom",
92 boundary_bottom_,
93 "Type of boundary condition enforced on the bottom side of the "
94 "domain (faces with normal in negative y direction)");
95 this->add_parameter(
96 "boundary condition top",
97 boundary_top_,
98 "Type of boundary condition enforced on the top side of the "
99 "domain (faces with normal in positive y direction)");
100 }
101
102 if constexpr (dim == 3) {
103 this->add_parameter("subdivisions z",
104 subdivisions_z_,
105 "number of subdivisions in z direction");
106 this->add_parameter(
107 "boundary condition back",
108 boundary_back_,
109 "Type of boundary condition enforced on the back side of the "
110 "domain (faces with normal in negative z direction)");
111 this->add_parameter(
112 "boundary condition front",
113 boundary_front_,
114 "Type of boundary condition enforced on the front side of the "
115 "domain (faces with normal in positive z direction)");
116 }
117
118 use_simplices_ = false;
119 this->add_parameter("simplex mesh",
120 use_simplices_,
121 "If set to true, the triangulation uses simplices "
122 "instead of quadrangles.");
123 }
124
125
126 RectangularDomain(const std::string subsection)
127 : RectangularDomain("rectangular domain", subsection)
128 {
129 }
130
131
133 dealii::Triangulation<dim> &triangulation) const override
134 {
135 /* create mesh: */
136
137 dealii::Triangulation<dim, dim> tria1;
138 tria1.set_mesh_smoothing(triangulation.get_mesh_smoothing());
139
140 if (use_simplices_) {
141 /* Mesh with simplices: */
142
143 if constexpr (dim == 1) {
144 dealii::GridGenerator::
145 subdivided_hyper_rectangle_with_simplices<dim, dim>(
146 tria1, {subdivisions_x_}, point_left_, point_right_);
147 } else if constexpr (dim == 2) {
148 dealii::GridGenerator::subdivided_hyper_rectangle_with_simplices(
149 tria1,
150 {subdivisions_x_, subdivisions_y_},
151 point_left_,
152 point_right_);
153 } else if constexpr (dim == 3) {
154 dealii::GridGenerator::subdivided_hyper_rectangle_with_simplices(
155 tria1,
156 {subdivisions_x_, subdivisions_y_, subdivisions_z_},
157 point_left_,
158 point_right_);
159 }
160
161 } else {
162 /* Regular mesh with quadrilaterals/hexahedra: */
163
164 if constexpr (dim == 1) {
165 dealii::GridGenerator::subdivided_hyper_rectangle<dim, dim>(
166 tria1, {subdivisions_x_}, point_left_, point_right_);
167 } else if constexpr (dim == 2) {
168 dealii::GridGenerator::subdivided_hyper_rectangle(
169 tria1,
170 {subdivisions_x_, subdivisions_y_},
171 point_left_,
172 point_right_);
173 } else if constexpr (dim == 3) {
174 dealii::GridGenerator::subdivided_hyper_rectangle(
175 tria1,
176 {subdivisions_x_, subdivisions_y_, subdivisions_z_},
177 point_left_,
178 point_right_);
179 }
180 }
181
182 triangulation.copy_triangulation(tria1);
183
184 /* create grading: */
185
186 // FIXME: We should ideally check the push forward and pull back
187 // for compatiblity.
188 if (grading_push_forward_ != grading_pull_back_) {
189 dealii::FunctionManifold<dim> grading(grading_push_forward_,
190 grading_pull_back_);
191 triangulation.set_all_manifold_ids(1);
192 triangulation.set_manifold(1, grading);
193 }
194
195 /* set boundary ids: */
196
197 for (auto cell : triangulation.active_cell_iterators()) {
198 for (auto f : cell->face_indices()) {
199 auto face = cell->face(f);
200 if (!face->at_boundary())
201 continue;
202 const auto position = face->center();
203
204 if (position[0] < point_left_[0] + 1.e-8)
205 face->set_boundary_id(boundary_left_);
206 if (position[0] > point_right_[0] - 1.e-8)
207 face->set_boundary_id(boundary_right_);
208
209 if constexpr (dim >= 2) {
210 if (position[1] < point_left_[1] + 1.e-8)
211 face->set_boundary_id(boundary_bottom_);
212 if (position[1] > point_right_[1] - 1.e-8)
213 face->set_boundary_id(boundary_top_);
214 }
215
216 if constexpr (dim == 3) {
217 if (position[2] < point_left_[2] + 1.e-8)
218 face->set_boundary_id(boundary_back_);
219 if (position[2] > point_right_[2] - 1.e-8)
220 face->set_boundary_id(boundary_front_);
221 }
222 } /*for*/
223 } /*for*/
224
225 std::vector<int> directions;
226
227 if (boundary_left_ == Boundary::periodic ||
228 boundary_right_ == Boundary::periodic) {
229 AssertThrow(boundary_left_ == boundary_right_,
230 dealii::ExcMessage(
231 "For prescribing periodic boundaries in x-direction, "
232 "both, the left and right boundary conditions must "
233 "be set to periodic"));
234 directions.push_back(0);
235 }
236
237 if (dim >= 2 && (boundary_bottom_ == Boundary::periodic ||
238 boundary_top_ == Boundary::periodic)) {
239 AssertThrow(boundary_bottom_ == boundary_top_,
240 dealii::ExcMessage(
241 "For prescribing periodic boundaries in y-direction, "
242 "both, the bottom and top boundary conditions must "
243 "be set to periodic"));
244 directions.push_back(1);
245 }
246
247 if (dim == 3 && (boundary_back_ == Boundary::periodic ||
248 boundary_front_ == Boundary::periodic)) {
249 AssertThrow(boundary_back_ == boundary_front_,
250 dealii::ExcMessage(
251 "For prescribing periodic boundaries in z-direction, "
252 "both, the back and front boundary conditions must "
253 "be set to periodic"));
254 directions.push_back(2);
255 }
256
257#ifndef BUG_COLLECT_PERIODIC_FACES_INSTANTIATION
258 if (!directions.empty()) {
259 std::vector<dealii::GridTools::PeriodicFacePair<
260 typename dealii::Triangulation<dim>::cell_iterator>>
261 periodic_faces;
262
263 for (const auto direction : directions)
264 dealii::GridTools::collect_periodic_faces(
265 triangulation,
266 /*b_id */ Boundary::periodic,
267 /*direction*/ direction,
268 periodic_faces);
269
270 triangulation.add_periodicity(periodic_faces);
271 }
272#endif
273 }
274
275
277 populate_hp_collections(const unsigned int /*fe_degree*/,
278 const bool /*have_discontinuous_ansatz*/,
280 & /*collection*/) const override
281 {
282 if (use_simplices_) {
284 } else {
286 }
287 }
288
289
290 private:
291 dealii::Point<dim> point_left_;
292 dealii::Point<dim> point_right_;
293
294 unsigned int subdivisions_x_;
295 unsigned int subdivisions_y_;
296 unsigned int subdivisions_z_;
297
298 std::string grading_push_forward_;
299 std::string grading_pull_back_;
300
301 Boundary boundary_back_;
302 Boundary boundary_bottom_;
303 Boundary boundary_front_;
304 Boundary boundary_left_;
305 Boundary boundary_right_;
306 Boundary boundary_top_;
307
308 bool use_simplices_;
309 };
310 } /* namespace Geometries */
311} /* namespace ryujin */
RectangularDomain(const std::string &name, const std::string &subsection)
void create_coarse_triangulation(dealii::Triangulation< dim > &triangulation) const override
Geometry< dim >::HP_Collection populate_hp_collections(const unsigned int, const bool, typename ryujin::Discretization< dim >::Collection &) const override