]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Minor adjustments to the pipe junction geometry.
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 27 Jan 2022 22:07:23 +0000 (15:07 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 28 Jan 2022 01:46:34 +0000 (18:46 -0700)
include/deal.II/grid/grid_generator.h
source/grid/grid_generator_pipe_junction.cc

index 66b9e6e3d864b00a08be69a479f5d8f4fc90343b..e42073c3c7007872815dc16fbe3e05d27cb6171f 100644 (file)
@@ -1050,7 +1050,7 @@ namespace GridGenerator
    * opening to match their index. All other boundary faces will be assigned
    * boundary ID 3.
    *
-   * <em>Manifold IDs</em> will be set on the mantles of each truncated cone in
+   * @ref GlossManifoldIndicator "Manifold IDs" will be set on the mantles of each truncated cone in
    * the same way. Each cone will have a special manifold object assigned, which
    * is based on the CylindricalManifold class. Further, all cells adjacent to
    * the mantle are given the manifold ID 3. If desired, you can assign an
@@ -1075,8 +1075,9 @@ namespace GridGenerator
    * @param aspect_ratio Aspect ratio of cells, specified as radial over z-extension.
    *                     Default ratio is $\Delta r/\Delta z = 1/2$.
    *
-   * Common configurations of tee fittings that can be generated with the
-   * following sets of parameters are:
+   * Common configurations of tee fittings (that is, "T" fittings, mimicking the
+   * geometry of the letter "T") can be generated with the
+   * following sets of parameters:
    * <div class="threecolumn" style="width: 80%; text-align: center;">
    *   <div>
    *     \htmlonly <style>div.image
index e2f1d675da29a14251eb9f4fe195545f79b140d5..47c6805dad892ed09d71c0ee49a3424a31215d2c 100644 (file)
@@ -288,7 +288,7 @@ namespace GridGenerator
   {
     constexpr unsigned int dim      = 3;
     constexpr unsigned int spacedim = 3;
-    using vector                    = Tensor<1, spacedim, double>;
+    using vector3d                  = Tensor<1, spacedim, double>;
 
     constexpr unsigned int n_pipes   = 3;
     constexpr double       tolerance = 1.e-12;
@@ -316,14 +316,14 @@ namespace GridGenerator
     };
 
     // Cartesian base represented by unit vectors.
-    constexpr std::array<vector, spacedim> directions = {
-      {vector({1., 0., 0.}), vector({0., 1., 0.}), vector({0., 0., 1.})}};
+    constexpr std::array<vector3d, spacedim> directions = {
+      {vector3d({1., 0., 0.}), vector3d({0., 1., 0.}), vector3d({0., 0., 1.})}};
 
     // The skeleton corresponds to the axis of symmetry in the center of each
     // pipe segment. Each skeleton vector points from the associated opening to
     // the common bifurcation point. For convenience, we also compute length and
     // unit vector of every skeleton vector here.
-    std::array<vector, n_pipes> skeleton;
+    std::array<vector3d, n_pipes> skeleton;
     for (unsigned int p = 0; p < n_pipes; ++p)
       skeleton[p] = bifurcation.first - openings[p].first;
 
@@ -339,7 +339,7 @@ namespace GridGenerator
       tolerance *
       *std::max_element(skeleton_length.begin(), skeleton_length.end());
 
-    std::array<vector, n_pipes> skeleton_unit;
+    std::array<vector3d, n_pipes> skeleton_unit;
     for (unsigned int p = 0; p < n_pipes; ++p)
       {
         Assert(skeleton_length[p] > tolerance_length,
@@ -356,8 +356,8 @@ namespace GridGenerator
     // which all pipe segments meet. If we would interpret the bifurcation as a
     // ball joint, the normal vector would correspond to the polar axis of the
     // ball.
-    vector normal = cross_product_3d(skeleton_unit[1] - skeleton_unit[0],
-                                     skeleton_unit[2] - skeleton_unit[0]);
+    vector3d normal = cross_product_3d(skeleton_unit[1] - skeleton_unit[0],
+                                       skeleton_unit[2] - skeleton_unit[0]);
     Assert(normal.norm() > tolerance_length,
            ExcMessage("Invalid input: all three openings "
                       "are located on one line."));
@@ -365,7 +365,7 @@ namespace GridGenerator
 
     // Projections of all skeleton vectors perpendicular to the normal vector,
     // or in other words, onto the plane described above.
-    std::array<vector, n_pipes> skeleton_plane;
+    std::array<vector3d, n_pipes> skeleton_plane;
     for (unsigned int p = 0; p < n_pipes; ++p)
       {
         skeleton_plane[p] = skeleton[p] - (skeleton[p] * normal) * normal;
@@ -447,7 +447,7 @@ namespace GridGenerator
         // Step 2: transform unit cylinder to pipe segment
         //
         // For the given cylinder, we will interpret the base in the xy-plane as
-        // the cross section of the opening, and the base at z=1 as the surface
+        // the cross section of the opening, and the top at z=1 as the surface
         // where all pipe segments meet. On the latter surface, we assign the
         // section in positive y-direction to face the next (right/cyclic) pipe
         // segment, and allocate the domain in negative y-direction to border
@@ -510,7 +510,8 @@ namespace GridGenerator
           std::cos(.5 * azimuth_angle_left) / std::sin(.5 * azimuth_angle_left);
 
         // Now transform the cylinder as described above.
-        const auto pipe_segment = [&](const Point<spacedim> &pt) {
+        const auto pipe_segment_transform =
+          [&](const Point<spacedim> &pt) -> Point<spacedim> {
           // We transform the cylinder in x- and y-direction to become a
           // truncated cone, similarly to GridGenerator::truncated_cone().
           const double r_factor =
@@ -528,15 +529,9 @@ namespace GridGenerator
                             "is not long enough in this configuration"));
           const double z_new = z_factor * pt[2];
 
-          // TODO: MSVC can't capture const or constexpr values in lambda
-          // functions (due to either a missing implementation or a bug).
-          // Instead, we duplicate the declaration here.
-          //   See also: https://developercommunity.visualstudio.com/t/
-          //             invalid-template-argument-expected-compile-time-co/187862
-          constexpr unsigned int spacedim = 3;
-          return Point<spacedim>(x_new, y_new, z_new);
+          return {x_new, y_new, z_new};
         };
-        GridTools::transform(pipe_segment, pipe);
+        GridTools::transform(pipe_segment_transform, pipe);
 
         //
         // Step 3: rotate pipe segment to match skeleton direction
@@ -548,8 +543,8 @@ namespace GridGenerator
         // cross product of both vectors.
         const double rotation_angle =
           Physics::VectorRelations::angle(directions[2], skeleton_unit[p]);
-        const vector rotation_axis = [&]() {
-          const vector rotation_axis =
+        const vector3d rotation_axis = [&]() {
+          const vector3d rotation_axis =
             cross_product_3d(directions[2], skeleton_unit[p]);
           const double norm = rotation_axis.norm();
           if (norm < tolerance)
@@ -576,11 +571,11 @@ namespace GridGenerator
         // With the latest rotation however, this is no longer the case. We
         // rotate the unit vector in x-direction in the same fashion, which
         // gives us the current direction of the projected edge.
-        const vector Rx = rotation_matrix * directions[0];
+        const vector3d Rx = rotation_matrix * directions[0];
 
         // To determine how far we need to rotate, we also need to project the
         // polar axis of the bifurcation ball joint into the same plane.
-        const vector normal_projected_on_opening =
+        const vector3d normal_projected_on_opening =
           normal - (normal * skeleton_unit[p]) * skeleton_unit[p];
 
         // Both the projected normal and Rx must be in the opening plane.
@@ -598,7 +593,8 @@ namespace GridGenerator
         GridTools::rotate(skeleton_unit[p], lateral_angle, pipe);
 
         //
-        // Step 5: shift to final position
+        // Step 5: shift to final position and merge this pipe into the entire
+        // assembly
         //
         GridTools::shift(openings[p].first, pipe);
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.