42 double skeleton_length;
44 double cosecant_polar;
45 double cotangent_polar;
46 double cotangent_azimuth_half_right;
47 double cotangent_azimuth_half_left;
57 compute_z_expansion(
const double x,
59 const AdditionalData &data)
66 + x * data.cotangent_polar
68 -
std::abs(y) * data.cosecant_polar *
69 ((y > 0) ? data.cotangent_azimuth_half_right :
70 data.cotangent_azimuth_half_left);
81 template <
int dim,
int spacedim = dim>
107 const AdditionalData & data,
108 const double tolerance = 1e-10);
113 virtual std::unique_ptr<::Manifold<dim, spacedim>>
114 clone()
const override;
133 push_forward(
const Point<3> &chart_point)
const override;
154 const AdditionalData data;
159 const double tolerance;
170 template <
int dim,
int spacedim>
171 Manifold<dim, spacedim>::Manifold(
175 const AdditionalData & data,
176 const double tolerance)
178 , normal_direction(normal_direction)
179 , direction(direction)
180 , point_on_axis(point_on_axis)
182 , tolerance(tolerance)
183 , dxn(cross_product_3d(direction, normal_direction))
187 "PipeSegment::Manifold can only be used for spacedim==3!"));
190 ExcMessage(
"Normal direction must be unit vector."));
192 ExcMessage(
"Direction must be unit vector."));
193 Assert(normal_direction * direction < tolerance,
195 "Direction and normal direction must be perpendicular."));
200 template <
int dim,
int spacedim>
201 std::unique_ptr<::Manifold<dim, spacedim>>
204 return std::make_unique<Manifold<dim, spacedim>>(*this);
209 template <
int dim,
int spacedim>
211 Manifold<dim, spacedim>::pull_back(
const Point<spacedim> &space_point)
const
215 double lambda = normalized_point * direction;
218 const double r = p_diff.
norm();
220 Assert(r > tolerance * data.skeleton_length,
222 "This class won't handle points on the direction axis."));
235 return {r, phi, lambda};
240 template <
int dim,
int spacedim>
242 Manifold<dim, spacedim>::push_forward(
const Point<3> &chart_point)
const
245 const double sine_r = chart_point(0) *
std::sin(chart_point(1));
246 const double cosine_r = chart_point(0) *
std::cos(chart_point(1));
249 normal_direction * cosine_r + dxn * sine_r;
252 const double lambda =
253 chart_point(2) * compute_z_expansion(cosine_r, sine_r, data);
256 return point_on_axis + direction * lambda + intermediate;
265 template <
int dim,
int spacedim>
283 const std::vector<std::pair<
Point<3>,
double>> &openings,
284 const std::pair<
Point<3>,
double> & bifurcation,
285 const double aspect_ratio)
287 constexpr unsigned int dim = 3;
288 constexpr unsigned int spacedim = 3;
291 constexpr unsigned int n_pipes = 3;
292 constexpr double tolerance = 1.e-12;
296 Assert(bifurcation.second > 0,
297 ExcMessage(
"Invalid input: negative radius."));
298 Assert(openings.size() == n_pipes,
299 ExcMessage(
"Invalid input: only 3 openings allowed."));
300 for (
const auto &opening : openings)
307 const auto cyclic = [](
const unsigned int i) ->
unsigned int {
308 constexpr unsigned int n_pipes = 3;
309 return (i < (n_pipes - 1)) ? i + 1 : 0;
311 const auto anticyclic = [](
const unsigned int i) ->
unsigned int {
312 constexpr unsigned int n_pipes = 3;
313 return (i > 0) ? i - 1 : n_pipes - 1;
317 const std::array<vector3d, spacedim> directions = {
318 {vector3d({1., 0., 0.}), vector3d({0., 1., 0.}), vector3d({0., 0., 1.})}};
324 std::array<vector3d, n_pipes> skeleton;
325 for (
unsigned int p = 0; p < n_pipes; ++p)
326 skeleton[p] = bifurcation.first - openings[p].first;
328 std::array<double, n_pipes> skeleton_length;
329 for (
unsigned int p = 0; p < n_pipes; ++p)
330 skeleton_length[p] = skeleton[p].
norm();
336 const double tolerance_length =
338 *std::max_element(skeleton_length.begin(), skeleton_length.end());
340 std::array<vector3d, n_pipes> skeleton_unit;
341 for (
unsigned int p = 0; p < n_pipes; ++p)
343 Assert(skeleton_length[p] > tolerance_length,
344 ExcMessage(
"Invalid input: bifurcation matches opening."));
345 skeleton_unit[p] = skeleton[p] / skeleton_length[p];
357 vector3d normal = cross_product_3d(skeleton_unit[1] - skeleton_unit[0],
358 skeleton_unit[2] - skeleton_unit[0]);
359 Assert(normal.norm() > tolerance_length,
360 ExcMessage(
"Invalid input: all three openings "
361 "are located on one line."));
362 normal /= normal.norm();
366 std::array<vector3d, n_pipes> skeleton_plane;
367 for (
unsigned int p = 0; p < n_pipes; ++p)
369 skeleton_plane[p] = skeleton[p] - (skeleton[p] * normal) * normal;
371 tolerance * skeleton_plane[p].
norm(),
373 Assert(skeleton_plane[p].
norm() > tolerance_length,
389 ExcMessage(
"The output triangulation object needs to be empty."));
391 std::vector<PipeSegment::Manifold<dim, spacedim>> manifolds;
392 for (
unsigned int p = 0; p < n_pipes; ++p)
405 const unsigned int n_slices =
406 1 +
static_cast<unsigned int>(std::ceil(
407 aspect_ratio * skeleton_length[p] /
418 for (
const auto &cell : pipe.active_cell_iterators())
420 cell->set_material_id(p);
422 for (
const auto &face : cell->face_iterators())
423 if (face->at_boundary())
425 const auto center_z = face->center()[2];
431 else if (
std::abs(center_z - 1.) < tolerance)
438 cell->set_all_manifold_ids(n_pipes);
439 face->set_all_manifold_ids(p);
474 const double polar_angle =
483 const double azimuth_angle_right =
485 skeleton_plane[cyclic(p)],
488 ExcMessage(
"Invalid input: at least two openings located "
489 "in same direction from bifurcation"));
491 const double azimuth_angle_left =
493 skeleton_plane[anticyclic(p)],
496 ExcMessage(
"Invalid input: at least two openings located "
497 "in same direction from bifurcation"));
501 PipeSegment::AdditionalData data;
502 data.skeleton_length = skeleton_length[p];
503 data.cosecant_polar = 1. /
std::sin(polar_angle);
504 data.cotangent_polar =
std::cos(polar_angle) * data.cosecant_polar;
505 data.cotangent_azimuth_half_right =
std::cos(.5 * azimuth_angle_right) /
507 data.cotangent_azimuth_half_left =
511 const auto pipe_segment_transform =
515 const double r_factor =
516 (bifurcation.second - openings[p].second) * pt[2] +
518 const double x_new = r_factor * pt[0];
519 const double y_new = r_factor * pt[1];
523 const double z_factor =
524 PipeSegment::compute_z_expansion(x_new, y_new, data);
526 ExcMessage(
"Invalid input: at least one pipe segment "
527 "is not long enough in this configuration"));
528 const double z_new = z_factor * pt[2];
530 return {x_new, y_new, z_new};
542 const double rotation_angle =
544 const vector3d rotation_axis = [&]() {
545 const vector3d rotation_axis =
546 cross_product_3d(directions[2], skeleton_unit[p]);
547 const double norm = rotation_axis.norm();
548 if (norm < tolerance)
549 return directions[1];
551 return rotation_axis /
norm;
555 rotation_axis, rotation_angle);
572 const vector3d Rx = rotation_matrix * directions[0];
576 const vector3d normal_projected_on_opening =
577 normal - (normal * skeleton_unit[p]) * skeleton_unit[p];
587 const double lateral_angle =
589 normal_projected_on_opening,
603 manifolds.emplace_back(
604 normal_projected_on_opening /
605 normal_projected_on_opening.norm(),
612 pipe,
tria,
tria, tolerance_length,
true);
615 for (
unsigned int p = 0; p < n_pipes; ++p)
622 for (
const auto &cell :
tria.active_cell_iterators())
623 for (const auto &face : cell->face_iterators())
624 if (face->at_boundary())
627 face->manifold_id() == n_pipes)
629 face->set_boundary_id(cell->material_id());
632 face->set_boundary_id(n_pipes);
642#include "grid_generator_pipe_junction.inst"
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const =0
numbers::NumberTraits< Number >::real_type norm() const
unsigned int n_cells() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
void hyper_ball_balanced(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.)
void pipe_junction(Triangulation< dim, spacedim > &tria, const std::vector< std::pair< Point< spacedim >, double > > &openings, const std::pair< Point< spacedim >, double > &bifurcation, const double aspect_ratio=0.5)
void extrude_triangulation(const Triangulation< 2, 2 > &input, const unsigned int n_slices, const double height, Triangulation< 3, 3 > &result, const bool copy_manifold_ids=false, const std::vector< types::manifold_id > &manifold_priorities={})
void merge_triangulations(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result, const double duplicated_vertex_tolerance=1.0e-12, const bool copy_manifold_ids=false, const bool copy_boundary_ids=false)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Number signed_angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b, const Tensor< 1, spacedim, Number > &axis)
Number angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b)
static constexpr double PI
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::Triangulation< dim, spacedim > & tria