8#include <compile_time_options.h>
43 this->add_parameter(
"position bottom left",
45 "Position of bottom left corner");
47 for (
unsigned int d = 0; d < dim; ++d)
48 point_right_[d] = 20.0;
50 "position top right", point_right_,
"Position of top right corner");
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");
57 grading_pull_back_ = grading_push_forward_;
58 this->add_parameter(
"grading pull back",
60 "pull back of grading manifold");
72 this->add_parameter(
"subdivisions x",
74 "number of subdivisions in x direction");
76 "boundary condition left",
78 "Type of boundary condition enforced on the left side of the "
79 "domain (faces with normal in negative x direction)");
81 "boundary condition right",
83 "Type of boundary condition enforced on the right side of the "
84 "domain (faces with normal in positive x direction)");
86 if constexpr (dim >= 2) {
87 this->add_parameter(
"subdivisions y",
89 "number of subdivisions in y direction");
91 "boundary condition bottom",
93 "Type of boundary condition enforced on the bottom side of the "
94 "domain (faces with normal in negative y direction)");
96 "boundary condition top",
98 "Type of boundary condition enforced on the top side of the "
99 "domain (faces with normal in positive y direction)");
102 if constexpr (dim == 3) {
103 this->add_parameter(
"subdivisions z",
105 "number of subdivisions in z direction");
107 "boundary condition back",
109 "Type of boundary condition enforced on the back side of the "
110 "domain (faces with normal in negative z direction)");
112 "boundary condition front",
114 "Type of boundary condition enforced on the front side of the "
115 "domain (faces with normal in positive z direction)");
118 use_simplices_ =
false;
119 this->add_parameter(
"simplex mesh",
121 "If set to true, the triangulation uses simplices "
122 "instead of quadrangles.");
133 dealii::Triangulation<dim> &triangulation)
const override
137 dealii::Triangulation<dim, dim> tria1;
138 tria1.set_mesh_smoothing(triangulation.get_mesh_smoothing());
140 if (use_simplices_) {
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(
150 {subdivisions_x_, subdivisions_y_},
153 }
else if constexpr (dim == 3) {
154 dealii::GridGenerator::subdivided_hyper_rectangle_with_simplices(
156 {subdivisions_x_, subdivisions_y_, subdivisions_z_},
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(
170 {subdivisions_x_, subdivisions_y_},
173 }
else if constexpr (dim == 3) {
174 dealii::GridGenerator::subdivided_hyper_rectangle(
176 {subdivisions_x_, subdivisions_y_, subdivisions_z_},
182 triangulation.copy_triangulation(tria1);
188 if (grading_push_forward_ != grading_pull_back_) {
189 dealii::FunctionManifold<dim> grading(grading_push_forward_,
191 triangulation.set_all_manifold_ids(1);
192 triangulation.set_manifold(1, grading);
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())
202 const auto position = face->center();
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_);
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_);
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_);
225 std::vector<int> directions;
229 AssertThrow(boundary_left_ == boundary_right_,
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);
239 AssertThrow(boundary_bottom_ == boundary_top_,
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);
249 AssertThrow(boundary_back_ == boundary_front_,
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);
257#ifndef BUG_COLLECT_PERIODIC_FACES_INSTANTIATION
258 if (!directions.empty()) {
259 std::vector<dealii::GridTools::PeriodicFacePair<
260 typename dealii::Triangulation<dim>::cell_iterator>>
263 for (
const auto direction : directions)
264 dealii::GridTools::collect_periodic_faces(
270 triangulation.add_periodicity(periodic_faces);
282 if (use_simplices_) {
291 dealii::Point<dim> point_left_;
292 dealii::Point<dim> point_right_;
294 unsigned int subdivisions_x_;
295 unsigned int subdivisions_y_;
296 unsigned int subdivisions_z_;
298 std::string grading_push_forward_;
299 std::string grading_pull_back_;
RectangularDomain(const std::string &name, const std::string &subsection)
RectangularDomain(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