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