]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Also colorize the edges, not just the faces, of a hypershell. 451/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 17 Jan 2015 16:33:04 +0000 (10:33 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 17 Jan 2015 20:32:12 +0000 (14:32 -0600)
doc/news/changes.h
include/deal.II/grid/grid_generator.h
source/grid/grid_generator.cc
tests/deal.II/grid_hyper_shell_06.cc [new file with mode: 0644]
tests/deal.II/grid_hyper_shell_06.output [new file with mode: 0644]

index 3a30356d57bab1bacd12c95c6dc624d8ca094afa..c6ce5b1279d1f58c75aa565adf75486ec6aabdac 100644 (file)
@@ -213,6 +213,14 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> Fixed: In 3d, when you set the <code>colorize</code> flag of
+  GridGenerator::hyper_shell(), the faces of the domain were colored but
+  the edges were not. This was an oversight because to refine correctly,
+  the edges also have to have the appropriate boundary indicator set.
+  <br>
+  (Wolfgang Bangerth, 2015/01/16)
+  </li>
+
   <li> New: dealii::multithread_info.n_cpus returns the correct number of CPU 
   on FreeBSD.
   <br>
index 8a4b2528d4118f20583852a5982525ea970f7bec..8cc20b629d81207d37a27a5db12b62077e1a2813 100644 (file)
@@ -428,9 +428,10 @@ namespace GridGenerator
    * of cells of the resulting triangulation, i.e., how many cells form the
    * ring (in 2d) or the shell (in 3d).
    *
-   * If the flag @p colorize is @p true, then the outer boundary will have the
-   * indicator 1, while the inner boundary has id zero. If the flag is @p
-   * false, both have indicator zero.
+   * If the flag @p colorize is @p true, then the outer boundary will
+   * have the indicator 1, while the inner boundary has id zero. In
+   * 3d, this applies to both the faces and the edges of these
+   * boundaries. If the flag is @p false, both have indicator zero.
    *
    * In 2d, the number <tt>n_cells</tt> of elements for this initial
    * triangulation can be chosen arbitrarily. If the number of initial cells
index caf976491b410ed28d2b636a9d8fef37db008bf8..845ffb15463cfc4a53f3774c659ddc9fc86669f0 100644 (file)
@@ -173,7 +173,7 @@ namespace GridGenerator
            cell != tria.end (); ++cell)
         {
           Assert(cell->face(2)->at_boundary(), ExcInternalError());
-          cell->face (2)->set_boundary_indicator (1);
+          cell->face (2)->set_all_boundary_indicators (1);
         }
     }
 
@@ -196,23 +196,28 @@ namespace GridGenerator
         {
           Triangulation<3>::cell_iterator cell = tria.begin();
 
-          cell->face(4)->set_boundary_indicator(1);
           Assert (cell->face(4)->at_boundary(), ExcInternalError());
+          cell->face(4)->set_all_boundary_indicators(1);
 
-          (++cell)->face(2)->set_boundary_indicator(1);
+          ++cell;
           Assert (cell->face(2)->at_boundary(), ExcInternalError());
+          cell->face(2)->set_all_boundary_indicators(1);
 
-          (++cell)->face(2)->set_boundary_indicator(1);
+          ++cell;
           Assert (cell->face(2)->at_boundary(), ExcInternalError());
+          cell->face(2)->set_all_boundary_indicators(1);
 
-          (++cell)->face(0)->set_boundary_indicator(1);
+          ++cell;
           Assert (cell->face(0)->at_boundary(), ExcInternalError());
+          cell->face(0)->set_all_boundary_indicators(1);
 
-          (++cell)->face(2)->set_boundary_indicator(1);
+          ++cell;
           Assert (cell->face(2)->at_boundary(), ExcInternalError());
+          cell->face(2)->set_all_boundary_indicators(1);
 
-          (++cell)->face(0)->set_boundary_indicator(1);
+          ++cell;
           Assert (cell->face(0)->at_boundary(), ExcInternalError());
+          cell->face(0)->set_all_boundary_indicators(1);
         }
       else if (tria.n_cells() == 12)
         {
@@ -222,7 +227,7 @@ namespace GridGenerator
                cell != tria.end(); ++cell)
             {
               Assert (cell->face(5)->at_boundary(), ExcInternalError());
-              cell->face(5)->set_boundary_indicator(1);
+              cell->face(5)->set_all_boundary_indicators(1);
             }
         }
       else if (tria.n_cells() == 96)
@@ -243,7 +248,7 @@ namespace GridGenerator
                cell != tria.end(); ++cell)
             if (cell->face(5)->at_boundary())
               {
-                cell->face(5)->set_boundary_indicator(1);
+                cell->face(5)->set_all_boundary_indicators(1);
                 ++count;
               }
           Assert (count == 48, ExcInternalError());
diff --git a/tests/deal.II/grid_hyper_shell_06.cc b/tests/deal.II/grid_hyper_shell_06.cc
new file mode 100644 (file)
index 0000000..fd1463b
--- /dev/null
@@ -0,0 +1,74 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 2015 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.
+//
+// ---------------------------------------------------------------------
+
+
+// GridGenerator::hyper_shell colorized the faces but forgot the
+// edges. This is not useful because the colorization is usually done
+// so that one can attach a boundary or manifold object to these parts
+// of the boundary
+
+#include "../tests.h"
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/base/logstream.h>
+#include <cmath>
+#include <cstdlib>
+
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+
+std::ofstream logfile("output");
+
+
+template<int dim>
+void check (double r1, double r2, unsigned int n)
+{
+  Point<dim> center;
+  Triangulation<dim> tria (Triangulation<dim>::none);
+  GridGenerator::hyper_shell (tria, center, r1, r2, n, true);
+  static const HyperShellBoundary<dim> boundary(center);
+  tria.set_boundary(0, boundary);
+
+  for (typename Triangulation<dim>::cell_iterator cell=tria.begin();
+       cell != tria.end(); ++cell)
+    for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
+      if (cell->face(f)->at_boundary())
+       for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_face; ++l)
+         Assert (cell->face(f)->line(l)->boundary_indicator()
+                 ==
+                 cell->face(f)->boundary_indicator(),
+                 ExcInternalError());
+
+  deallog << "OK" << std::endl;
+}
+
+
+int main()
+{
+  deallog << std::setprecision(3);
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  check<3> (.5, 1, 6);
+  check<3> (.5, 1, 12);
+  check<3> (.5, 1, 96);
+}
diff --git a/tests/deal.II/grid_hyper_shell_06.output b/tests/deal.II/grid_hyper_shell_06.output
new file mode 100644 (file)
index 0000000..fb71de2
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::OK
+DEAL::OK
+DEAL::OK

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.