]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Update step-49. 6370/head
authorDavid Wells <wellsd2@rpi.edu>
Sat, 28 Apr 2018 03:01:40 +0000 (23:01 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Sat, 28 Apr 2018 03:01:40 +0000 (23:01 -0400)
This clarifies the use of manifolds since GridGenerator now attaches manifolds
by default.

examples/step-49/doc/intro.dox
examples/step-49/doc/results.dox
examples/step-49/step-49.cc

index 0c6a3c8bc80f3840f4fbbf504a61cddc3642f0f6..381203457870f5e19097b89c7d82f336a936d878 100644 (file)
@@ -161,7 +161,7 @@ coordinate of a mesh with a sine curve:
 
 Similarly, we can transform a regularly refined
 unit square to a wall-adapted mesh in y direction using the formula
-$(x,y) \mapsto (x,\tanh(2*y)/\tanh(2))$. This is done in <code>grid_6()</code>
+$(x,y) \mapsto (x,\tanh(2 y)/\tanh(2))$. This is done in <code>grid_6()</code>
 of this tutorial:
 <TABLE WIDTH="60%" ALIGN="center">
   <tr>
index f1b381bf1a4d35fc3cdcc1ec0f3a7a52c086b581..50d3f8f014ff8900ffcd177c3d856caeafc88134 100644 (file)
@@ -40,16 +40,18 @@ to use a curved geometry. The way to do this requires three steps:
   <code>manifold_id</code> with a Manifold object. For more information on this
   see the @ref GlossManifoldIndicator "glossary entry on this topic".
 - Finally, you must mark cells and cell faces with the correct
-  <code>manifold_id</code>. For example, you could get an annular domain with
+  <code>manifold_id</code>. For example, you could get an annular sector with
   curved cells in Cartesian coordinates (but rectangles in polar coordinates)
   by doing the following:
   @code
-  PolarManifold<2> polar_manifold;
-  Triangulation<2> triangulation;
-  const types::manifold_id polar_id = 42;
-  GridGenerator::hyper_shell(triangulation, Point<2>(), 0.5, 1.0, 10);
-  triangulation.set_manifold(polar_id, polar_manifold);
-  triangulation.set_all_manifold_ids(polar_id);
+  Triangulation<2> tria;
+  GridGenerator::hyper_cube(tria);
+  const auto cell = tria.begin_active();
+  cell->vertex(2) = Point<2>(-0.5, 1.0);
+  cell->vertex(3) = Point<2>(1.5, 1.0);
+  tria.set_all_manifold_ids(42);
+  tria.set_manifold(42, PolarManifold<2>(Point<2>(0.5, -1.0)));
+  tria.refine_global(3);
   @endcode
   Now, when the grid is refined, all cell splitting calculations will be done in
   polar coordinates.
index 80d777c67abd8103d88faa2cedc14295dc4aa626..cf2dd9a620dc93e98e4d4316905f0202249d5a5d 100644 (file)
@@ -181,30 +181,15 @@ void grid_3 ()
         }
     }
 
-  // In the second step we will refine the mesh twice. To do this
-  // correctly, we have to associate a geometry object with the
-  // boundary of the hole; since the boundary of the hole has boundary
-  // indicator 1 (see the documentation of the function that generates
-  // the mesh), we need to create an object that describes a spherical
-  // manifold (i.e., a hyper ball) with appropriate center and assign
-  // it to the triangulation. Notice that the function that generates
-  // the triangulation sets the boundary indicators of the inner mesh,
-  // but leaves unchanged the manifold indicator. We copy the boundary
-  // indicator to the manifold indicators in order for the object to
-  // be refined accordingly.
-  // We can then refine twice:
-  GridTools::copy_boundary_to_manifold_id(triangulation);
-  const SphericalManifold<2> boundary_description(Point<2>(0,0));
-  triangulation.set_manifold (1, boundary_description);
+  // In the second step we will refine the mesh twice. To do this correctly,
+  // we should place new points on the interior boundary along the surface of
+  // a circle centered at the origin. Fortunately,
+  // GridGenerator::hyper_cube_with_cylindrical_hole already attaches a
+  // Manifold object to the interior boundary, so we do not need to do
+  // anything but refine the mesh (see the @ref Results results section for a
+  // fully worked example where we <em>do</em> attach a Manifold object).
   triangulation.refine_global(2);
-
-  // The mesh so generated is then passed to the function that generates
-  // output. In a final step we remove the boundary object again so that it is
-  // no longer in use by the triangulation when it is destroyed (the boundary
-  // object is destroyed first in this function since it was declared after
-  // the triangulation).
   print_mesh_info (triangulation, "grid-3.eps");
-  triangulation.reset_manifold(1);
 }
 
 // There is one snag to doing things as shown above: If one moves the nodes on

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.