8#include <compile_time_options.h>
14 namespace GridGenerator
30 template <
int dim,
int spacedim,
template <
int,
int>
class Triangulation>
37 AssertThrow(
false, dealii::ExcNotImplemented());
43 template <
template <
int,
int>
class Triangulation>
44 void cylinder(Triangulation<2, 2> &triangulation,
47 const double cylinder_position,
48 const double cylinder_diameter)
50 constexpr int dim = 2;
52 using namespace dealii;
54 dealii::Triangulation<dim, dim> tria1, tria2, tria3, tria4, tria5, tria6,
57 GridGenerator::hyper_cube_with_cylindrical_hole(
58 tria1, cylinder_diameter / 2., cylinder_diameter, 0.5, 1,
false);
60 GridGenerator::subdivided_hyper_rectangle(
63 Point<2>(-cylinder_diameter, -cylinder_diameter),
64 Point<2>(cylinder_diameter, -height / 2.));
66 GridGenerator::subdivided_hyper_rectangle(
69 Point<2>(-cylinder_diameter, cylinder_diameter),
70 Point<2>(cylinder_diameter, height / 2.));
72 GridGenerator::subdivided_hyper_rectangle(
75 Point<2>(cylinder_diameter, -cylinder_diameter),
76 Point<2>(length - cylinder_position, cylinder_diameter));
78 GridGenerator::subdivided_hyper_rectangle(
81 Point<2>(cylinder_diameter, cylinder_diameter),
82 Point<2>(length - cylinder_position, height / 2.));
84 GridGenerator::subdivided_hyper_rectangle(
87 Point<2>(cylinder_diameter, -height / 2.),
88 Point<2>(length - cylinder_position, -cylinder_diameter));
90 tria7.set_mesh_smoothing(triangulation.get_mesh_smoothing());
91 GridGenerator::merge_triangulations(
92 {&tria1, &tria2, &tria3, &tria4, &tria5, &tria6},
96 triangulation.copy_triangulation(tria7);
100 triangulation.set_manifold(0, PolarManifold<2>(Point<2>()));
104 for (
auto cell : triangulation.active_cell_iterators())
105 for (
unsigned int v : cell->vertex_indices()) {
106 auto &vertex = cell->vertex(v);
107 if (vertex[0] <= -cylinder_diameter + 1.e-6)
108 vertex[0] = -cylinder_position;
115 for (
auto cell : triangulation.active_cell_iterators()) {
116 for (
auto f : cell->face_indices()) {
117 const auto face = cell->face(f);
119 if (!face->at_boundary())
129 const auto center = face->center();
131 if (center[0] > length - cylinder_position - 1.e-6) {
136 if (center[0] < -cylinder_position + 1.e-6) {
148 template <
template <
int,
int>
class Triangulation>
149 void cylinder(Triangulation<3, 3> &triangulation,
152 const double cylinder_position,
153 const double cylinder_diameter)
155 using namespace dealii;
157 dealii::Triangulation<2, 2> tria1;
159 cylinder(tria1, length, height, cylinder_position, cylinder_diameter);
161 dealii::Triangulation<3, 3> tria2;
162 tria2.set_mesh_smoothing(triangulation.get_mesh_smoothing());
164 GridGenerator::extrude_triangulation(tria1, 4, height, tria2,
true);
165 GridTools::transform(
166 [height](
auto point) {
167 return point - dealii::Tensor<1, 3>{{0, 0, height / 2.}};
171 triangulation.copy_triangulation(tria2);
177 triangulation.set_manifold(
178 0, CylindricalManifold<3>(Tensor<1, 3>{{0., 0., 1.}}, Point<3>()));
184 for (
auto cell : triangulation.active_cell_iterators()) {
185 for (
auto f : cell->face_indices()) {
186 const auto face = cell->face(f);
188 if (!face->at_boundary())
198 const auto center = face->center();
200 if (center[0] > length - cylinder_position - 1.e-6) {
205 if (center[0] < -cylinder_position + 1.e-6) {
232 :
Geometry<dim>(
"cylinder", subsection)
236 "length", length_,
"length of computational domain");
240 "height", height_,
"height of computational domain");
242 object_position_ = 0.6;
243 this->add_parameter(
"object position",
245 "x position of immersed cylinder center point");
247 object_diameter_ = 0.5;
248 this->add_parameter(
"object diameter",
250 "diameter of immersed cylinder");
254 dealii::Triangulation<dim> &triangulation)
const final
266 double object_position_;
267 double object_diameter_;
void create_coarse_triangulation(dealii::Triangulation< dim > &triangulation) const final
Cylinder(const std::string &subsection)
void cylinder(Triangulation< dim, spacedim > &, const double, const double, const double, const double)