ryujin 2.1.1 revision dbf0e3ba7acdb60b6d558e4257815df4a8f8daf9
geometry_step.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 GridGenerator
15 {
24 template <int dim, template <int, int> class Triangulation>
25 void step(Triangulation<dim, dim> &,
26 const double /*length*/,
27 const double /*height*/,
28 const double /*step_position*/,
29 const double /*step_height*/)
30 {
31 AssertThrow(false, dealii::ExcNotImplemented());
32 __builtin_trap();
33 }
34
35
36#ifndef DOXYGEN
37 template <template <int, int> class Triangulation>
38 void step(Triangulation<2, 2> &triangulation,
39 const double length,
40 const double height,
41 const double step_position,
42 const double step_height)
43 {
44 using namespace dealii;
45
46 dealii::Triangulation<2, 2> tria1, tria2, tria3;
47 tria3.set_mesh_smoothing(triangulation.get_mesh_smoothing());
48
49 GridGenerator::subdivided_hyper_rectangle(
50 tria1, {15, 4}, Point<2>(0., step_height), Point<2>(length, height));
51
52 GridGenerator::subdivided_hyper_rectangle(
53 tria2,
54 {3, 1},
55 Point<2>(0., 0.),
56 Point<2>(step_position, step_height));
57
58 GridGenerator::merge_triangulations(tria1, tria2, tria3);
59 triangulation.copy_triangulation(tria3);
60
61 /*
62 * Set boundary ids:
63 */
64
65 for (auto cell : triangulation.active_cell_iterators()) {
66 for (auto f : cell->face_indices()) {
67 const auto face = cell->face(f);
68
69 if (!face->at_boundary())
70 continue;
71
72 /*
73 * We want slip boundary conditions (i.e. indicator 1) at top
74 * and bottom of the rectangle. On the left we set inflow
75 * conditions (i.e. indicator 2), and we do nothing on the right
76 * side (i.e. indicator 0).
77 */
78
79 const auto center = face->center();
80
81 if (center[0] > 0. + 1.e-6 && center[0] < length - 1.e-6)
82 face->set_boundary_id(Boundary::slip);
83
84 if (center[0] < 0. + 1.e-06)
85 face->set_boundary_id(Boundary::dirichlet);
86 }
87 }
88
89 /*
90 * Refine four times and round off corner with radius 0.0125:
91 */
92
93 triangulation.refine_global(4);
94
95 Point<2> point(step_position + 0.0125, step_height - 0.0125);
96 triangulation.set_manifold(1, SphericalManifold<2>(point));
97
98 for (auto cell : triangulation.active_cell_iterators())
99 for (unsigned int v : cell->vertex_indices()) {
100 double distance =
101 (cell->vertex(v) - Point<2>(step_position, step_height)).norm();
102 if (distance < 1.e-6) {
103 for (auto f : cell->face_indices()) {
104 const auto face = cell->face(f);
105 if (face->at_boundary())
106 face->set_manifold_id(1);
107 cell->set_manifold_id(1); // temporarily for second loop
108 }
109 }
110 }
111
112 for (auto cell : triangulation.active_cell_iterators()) {
113 if (cell->manifold_id() != 1)
114 continue;
115
116 cell->set_manifold_id(0); // reset manifold id again
117
118 for (unsigned int v : cell->vertex_indices()) {
119 auto &vertex = cell->vertex(v);
120
121 if (std::abs(vertex[0] - step_position) < 1.e-6 &&
122 vertex[1] > step_height - 1.e-6)
123 vertex[0] = step_position + 0.0125 * (1 - std::sqrt(1. / 2.));
124
125 if (std::abs(vertex[1] - step_height) < 1.e-6 &&
126 vertex[0] < step_position + 0.005)
127 vertex[1] = step_height - 0.0125 * (1 - std::sqrt(1. / 2.));
128 }
129 }
130 }
131#endif
132 } /* namespace GridGenerator */
133
134
135 namespace Geometries
136 {
142 template <int dim>
143 class Step : public Geometry<dim>
144 {
145 public:
146 Step(const std::string &subsection)
147 : Geometry<dim>("step", subsection)
148 {
149 length_ = 3.;
150 this->add_parameter(
151 "length", length_, "length of computational domain");
152
153 height_ = 1.;
154 this->add_parameter(
155 "height", height_, "height of computational domain");
156
157 step_position_ = 0.6;
158 this->add_parameter(
159 "step position", step_position_, "x position of step");
160
161 step_height_ = 0.2;
162 this->add_parameter("step height", step_height_, "height of step");
163 }
164
166 dealii::Triangulation<dim> &triangulation) const final
167 {
169 triangulation, length_, height_, step_position_, step_height_);
170 }
171
172 private:
173 double length_;
174 double height_;
175 double step_position_;
176 double step_height_;
177 };
178 } /* namespace Geometries */
179
180} /* namespace ryujin */
Step(const std::string &subsection)
void create_coarse_triangulation(dealii::Triangulation< dim > &triangulation) const final
void step(Triangulation< dim, dim > &, const double, const double, const double, const double)
Definition: geometry_step.h:25