]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Generalized tria extruding and reimplementing existing extrude function #6158 6161/head
authorWeixiong Zheng <weixiong.zheng@berkeley.edu>
Thu, 5 Apr 2018 01:17:34 +0000 (18:17 -0700)
committerWeixiong Zheng <weixiong.zheng@berkeley.edu>
Thu, 5 Apr 2018 11:01:05 +0000 (04:01 -0700)
1. Added extrude_triangulation overload with doxygen documentation
   as argument-generalized version to the existing function with testing
2. Re-implemented existing extruding by using the newly developed one
3. A new entry was added for doc/news/changes/minor/

doc/news/changes/minor/20180405WeixiongZheng [new file with mode: 0644]
include/deal.II/grid/grid_generator.h
source/grid/grid_generator.cc
tests/grid/extrude_orientation_02.cc
tests/grid/extrude_orientation_03.cc [new file with mode: 0644]
tests/grid/extrude_orientation_03.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20180405WeixiongZheng b/doc/news/changes/minor/20180405WeixiongZheng
new file mode 100644 (file)
index 0000000..02d9d14
--- /dev/null
@@ -0,0 +1,5 @@
+New: Extend GridGenerator::extrude_triangulation to using exact slicing z-axis
+values and reimplement the exsiting GridGenerator::extrude_triangulation using
+the newly developed overload.
+<br>
+(Weixiong Zheng, 2018/04/05)
\ No newline at end of file
index 175cdf55b06ee1b04f68078dcecd8be5fb92810b..df6f2ec0788a405dc5414b1ae110708650e398fa 100644 (file)
@@ -1038,7 +1038,6 @@ namespace GridGenerator
                                            const std::set<typename Triangulation<dim, spacedim>::active_cell_iterator> &cells_to_remove,
                                            Triangulation<dim, spacedim>       &result);
 
