41 double skeleton_length;
43 double cosecant_polar;
44 double cotangent_polar;
45 double cotangent_azimuth_half_right;
46 double cotangent_azimuth_half_left;
56 compute_z_expansion(
const double x,
58 const AdditionalData &data)
65 + x * data.cotangent_polar
67 -
std::abs(y) * data.cosecant_polar *
68 ((y > 0) ? data.cotangent_azimuth_half_right :
69 data.cotangent_azimuth_half_left);
80 template <
int dim,
int spacedim = dim>
106 const AdditionalData &data,
107 const double tolerance = 1e-10);
112 virtual std::unique_ptr<::Manifold<dim, spacedim>>
113 clone()
const override;
132 push_forward(
const Point<3> &chart_point)
const override;
153 const AdditionalData data;
158 const double tolerance;
169 template <
int dim,
int spacedim>
170 Manifold<dim, spacedim>::Manifold(
174 const AdditionalData &data,
175 const double tolerance)
177 , normal_direction(normal_direction)
178 , direction(direction)
179 , point_on_axis(point_on_axis)
181 , tolerance(tolerance)
182 , dxn(cross_product_3d(direction, normal_direction))
186 "PipeSegment::Manifold can only be used for spacedim==3!"));
189 ExcMessage(
"Normal direction must be unit vector."));
191 ExcMessage(
"Direction must be unit vector."));
192 Assert(normal_direction * direction < tolerance,
194 "Direction and normal direction must be perpendicular."));
199 template <
int dim,
int spacedim>
200 std::unique_ptr<::Manifold<dim, spacedim>>
203 return std::make_unique<Manifold<dim, spacedim>>(*this);
208 template <
int dim,
int spacedim>
210 Manifold<dim, spacedim>::pull_back(
const Point<spacedim> &space_point)
const
214 double lambda = normalized_point * direction;
217 const double r = p_diff.
norm();
219 Assert(r > tolerance * data.skeleton_length,
221 "This class won't handle points on the direction axis."));
234 return {r, phi, lambda};
239 template <
int dim,
int spacedim>
241 Manifold<dim, spacedim>::push_forward(
const Point<3> &chart_point)
const
244 const double sine_r = chart_point[0] *
std::sin(chart_point[1]);
245 const double cosine_r = chart_point[0] *
std::cos(chart_point[1]);
248 normal_direction * cosine_r + dxn * sine_r;
251 const double lambda =
252 chart_point[2] * compute_z_expansion(cosine_r, sine_r, data);
255 return point_on_axis + direction * lambda + intermediate;
264 template <
int dim,
int spacedim>
282 const std::vector<std::pair<
Point<3>,
double>> &openings,
283 const std::pair<
Point<3>,
double> &bifurcation,
284 const double aspect_ratio)
286 constexpr unsigned int dim = 3;
287 constexpr unsigned int spacedim = 3;
290 constexpr unsigned int n_pipes = 3;
291 constexpr double tolerance = 1.e-12;
295 Assert(bifurcation.second > 0,
296 ExcMessage(
"Invalid input: negative radius."));
297 Assert(openings.size() == n_pipes,
298 ExcMessage(
"Invalid input: only 3 openings allowed."));
299 for (
const auto &opening : openings)
306 const auto cyclic = [](
const unsigned int i) ->
unsigned int {
307 constexpr unsigned int n_pipes = 3;
308 return (i < (n_pipes - 1)) ? i + 1 : 0;
310 const auto anticyclic = [](
const unsigned int i) ->
unsigned int {
311 constexpr unsigned int n_pipes = 3;
312 return (i > 0) ? i - 1 : n_pipes - 1;
316 const std::array<vector3d, spacedim> directions = {
317 {vector3d({1., 0., 0.}), vector3d({0., 1., 0.}), vector3d({0., 0., 1.})}};
323 std::array<vector3d, n_pipes> skeleton;
324 for (
unsigned int p = 0; p < n_pipes; ++p)
325 skeleton[p] = bifurcation.first - openings[p].first;
327 std::array<double, n_pipes> skeleton_length;
328 for (
unsigned int p = 0; p < n_pipes; ++p)
329 skeleton_length[p] = skeleton[p].
norm();
335 const double tolerance_length =
337 *std::max_element(skeleton_length.begin(), skeleton_length.end());
339 std::array<vector3d, n_pipes> skeleton_unit;
340 for (
unsigned int p = 0; p < n_pipes; ++p)
342 Assert(skeleton_length[p] > tolerance_length,
343 ExcMessage(
"Invalid input: bifurcation matches opening."));
344 skeleton_unit[p] = skeleton[p] / skeleton_length[p];
356 vector3d normal = cross_product_3d(skeleton_unit[1] - skeleton_unit[0],
357 skeleton_unit[2] - skeleton_unit[0]);
358 Assert(normal.norm() > tolerance_length,
359 ExcMessage(
"Invalid input: all three openings "
360 "are located on one line."));
361 normal /= normal.norm();
365 std::array<vector3d, n_pipes> skeleton_plane;
366 for (
unsigned int p = 0; p < n_pipes; ++p)
368 skeleton_plane[p] = skeleton[p] - (skeleton[p] * normal) * normal;
370 tolerance * skeleton_plane[p].
norm(),
372 Assert(skeleton_plane[p].
norm() > tolerance_length,
388 ExcMessage(
"The output triangulation object needs to be empty."));
390 std::vector<PipeSegment::Manifold<dim, spacedim>> manifolds;
391 for (
unsigned int p = 0; p < n_pipes; ++p)
404 const unsigned int n_slices =
405 1 +
static_cast<unsigned int>(std::ceil(
406 aspect_ratio * skeleton_length[p] /
417 for (
const auto &cell : pipe.active_cell_iterators())
419 cell->set_material_id(p);
421 for (
const auto &face : cell->face_iterators())
422 if (face->at_boundary())
424 const auto center_z = face->center()[2];
430 else if (
std::abs(center_z - 1.) < tolerance)
437 cell->set_all_manifold_ids(n_pipes);
438 face->set_all_manifold_ids(p);
473 const double polar_angle =
482 const double azimuth_angle_right =
484 skeleton_plane[cyclic(p)],
487 ExcMessage(
"Invalid input: at least two openings located "
488 "in same direction from bifurcation"));
490 const double azimuth_angle_left =
492 skeleton_plane[anticyclic(p)],
495 ExcMessage(
"Invalid input: at least two openings located "
496 "in same direction from bifurcation"));
500 PipeSegment::AdditionalData data;
501 data.skeleton_length = skeleton_length[p];
502 data.cosecant_polar = 1. /
std::sin(polar_angle);
503 data.cotangent_polar =
std::cos(polar_angle) * data.cosecant_polar;
504 data.cotangent_azimuth_half_right =
std::cos(.5 * azimuth_angle_right) /
506 data.cotangent_azimuth_half_left =
510 const auto pipe_segment_transform =
514 const double r_factor =
515 (bifurcation.second - openings[p].second) * pt[2] +
517 const double x_new = r_factor * pt[0];
518 const double y_new = r_factor * pt[1];
522 const double z_factor =
523 PipeSegment::compute_z_expansion(x_new, y_new, data);
525 ExcMessage(
"Invalid input: at least one pipe segment "
526 "is not long enough in this configuration"));
527 const double z_new = z_factor * pt[2];
529 return {x_new, y_new, z_new};
541 const double rotation_angle =
543 const vector3d rotation_axis = [&]() {
544 const vector3d rotation_axis =
545 cross_product_3d(directions[2], skeleton_unit[p]);
546 const double norm = rotation_axis.norm();
547 if (norm < tolerance)
548 return directions[1];
550 return rotation_axis /
norm;
554 rotation_axis, rotation_angle);
571 const vector3d Rx = rotation_matrix * directions[0];
575 const vector3d normal_projected_on_opening =
576 normal - (normal * skeleton_unit[p]) * skeleton_unit[p];
586 const double lateral_angle =
588 normal_projected_on_opening,
602 manifolds.emplace_back(
603 normal_projected_on_opening /
604 normal_projected_on_opening.norm(),
611 pipe, tria, tria, tolerance_length,
true);
614 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
#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)
#define DEAL_II_NOT_IMPLEMENTED()
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 > &)