]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added the colorize option to subdivided_hyper_cube grid generator 9142/head
authorBruno Blais <blais.bruno@gmail.com>
Thu, 5 Dec 2019 14:45:58 +0000 (09:45 -0500)
committerBruno Blais <blais.bruno@gmail.com>
Thu, 5 Dec 2019 23:40:32 +0000 (18:40 -0500)
The subdivided_hyper_cube is the one of the only grid generator which
does not allow the boundaries to be colorized. This is a non-standard
behavior, because all of the other generators (hyper_cube,
hyper_rectangle, subdivided_hyper_rectangle) allow one to colorize the
boundaries. This small pull request addresses this by adding the
colorize option to the subdivided_hyper_cube grid. It also fixes the
grid_generator_from_name to make it coherent with the new
implementaiton.

doc/news/changes/minor/20191205BrunoBlais [new file with mode: 0644]
include/deal.II/grid/grid_generator.h
source/grid/grid_generator.cc
source/grid/grid_generator.inst.in
source/grid/grid_generator_from_name.cc
tests/grid/grid_generator_subdivided_hyper_cube_01.cc [new file with mode: 0644]
tests/grid/grid_generator_subdivided_hyper_cube_01.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20191205BrunoBlais b/doc/news/changes/minor/20191205BrunoBlais
new file mode 100644 (file)
index 0000000..72f61b3
--- /dev/null
@@ -0,0 +1,3 @@
+Improved: GridGenerator::subdivided_hyper_cube() gained an option to also set boundary ids.
+<br>
+(Bruno Blais, 2019/12/05)
index 616b3eaa4888d201d86249aa29eb70996facf935..4a50c738bc20b828d3bd111e506daa6f352baf51 100644 (file)
@@ -147,13 +147,16 @@ namespace GridGenerator
    * @param left Lower bound for the interval used to create the hyper cube.
    *
    * @param right Upper bound for the interval used to create the hyper cube.
+   *
+   * @param colorize Assign different boundary ids if set to true.
    */
   template <int dim, int spacedim>
   void
   subdivided_hyper_cube(Triangulation<dim, spacedim> &tria,
                         const unsigned int            repetitions,
-                        const double                  left  = 0.,
-                        const double                  right = 1.);
+                        const double                  left     = 0.,
+                        const double                  right    = 1.,
+                        const bool                    colorize = false);
 
   /**
    * Create a coordinate-parallel brick from the two diagonally opposite
index 4309dfb59bbd21bbe50997942ffa306538d6c5ba..dc938b892de767b77da8096a7200ddc6434ed444 100644 (file)
@@ -1234,7 +1234,8 @@ namespace GridGenerator
   subdivided_hyper_cube(Triangulation<dim, spacedim> &tria,
                         const unsigned int            repetitions,
                         const double                  left,
-                        const double                  right)
+                        const double                  right,
+                        const bool                    colorize)
   {
     Assert(repetitions >= 1, ExcInvalidRepetitions(repetitions));
     Assert(left < right,
@@ -1248,7 +1249,7 @@ namespace GridGenerator
       }
 
     std::vector<unsigned int> reps(dim, repetitions);
-    subdivided_hyper_rectangle(tria, reps, p0, p1);
+    subdivided_hyper_rectangle(tria, reps, p0, p1, colorize);
   }
 
 
index 8b0d1ebdc39b31261dda6556bb8d871dd4fd787f..76522efd9f1698d0be5cba6a740a8ce828480002 100644 (file)
@@ -38,7 +38,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
         Triangulation<deal_II_dimension, deal_II_space_dimension> &,
         const unsigned int,
         const double,
-        const double);
+        const double,
+        const bool);
 
 
       template void
index 321824e26887701619ffaafb20cd36f0e8dccfdc..7cf00ec86a1572f810ec87ce8d59490689533055 100644 (file)
@@ -322,7 +322,7 @@ namespace GridGenerator
                                                             arguments,
                                                             tria);
     else if (name == "subdivided_hyper_cube")
-      parse_and_create<dim, spacedim, unsigned int, double, double>(
+      parse_and_create<dim, spacedim, unsigned int, double, double, bool>(
         subdivided_hyper_cube, arguments, tria);
     else if (name == "hyper_rectangle")
       parse_and_create<dim,
diff --git a/tests/grid/grid_generator_subdivided_hyper_cube_01.cc b/tests/grid/grid_generator_subdivided_hyper_cube_01.cc
new file mode 100644 (file)
index 0000000..11d8714
--- /dev/null
@@ -0,0 +1,53 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Test for GridGenerator::subdivided_hyper_cube() with colorized boundary
+// conditions
+
+#include <deal.II/base/tensor.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+template <int dim>
+void
+test()
+{
+  deallog << "dim = " << dim << std::endl;
+  Triangulation<dim> tr;
+  GridGenerator::subdivided_hyper_cube(tr, 1, -1, 1, true);
+  for (const auto &cell : tr.active_cell_iterators())
+    {
+      deallog << "cell:" << std::endl;
+
+      for (const auto &face : cell->face_iterators())
+        {
+          if (face->at_boundary())
+            deallog << "boundary id = " << face->boundary_id()
+                    << " center = " << face->center()
+                    << " faceidx = " << face->index() << std::endl;
+        }
+    }
+}
+
+int
+main()
+{
+  initlog();
+  test<2>();
+  test<3>();
+}
diff --git a/tests/grid/grid_generator_subdivided_hyper_cube_01.output b/tests/grid/grid_generator_subdivided_hyper_cube_01.output
new file mode 100644 (file)
index 0000000..49679ef
--- /dev/null
@@ -0,0 +1,15 @@
+
+DEAL::dim = 2
+DEAL::cell:
+DEAL::boundary id = 0 center = -1.00000 0.00000 faceidx = 1
+DEAL::boundary id = 1 center = 1.00000 0.00000 faceidx = 2
+DEAL::boundary id = 2 center = 0.00000 -1.00000 faceidx = 0
+DEAL::boundary id = 3 center = 0.00000 1.00000 faceidx = 3
+DEAL::dim = 3
+DEAL::cell:
+DEAL::boundary id = 0 center = -1.00000 0.00000 0.00000 faceidx = 2
+DEAL::boundary id = 1 center = 1.00000 0.00000 0.00000 faceidx = 3
+DEAL::boundary id = 2 center = 0.00000 -1.00000 0.00000 faceidx = 0
+DEAL::boundary id = 3 center = 0.00000 1.00000 0.00000 faceidx = 4
+DEAL::boundary id = 4 center = 0.00000 0.00000 -1.00000 faceidx = 1
+DEAL::boundary id = 5 center = 0.00000 0.00000 1.00000 faceidx = 5

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.