-
   /**
    * Take a 2d Triangulation that is being extruded in z direction by the
    * total height of @p height using @p n_slices slices (minimum is 2). The
@@ -1055,6 +1054,26 @@ namespace GridGenerator
                          const double               height,
                          Triangulation<3,3>        &result);
 
+  /**
+   * Overload of the previous function. Take a 2d Triangulation that is being
+   * extruded. Differing from the previous function taking height and number of
+   * slices for uniform extrusion, this function takes z-axis values
+   * @p slice_coordinates where the slicing will happen. The boundary indicators
+   * of the faces of @p input are going to be assigned to the corresponding side
+   * walls in z direction. The bottom and top get the next two free boundary
+   * indicators.
+   *
+   * @note The 2d input triangulation @p input must be a coarse mesh that has
+   * no refined cells.
+   *
+   * @author Weixiong Zheng, 2018
+   */
+  void
+  extrude_triangulation (const Triangulation<2, 2>             &input,
+                         const std::vector<double> &slice_coordinates,
+                         Triangulation<3,3>                   &result);
+
+
   /**
    * Given an input triangulation @p in_tria, this function makes a new flat
    * triangulation @p out_tria which contains a single level with all active
index aef2264ec975eafd4f1ab0dcc2cbdbc2b52f63dd..cd44e3b4c9900bf9700a75c1cfab527da2da3c73 100644 (file)
@@ -4044,12 +4044,34 @@ namespace GridGenerator
     Assert(n_slices>=2,
            ExcMessage("The number of slices for extrusion must be at least 2."));
 
+    const double delta_h = height / (n_slices-1);
+    std::vector<double> slices_z_values;
+    for (unsigned int i=0; i<n_slices; ++i)
+      slices_z_values.push_back(i*delta_h);
+    extrude_triangulation(input, slices_z_values, result);
+  }
+
+  void
+  extrude_triangulation (const Triangulation<2,2> &input,
+                         const std::vector<double> &slice_coordinates,
+                         Triangulation<3,3> &result)
+  {
+    Assert (input.n_levels() == 1,
+            ExcMessage ("The input triangulation must be a coarse mesh, i.e., it must "
+                        "not have been refined."));
+    Assert(result.n_cells()==0,
+           ExcMessage("The output triangulation object needs to be empty."));
+    Assert(slice_coordinates.size()>=2,
+           ExcMessage("The number of slices for extrusion must be at least 2."));
+    const unsigned int n_slices = slice_coordinates.size();
+    Assert(std::is_sorted(slice_coordinates.begin(), slice_coordinates.end()),
+           ExcMessage("Slice z-coordinates should be in ascending order"));
     std::vector<Point<3> > points(n_slices*input.n_vertices());
     std::vector<CellData<3> > cells;
     cells.reserve((n_slices-1)*input.n_active_cells());
 
     // copy the array of points as many times as there will be slices,
-    // one slice at a time
+    // one slice at a time. The z-axis value are defined in slices_coordinates
     for (unsigned int slice=0; slice<n_slices; ++slice)
       {
         for (unsigned int i=0; i<input.n_vertices(); ++i)
@@ -4057,7 +4079,7 @@ namespace GridGenerator
             const Point<2> &v = input.get_vertices()[i];
             points[slice*input.n_vertices()+i](0) = v(0);
             points[slice*input.n_vertices()+i](1) = v(1);
-            points[slice*input.n_vertices()+i](2) = height * slice / (n_slices-1);
+            points[slice*input.n_vertices()+i](2) = slice_coordinates[slice];
           }
       }
 
@@ -4150,6 +4172,7 @@ namespace GridGenerator
     result.create_triangulation (points,
                                  cells,
                                  s);
+
   }
 
 
index d20cfad2c48a53925a3ed717cdd461426644cc5d..4877aa83556935827c8173eba45566b6049a84d8 100644 (file)
 // test this for a circle extruded to a cylinder
 
 #include "../tests.h"
-#include <deal.II/base/tensor.h>
 #include <deal.II/grid/tria.h>
 #include <deal.II/grid/grid_generator.h>
-#include <deal.II/grid/grid_out.h>
-#include <deal.II/fe/fe_q.h>
 
-
-
-void test(std::ostream &out)
+void test()
 {
   Triangulation<2> tr;
   GridGenerator::hyper_ball (tr);
@@ -71,5 +66,5 @@ int main()
 {
   initlog();
 
-  test(deallog.get_file_stream());
+  test();
 }
diff --git a/tests/grid/extrude_orientation_03.cc b/tests/grid/extrude_orientation_03.cc
new file mode 100644 (file)
index 0000000..306ad1d
--- /dev/null
@@ -0,0 +1,69 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// Test GridGenerator::extrude_triangulation taking slice z-coordinate values.
+// This test is just a replicate of the one in extrude_orientation_02.cc except
+// that the newly created overload of GridGenerator::extrude_triangulation is
+// tested for face and edge orientations for the resulting 3d mesh.
+// See https://github.com/dealii/dealii/issues/6158.
+//
+// test this for a circle extruded to a cylinder
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+
+void test()
+{
+  Triangulation<2> tr;
+  GridGenerator::hyper_ball (tr);
+
+  for (Triangulation<2>::active_cell_iterator c=tr.begin_active();
+       c!=tr.end(); ++c)
+    {
+      deallog << "2d cell " << c << " has the following face orientations:"
+              << std::endl;
+      for (unsigned int l=0; l<GeometryInfo<2>::faces_per_cell; ++l)
+        deallog << "    " << (c->face_orientation(l) ? "true" : "false")
+                << std::endl;
+    }
+
+  Triangulation<3> tr3;
+  std::vector<double> slice_points = {0, 0.1, 0.5};
+  GridGenerator::extrude_triangulation(tr, slice_points, tr3);
+
+  for (Triangulation<3>::active_cell_iterator c=tr3.begin_active();
+       c!=tr3.end(); ++c)
+    {
+      deallog << "3d cell " << c << " has the following face orientation/flips and edge orientations:"
+              << std::endl;
+      for (unsigned int f=0; f<GeometryInfo<3>::faces_per_cell; ++f)
+        deallog << "    face=" << f
+                << (c->face_orientation(f) ? " -> true" : " -> false")
+                << (c->face_flip(f) ? "/true" : "/false")
+                << std::endl;
+      for (unsigned int e=0; e<GeometryInfo<3>::lines_per_cell; ++e)
+        deallog << "    edge=" << e << (c->line_orientation(e) ? " -> true" : " -> false")
+                << std::endl;
+    }
+}
+
+
+int main()
+{
+  initlog();
+
+  test();
+}
\ No newline at end of file
diff --git a/tests/grid/extrude_orientation_03.output b/tests/grid/extrude_orientation_03.output
new file mode 100644 (file)
index 0000000..040a334
--- /dev/null
@@ -0,0 +1,216 @@
+
+DEAL::2d cell 0.0 has the following face orientations:
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::2d cell 0.1 has the following face orientations:
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::2d cell 0.2 has the following face orientations:
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::2d cell 0.3 has the following face orientations:
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::2d cell 0.4 has the following face orientations:
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::    true
+DEAL::3d cell 0.0 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.1 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.2 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> false/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.3 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> false/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.4 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.5 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> true/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.6 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> false/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.7 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> true/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> false/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.8 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> false/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> false/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true
+DEAL::3d cell 0.9 has the following face orientation/flips and edge orientations:
+DEAL::    face=0 -> true/false
+DEAL::    face=1 -> false/false
+DEAL::    face=2 -> true/false
+DEAL::    face=3 -> false/false
+DEAL::    face=4 -> true/false
+DEAL::    face=5 -> true/false
+DEAL::    edge=0 -> true
+DEAL::    edge=1 -> true
+DEAL::    edge=2 -> true
+DEAL::    edge=3 -> true
+DEAL::    edge=4 -> true
+DEAL::    edge=5 -> true
+DEAL::    edge=6 -> true
+DEAL::    edge=7 -> true
+DEAL::    edge=8 -> true
+DEAL::    edge=9 -> true
+DEAL::    edge=10 -> true
+DEAL::    edge=11 -> true

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.