12 namespace GridGenerator
27 template <
int dim,
int spacedim,
template <
int,
int>
class Triangulation>
34 AssertThrow(
false, dealii::ExcNotImplemented());
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)
47 using namespace dealii;
49 dealii::Triangulation<2, 2> res1, res2, res3, flume,
final;
51 const double tolerance = 1.e-8;
53 Assert(reservoir_width - flume_width > tolerance,
54 dealii::ExcInternalError());
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.));
63 GridGenerator::subdivided_hyper_rectangle(
66 Point<2>(-reservoir_length, -reservoir_width / 2.),
67 Point<2>(0, -flume_width / 2.));
69 GridGenerator::subdivided_hyper_rectangle(
72 Point<2>(-reservoir_length, flume_width / 2.),
73 Point<2>(0, reservoir_width / 2.));
75 sub_y = int(std::round(flume_width * 100.));
77 GridGenerator::subdivided_hyper_rectangle(
80 Point<2>(-reservoir_length, -flume_width / 2.),
81 Point<2>(0, flume_width / 2.));
83 sub_x = int(std::round(flume_length * 100.));
85 GridGenerator::subdivided_hyper_rectangle(
88 Point<2>(0., -flume_width / 2.),
89 Point<2>(flume_length, flume_width / 2.));
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);
100 for (
auto cell : triangulation.active_cell_iterators()) {
101 for (
auto f : cell->face_indices()) {
102 const auto face = cell->face(f);
104 if (!face->at_boundary())
113 const auto center = face->center();
114 if (center[0] > flume_length - tolerance)
136 :
Geometry<dim>(
"wave tank", subsection)
138 reservoir_length_ = 157. / 100.;
139 this->add_parameter(
"reservoir length",
141 "length of water reservoir [meters]");
143 reservoir_width_ = 8.1 / 100.;
144 this->add_parameter(
"reservoir width",
146 "width of water reservoir [meters]");
148 flume_length_ = 600.78 / 100.;
150 "flume length", flume_length_,
"length of flume [meters]");
152 flume_width_ = 24 / 100.;
154 "flume width", flume_width_,
"width of flume [meters]");
168 double reservoir_length_;
169 double reservoir_width_;
170 double flume_length_;
WaveTank(const std::string subsection)
void create_triangulation(typename Geometry< dim >::Triangulation &triangulation) final
typename Discretization< dim >::Triangulation Triangulation
void wavetank(Triangulation< dim, spacedim > &, const double, const double, const double, const double)