ryujin 2.1.1 revision dbf0e3ba7acdb60b6d558e4257815df4a8f8daf9
geometry_annulus.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
11
12namespace ryujin
13{
14 namespace GridGenerator
15 {
24 template <int dim, int spacedim, template <int, int> class Triangulation>
25 void annulus(Triangulation<dim, spacedim> &,
26 const double /*length*/,
27 const double /*inner_radius*/,
28 const double /*outer_radius*/,
29 const double /*angle*/)
30 {
31 AssertThrow(false, dealii::ExcNotImplemented());
32 __builtin_trap();
33 }
34
35
36#ifndef DOXYGEN
37 template <template <int, int> class Triangulation>
38 void annulus(Triangulation<2, 2> &triangulation,
39 const double length,
40 const double inner_radius,
41 const double outer_radius,
42 const double angle)
43 {
44 constexpr int dim = 2;
45 constexpr double eps = 1.0e-10;
46
47 using namespace dealii;
48
49 /*
50 * A lambda that conveniently sets and assigns manifold IDs and that
51 * we need a couple of times during the construction process of the
52 * triangulation:
53 */
54 const auto assign_manifolds = [&](auto &tria) {
55 for (const auto &cell : tria.cell_iterators()) {
56 /*
57 * Mark all cells that will comprise the annulus. That is, all
58 * cells whose vertices lie betwenn the two radii.
59 */
60 bool cell_on_annulus = true;
61 for (unsigned int v : cell->vertex_indices()) {
62 const auto vertex = cell->vertex(v);
63 const auto distance = vertex.norm();
64 if (!(inner_radius - eps <= distance && distance <= outer_radius)) {
65 cell_on_annulus = false;
66 break;
67 }
68 }
69 if (cell_on_annulus)
70 cell->set_all_manifold_ids(1);
71
72 /*
73 * Separately, mark all faces that touch the annulus.
74 */
75 for (const unsigned int f : cell->face_indices()) {
76 const auto face = cell->face(f);
77
78 bool face_on_annulus = true;
79 for (unsigned int v : face->vertex_indices()) {
80 const auto vertex = face->vertex(v);
81 const auto distance = vertex.norm();
82 if (!(inner_radius - eps <= distance &&
83 distance <= outer_radius)) {
84 face_on_annulus = false;
85 break;
86 }
87 }
88 if (face_on_annulus)
89 face->set_all_manifold_ids(1);
90 }
91 }
92
93 tria.set_manifold(1, SphericalManifold<dim>());
94 dealii::TransfiniteInterpolationManifold<dim> transfinite_manifold;
95 transfinite_manifold.initialize(tria);
96 tria.set_manifold(0, transfinite_manifold);
97 };
98
99 /* Create inner ball with radius=inner_radius: */
100 dealii::Triangulation<dim> tria_inner;
101 {
102 dealii::Triangulation<dim> temp;
103 GridGenerator::hyper_ball_balanced(
104 temp, dealii::Point<dim>(), inner_radius);
105 temp.refine_global(2);
106 GridGenerator::flatten_triangulation(temp, tria_inner);
107 }
108
109 /* Create outer annulus. Note part of this will be removed. */
110 dealii::Triangulation<dim> annulus;
111 GridGenerator::hyper_shell(
112 annulus, dealii::Point<dim>(), inner_radius, outer_radius, 32);
113
114 /* Create outside shell */
115 dealii::Triangulation<dim> tria_outer;
116 {
117 dealii::Triangulation<dim> temp;
118 GridGenerator::hyper_shell(temp,
119 dealii::Point<dim>(),
120 outer_radius,
121 length / 2. * std::sqrt(2),
122 8);
123 /* Fix up vertices so that we get back a unit square: */
124 for (const auto &cell : temp.cell_iterators()) {
125 static_assert(dim == 2, "not implemented");
126 for (unsigned int i = 0; i < 4; ++i) {
127 auto &vertex = cell->vertex(i);
128 if (std::abs(vertex[0]) < eps && std::abs(vertex[1]) > length / 2.)
129 vertex[1] = std::copysign(length / 2., vertex[1]);
130 if (std::abs(vertex[1]) < eps && std::abs(vertex[0]) > length / 2.)
131 vertex[0] = std::copysign(length / 2., vertex[0]);
132 }
133 }
134 assign_manifolds(temp);
135 temp.refine_global(2);
136 GridGenerator::flatten_triangulation(temp, tria_outer);
137 }
138
139 /* Create triangulation to merge: */
140 dealii::Triangulation<dim, dim> coarse_triangulation;
141 coarse_triangulation.set_mesh_smoothing(
142 triangulation.get_mesh_smoothing());
143 GridGenerator::merge_triangulations(
144 {&tria_inner, &annulus, &tria_outer}, coarse_triangulation, 1.e-12);
145
146 /*
147 * Set manifold IDs:
148 */
149
150 coarse_triangulation.reset_all_manifolds();
151 coarse_triangulation.set_all_manifold_ids(0);
152
153 assign_manifolds(coarse_triangulation);
154 coarse_triangulation.refine_global(2);
155
156 /* Remove mesh cells in the annulus */
157
158 std::set<typename dealii::Triangulation<dim>::active_cell_iterator>
159 cells_to_remove;
160
161 for (const auto &cell : coarse_triangulation.active_cell_iterators()) {
162 for (auto f : cell->face_indices()) {
163 auto face = cell->face(f);
164 const auto position = face->center();
165 const auto radius = position.norm();
166 const auto inner_value = inner_radius;
167 const auto outer_value = outer_radius;
168
169 bool in_anulus =
170 radius - inner_value > 1.e-8 && outer_value - radius > 1.e-3;
171
172 bool partial_annulus =
173 std::abs(position[1]) -
174 std::abs(position[0]) *
175 std::tan(dealii::numbers::PI / 180. * angle) <
176 1.e-8;
177
178 if (in_anulus && partial_annulus) {
179 cells_to_remove.insert(cell);
180 }
181 }
182 }
183
184 GridGenerator::create_triangulation_with_removed_cells(
185 coarse_triangulation, cells_to_remove, coarse_triangulation);
186
187 /*
188 * Flatten triangulation and copy over to distributed triangulation:
189 */
190 dealii::Triangulation<dim> flattened_triangulation;
191 flattened_triangulation.set_mesh_smoothing(
192 triangulation.get_mesh_smoothing());
193 GridGenerator::flatten_triangulation(coarse_triangulation,
194 flattened_triangulation);
195 triangulation.copy_triangulation(flattened_triangulation);
196 assign_manifolds(triangulation);
197
198
199 /*
200 * Set boundary ids:
201 */
202
203 for (auto cell : triangulation.active_cell_iterators()) {
204 for (auto f : cell->face_indices()) {
205 const auto face = cell->face(f);
206
207 if (!face->at_boundary())
208 continue;
209
210 /*
211 * We want slip boundary conditions everywhere.
212 */
213 face->set_boundary_id(Boundary::slip);
214 }
215 }
216 }
217
218
219 template <template <int, int> class Triangulation>
220 void annulus(Triangulation<3, 3> & /* triangulation */,
221 const double /* length */,
222 const double /* inner_radius */,
223 const double /* outer_radius */,
224 const double /* angle */)
225 {
226 using namespace dealii;
227 AssertThrow(false, dealii::ExcNotImplemented());
228 __builtin_trap();
229 }
230#endif
231 } /* namespace GridGenerator */
232
233
234 namespace Geometries
235 {
242 template <int dim>
243 class Annulus : public Geometry<dim>
244 {
245 public:
246 Annulus(const std::string &subsection)
247 : Geometry<dim>("annulus", subsection)
248 {
249 length_ = 2.;
250 this->add_parameter(
251 "length", length_, "length of computational domain [-L/2,L/2]^d");
252
253 inner_radius_ = 0.6;
254 this->add_parameter(
255 "inner radius", inner_radius_, "inner radius of partial annulus");
256
257 outer_radius_ = 0.7;
258 this->add_parameter(
259 "outer radius", outer_radius_, "outer radius of partial annulus");
260
261 angle_ = 45.;
262 this->add_parameter("coverage angle",
263 angle_,
264 "angle coverage of partial annulus above y-axis");
265 }
266
268 dealii::Triangulation<dim> &triangulation) const final
269 {
271 triangulation, length_, inner_radius_, outer_radius_, angle_);
272 }
273
274 private:
275 double length_;
276 double inner_radius_;
277 double outer_radius_;
278 double angle_;
279 };
280 } /* namespace Geometries */
281} /* namespace ryujin */
Annulus(const std::string &subsection)
void create_coarse_triangulation(dealii::Triangulation< dim > &triangulation) const final
void annulus(Triangulation< dim, spacedim > &, const double, const double, const double, const double)