8#include <compile_time_options.h>
14 namespace GridGenerator
24 template <
int dim,
int spacedim,
template <
int,
int>
class Triangulation>
25 void annulus(Triangulation<dim, spacedim> &,
31 AssertThrow(
false, dealii::ExcNotImplemented());
37 template <
template <
int,
int>
class Triangulation>
38 void annulus(Triangulation<2, 2> &triangulation,
40 const double inner_radius,
41 const double outer_radius,
44 constexpr int dim = 2;
45 constexpr double eps = 1.0e-10;
47 using namespace dealii;
54 const auto assign_manifolds = [&](
auto &tria) {
55 for (
const auto &cell : tria.cell_iterators()) {
60 bool cell_on_annulus =
true;
61 for (
unsigned int v : cell->vertex_indices()) {
62 const auto vertex = cell->vertex(v);
63 const auto distance = vertex.norm();
64 if (!(inner_radius - eps <= distance && distance <= outer_radius)) {
65 cell_on_annulus =
false;
70 cell->set_all_manifold_ids(1);
75 for (
const unsigned int f : cell->face_indices()) {
76 const auto face = cell->face(f);
78 bool face_on_annulus =
true;
79 for (
unsigned int v : face->vertex_indices()) {
80 const auto vertex = face->vertex(v);
81 const auto distance = vertex.norm();
82 if (!(inner_radius - eps <= distance &&
83 distance <= outer_radius)) {
84 face_on_annulus =
false;
89 face->set_all_manifold_ids(1);
93 tria.set_manifold(1, SphericalManifold<dim>());
94 dealii::TransfiniteInterpolationManifold<dim> transfinite_manifold;
95 transfinite_manifold.initialize(tria);
96 tria.set_manifold(0, transfinite_manifold);
100 dealii::Triangulation<dim> tria_inner;
102 dealii::Triangulation<dim> temp;
103 GridGenerator::hyper_ball_balanced(
104 temp, dealii::Point<dim>(), inner_radius);
105 temp.refine_global(2);
106 GridGenerator::flatten_triangulation(temp, tria_inner);
110 dealii::Triangulation<dim>
annulus;
111 GridGenerator::hyper_shell(
112 annulus, dealii::Point<dim>(), inner_radius, outer_radius, 32);
115 dealii::Triangulation<dim> tria_outer;
117 dealii::Triangulation<dim> temp;
118 GridGenerator::hyper_shell(temp,
119 dealii::Point<dim>(),
121 length / 2. * std::sqrt(2),
124 for (
const auto &cell : temp.cell_iterators()) {
125 static_assert(dim == 2,
"not implemented");
126 for (
unsigned int i = 0; i < 4; ++i) {
127 auto &vertex = cell->vertex(i);
128 if (std::abs(vertex[0]) < eps && std::abs(vertex[1]) > length / 2.)
129 vertex[1] = std::copysign(length / 2., vertex[1]);
130 if (std::abs(vertex[1]) < eps && std::abs(vertex[0]) > length / 2.)
131 vertex[0] = std::copysign(length / 2., vertex[0]);
134 assign_manifolds(temp);
135 temp.refine_global(2);
136 GridGenerator::flatten_triangulation(temp, tria_outer);
140 dealii::Triangulation<dim, dim> coarse_triangulation;
141 coarse_triangulation.set_mesh_smoothing(
142 triangulation.get_mesh_smoothing());
143 GridGenerator::merge_triangulations(
144 {&tria_inner, &
annulus, &tria_outer}, coarse_triangulation, 1.e-12);
150 coarse_triangulation.reset_all_manifolds();
151 coarse_triangulation.set_all_manifold_ids(0);
153 assign_manifolds(coarse_triangulation);
154 coarse_triangulation.refine_global(2);
158 std::set<typename dealii::Triangulation<dim>::active_cell_iterator>
161 for (
const auto &cell : coarse_triangulation.active_cell_iterators()) {
162 for (
auto f : cell->face_indices()) {
163 auto face = cell->face(f);
164 const auto position = face->center();
165 const auto radius = position.norm();
166 const auto inner_value = inner_radius;
167 const auto outer_value = outer_radius;
170 radius - inner_value > 1.e-8 && outer_value - radius > 1.e-3;
172 bool partial_annulus =
173 std::abs(position[1]) -
174 std::abs(position[0]) *
175 std::tan(dealii::numbers::PI / 180. * angle) <
178 if (in_anulus && partial_annulus) {
179 cells_to_remove.insert(cell);
184 GridGenerator::create_triangulation_with_removed_cells(
185 coarse_triangulation, cells_to_remove, coarse_triangulation);
190 dealii::Triangulation<dim> flattened_triangulation;
191 flattened_triangulation.set_mesh_smoothing(
192 triangulation.get_mesh_smoothing());
193 GridGenerator::flatten_triangulation(coarse_triangulation,
194 flattened_triangulation);
195 triangulation.copy_triangulation(flattened_triangulation);
196 assign_manifolds(triangulation);
203 for (
auto cell : triangulation.active_cell_iterators()) {
204 for (
auto f : cell->face_indices()) {
205 const auto face = cell->face(f);
207 if (!face->at_boundary())
219 template <
template <
int,
int>
class Triangulation>
220 void annulus(Triangulation<3, 3> & ,
226 using namespace dealii;
227 AssertThrow(
false, dealii::ExcNotImplemented());
247 :
Geometry<dim>(
"annulus", subsection)
251 "length", length_,
"length of computational domain [-L/2,L/2]^d");
255 "inner radius", inner_radius_,
"inner radius of partial annulus");
259 "outer radius", outer_radius_,
"outer radius of partial annulus");
262 this->add_parameter(
"coverage angle",
264 "angle coverage of partial annulus above y-axis");
268 dealii::Triangulation<dim> &triangulation)
const final
271 triangulation, length_, inner_radius_, outer_radius_, angle_);
276 double inner_radius_;
277 double outer_radius_;
Annulus(const std::string &subsection)
void create_coarse_triangulation(dealii::Triangulation< dim > &triangulation) const final
void annulus(Triangulation< dim, spacedim > &, const double, const double, const double, const double)