ryujin 2.1.1 revision a15074459a388761bd8df6bd4ef7e6abe9d8077b
geometry_tank.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2024 by the ryujin authors
4//
5
6#pragma once
7
9
10namespace ryujin
11{
12 namespace GridGenerator
13 {
27 template <int dim, int spacedim, template <int, int> class Triangulation>
28 void wavetank(Triangulation<dim, spacedim> &,
29 const double /*reservoir_length*/,
30 const double /*reservoir_width*/,
31 const double /*flume_length*/,
32 const double /*flume_width*/)
33 {
34 AssertThrow(false, dealii::ExcNotImplemented());
35 __builtin_trap();
36 }
37
38
39#ifndef DOXYGEN
40 template <template <int, int> class Triangulation>
41 void wavetank(Triangulation<2, 2> &triangulation,
42 const double reservoir_length,
43 const double reservoir_width,
44 const double flume_length,
45 const double flume_width)
46 {
47 using namespace dealii;
48
49 dealii::Triangulation<2, 2> res1, res2, res3, flume, final;
50
51 const double tolerance = 1.e-8;
52
53 Assert(reservoir_width - flume_width > tolerance,
54 dealii::ExcInternalError());
55
56 /* We split the reservoir into three triangulations and subdivide to
57 * get close to uniform refinement */
58
59 const double diff = (reservoir_width - flume_width) / 2.;
60 unsigned int sub_x = int(std::round(reservoir_length * 100.));
61 unsigned int sub_y = int(std::round(diff * 100.));
62
63 GridGenerator::subdivided_hyper_rectangle(
64 res1,
65 {sub_x, sub_y},
66 Point<2>(-reservoir_length, -reservoir_width / 2.),
67 Point<2>(0, -flume_width / 2.));
68
69 GridGenerator::subdivided_hyper_rectangle(
70 res3,
71 {sub_x, sub_y},
72 Point<2>(-reservoir_length, flume_width / 2.),
73 Point<2>(0, reservoir_width / 2.));
74
75 sub_y = int(std::round(flume_width * 100.));
76
77 GridGenerator::subdivided_hyper_rectangle(
78 res2,
79 {sub_x, sub_y},
80 Point<2>(-reservoir_length, -flume_width / 2.),
81 Point<2>(0, flume_width / 2.));
82
83 sub_x = int(std::round(flume_length * 100.));
84
85 GridGenerator::subdivided_hyper_rectangle(
86 flume,
87 {sub_x, sub_y},
88 Point<2>(0., -flume_width / 2.),
89 Point<2>(flume_length, flume_width / 2.));
90
91 final.set_mesh_smoothing(triangulation.get_mesh_smoothing());
92 GridGenerator::merge_triangulations(
93 {&res1, &res2, &res3, &flume}, final, tolerance);
94 triangulation.copy_triangulation(final);
95
96 /*
97 * Set boundary ids:
98 */
99
100 for (auto cell : triangulation.active_cell_iterators()) {
101 for (auto f : GeometryInfo<2>::face_indices()) {
102 const auto face = cell->face(f);
103
104 if (!face->at_boundary())
105 continue;
106
107 /*
108 * We want slip everywhere except right edge of tank.
109 */
110
111 face->set_boundary_id(Boundary::slip);
112
113 const auto center = face->center();
114 if (center[0] > flume_length - tolerance)
115 face->set_boundary_id(Boundary::dynamic);
116
117 } /*f*/
118 } /*cell*/
119 }
120#endif
121 } /* namespace GridGenerator */
122
123
124 namespace Geometries
125 {
131 template <int dim>
132 class WaveTank : public Geometry<dim>
133 {
134 public:
135 WaveTank(const std::string subsection)
136 : Geometry<dim>("wave tank", subsection)
137 {
138 reservoir_length_ = 157. / 100.;
139 this->add_parameter("reservoir length",
140 reservoir_length_,
141 "length of water reservoir [meters]");
142
143 reservoir_width_ = 8.1 / 100.;
144 this->add_parameter("reservoir width",
145 reservoir_width_,
146 "width of water reservoir [meters]");
147
148 flume_length_ = 600.78 / 100.;
149 this->add_parameter(
150 "flume length", flume_length_, "length of flume [meters]");
151
152 flume_width_ = 24 / 100.;
153 this->add_parameter(
154 "flume width", flume_width_, "width of flume [meters]");
155 }
156
158 typename Geometry<dim>::Triangulation &triangulation) final
159 {
160 GridGenerator::wavetank(triangulation,
161 reservoir_length_,
162 reservoir_width_,
163 flume_length_,
164 flume_width_);
165 }
166
167 private:
168 double reservoir_length_;
169 double reservoir_width_;
170 double flume_length_;
171 double flume_width_;
172 };
173 } /* namespace Geometries */
174} /* namespace ryujin */
WaveTank(const std::string subsection)
void create_triangulation(typename Geometry< dim >::Triangulation &triangulation) final
typename Discretization< dim >::Triangulation Triangulation
Definition: geometry.h:38
void wavetank(Triangulation< dim, spacedim > &, const double, const double, const double, const double)
Definition: geometry_tank.h:28