]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added TriaAccessor::point function, which returns points on the manifold given their...
authorLuca Heltai <luca.heltai@sissa.it>
Tue, 26 Aug 2014 10:25:10 +0000 (12:25 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Sat, 13 Sep 2014 12:15:18 +0000 (14:15 +0200)
Simplified Manifolds::get_default_quadrature.

14 files changed:
doc/news/changes.h
include/deal.II/grid/manifold.h
include/deal.II/grid/tria_accessor.h
include/deal.II/grid/tria_accessor.templates.h
source/grid/manifold.cc
source/grid/manifold_lib.cc
source/grid/manifold_lib.inst.in
source/grid/tria_accessor.cc
tests/manifold/tria_accessor_point_01.cc [new file with mode: 0644]
tests/manifold/tria_accessor_point_01.output [new file with mode: 0644]
tests/manifold/tria_accessor_point_02.cc [new file with mode: 0644]
tests/manifold/tria_accessor_point_02.output [new file with mode: 0644]
tests/manifold/tria_accessor_point_03.cc [new file with mode: 0644]
tests/manifold/tria_accessor_point_03.output [new file with mode: 0644]

index b5cadb98a32b8431d12adf577e53e01428c6c097..1dbfb4ad744339a05008d395d0be08e8de133d17 100644 (file)
@@ -97,11 +97,22 @@ inconvenience this causes.
 
 
 <ol>
+  
+  <li> New: Added two optional parameters to TriaAccessor::center()
+  and a new method TriaAccessor::point(). They allow to query for a 
+  geometrically coherent center, or ask for arbitrary points on the 
+  underlying Manifold, given the dim coordinates in the reference 
+  element.
+  <br>
+  (Luca Heltai, 2014/09/12)
+
+
   <li> New: The new tutorial program step-52 explains how to use the 
   new time stepping methods.
   <br>
   (Bruno Turcksin, Damien Lebrun-Grandie, 2014/09/12)
 
+
   <li> New: The new tutorial program step-53 explains how to deal with
   complicated geometries.
   <br>
index 02dd9db50d39d7c28dc217995e79477cda7067cf..bfb66354778ed65bb98e4a93caa19d153de0f5f7 100644 (file)
@@ -38,15 +38,6 @@ template <int dim, int space_dim> class Triangulation;
  */
 namespace Manifolds
 {
-  /**
-   * Given a hex iterator, construct a quadrature with the Laplace
-   * weigths, and all relevant points of the hex: vertices, line
-   * centers and face centers, which can be called when creating
-   * middle vertices in the manifold routines.
-   */
-  Quadrature<3>
-  get_default_quadrature(const TriaIterator<CellAccessor<3, 3> > &hex);
-
   /**
     * Given a general mesh iterator, construct a quadrature with the
     * Laplace weights or with uniform weights according the parameter
@@ -54,8 +45,8 @@ namespace Manifolds
     * vertices, line centers and/or face centers, which can be called
     * when creating new vertices in the manifold routines.
     */
-  template <typename OBJECT, int spacedim>
-  Quadrature<spacedim>
+  template <typename OBJECT>
+  Quadrature<OBJECT::AccessorType::space_dimension>
   get_default_quadrature(const OBJECT &obj, bool with_laplace = false);
 }
 
@@ -527,20 +518,21 @@ get_new_point_on_hex (const Triangulation<3,3>::hex_iterator &) const;
 namespace Manifolds
 {
 
-  template <typename OBJECT, int spacedim>
-  Quadrature<spacedim>
+  template <typename OBJECT>
+  Quadrature<OBJECT::AccessorType::space_dimension>
   get_default_quadrature(const OBJECT &obj, bool with_laplace)
   {
+    const int spacedim = OBJECT::AccessorType::space_dimension;
+    const int dim = OBJECT::AccessorType::structure_dimension;
+
     std::vector<Point<spacedim> > sp;
     std::vector<double> wp;
 
-    const int dim = OBJECT::AccessorType::structure_dimension;
 
     // note that the exact weights are chosen such as to minimize the
     // distortion of the four new quads from the optimal shape; their
     // derivation and values is copied over from the
     // @p{MappingQ::set_laplace_on_vector} function
-    AssertDimension(spacedim, OBJECT::AccessorType::space_dimension);
     switch (dim)
       {
       case 1:
@@ -554,23 +546,15 @@ namespace Manifolds
       case 2:
         sp.resize(8);
         wp.resize(8);
-        sp[0] = obj->vertex(0);
-        sp[1] = obj->vertex(1);
-        sp[2] = obj->vertex(2);
-        sp[3] = obj->vertex(3);
-
-        sp[4] = obj->line(0)->has_children() ?
-                obj->line(0)->child(0)->vertex(1) :
-                obj->line(0)->get_manifold().get_new_point_on_line(obj->line(0));
-        sp[5] = obj->line(1)->has_children() ?
-                obj->line(1)->child(0)->vertex(1) :
-                obj->line(1)->get_manifold().get_new_point_on_line(obj->line(1));
-        sp[6] = obj->line(2)->has_children() ?
-                obj->line(2)->child(0)->vertex(1) :
-                obj->line(2)->get_manifold().get_new_point_on_line(obj->line(2));
-        sp[7] = obj->line(3)->has_children() ?
-                obj->line(3)->child(0)->vertex(1) :
-                obj->line(3)->get_manifold().get_new_point_on_line(obj->line(3));
+
+        for (unsigned int i=0; i<4; ++i)
+          {
+            sp[i] = obj->vertex(i);
+            sp[4+i] = ( obj->line(i)->has_children() ?
+                        obj->line(i)->child(0)->vertex(1) :
+                        obj->line(i)->get_manifold().get_new_point_on_line(obj->line(i)) );
+          }
+
         if (with_laplace)
           {
             std::fill(wp.begin(), wp.begin()+4, 1.0/16.0);
@@ -579,6 +563,49 @@ namespace Manifolds
         else
           std::fill(wp.begin(), wp.end(), 1.0/8.0);
         break;
+      case 3:
+      {
+        TriaIterator<TriaAccessor<3, 3, 3> > hex
+          = static_cast<TriaIterator<TriaAccessor<3, 3, 3> > >(obj);
+        const unsigned int np =
+          GeometryInfo<dim>::vertices_per_cell+
+          GeometryInfo<dim>::lines_per_cell+
+          GeometryInfo<dim>::faces_per_cell;
+        sp.resize(np);
+        wp.resize(np);
+        std::vector<Point<3> > *sp3 = reinterpret_cast<std::vector<Point<3> > *>(&sp);
+
+        unsigned int j=0;
+
+        // note that the exact weights are chosen such as to minimize the
+        // distortion of the eight new hexes from the optimal shape; their
+        // derivation and values is copied over from the
+        // @p{MappingQ::set_laplace_on_vector} function
+        for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i, ++j)
+          {
+            (*sp3)[j] = hex->vertex(i);
+            wp[j] = 1.0/128.0;
+          }
+        for (unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i, ++j)
+          {
+            (*sp3)[j] = (hex->line(i)->has_children() ?
+                         hex->line(i)->child(0)->vertex(1) :
+                         hex->line(i)->get_manifold().get_new_point_on_line(hex->line(i)));
+            wp[j] = 7.0/192.0;
+          }
+        for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i, ++j)
+          {
+            (*sp3)[j] = (hex->quad(i)->has_children() ?
+                         hex->quad(i)->isotropic_child(0)->vertex(3) :
+                         hex->quad(i)->get_manifold().get_new_point_on_quad(hex->quad(i)));
+            wp[j] = 1.0/12.0;
+          }
+        // Overwrited the weights with 1/np if we don't want to use
+        // laplace vectors.
+        if (with_laplace == false)
+          std::fill(wp.begin(), wp.end(), 1.0/np);
+      }
+      break;
       default:
         Assert(false, ExcInternalError());
         break;
index 6aea7debe7df7514defdf5a0c0a4a4635025b635..b5433347b472ad454853f26f02ad6e653bcc6812 100644 (file)
@@ -1539,18 +1539,35 @@ public:
   double minimum_vertex_distance () const;
 
   /**
-   * Center of the object. The center of an
-   * object is defined to be the average of
-   * the locations of the vertices, which
-   * is also where the (dim-)linear mapping
-   * places the midpoint of the unit cell
-   * in real space.  However, this may not
-   * be the barycenter of the object and it
-   * may also not be the true center of an
-   * object if higher order mappings are
-   * used.
-   */
-  Point<spacedim> center () const;
+   * Returns a point belonging to the Manifold<dim,spacedim> where
+   * this object lives, given its parametric coordinates on the
+   * reference #strucdim cell. This function queries the underlying
+   * manifold object, and can be used to obtain the exact geometrical
+   * location of arbitrary points on this object.
+   */
+  Point<spacedim> point(const Point<structdim> &coordinates) const;
+
+  /**
+   * Center of the object. The center of an object is defined to be
+   * the average of the locations of the vertices. If required, the
+   * user may ask this function to return the average of the point
+   * according to the underlyinging Manifold object, by setting to
+   * true the optional parameter @p respect_manifold.
+   *
+   * When the geometry of a TriaAccessor is not flat, or when part of
+   * the bounding objects of this TriaAccessor are not flat, the
+   * result given by the TriaAccessor::center() function may not be
+   * accurate enough, even when parameter @p respect_manifold is set
+   * to true. If you find this to be case, than you can further refine
+   * the computation of the center by setting to true the second
+   * additional parameter @p use_laplace_transformation, which will
+   * force this function to compute the center performing a Laplace
+   * transformation on all bounding support points. The Laplace
+   * transformation is computed similarly to what happens in
+   * @p{MappingQ::set_laplace_on_vector}.
+   */
+  Point<spacedim> center (const bool respect_manifold=false,
+                          const bool use_laplace_transformation=false) const;
 
   /**
    * Barycenter of the object.
index 842bd1633774844ca7030dbd89b3779bfa2e0243..aad3196b533970e8629eff5b16f240501295fac2 100644 (file)
@@ -2060,18 +2060,6 @@ TriaAccessor<structdim, dim, spacedim>::minimum_vertex_distance () const
 }
 
 
-
-template <int structdim, int dim, int spacedim>
-Point<spacedim>
-TriaAccessor<structdim, dim, spacedim>::center () const
-{
-  Point<spacedim> p;
-  for (unsigned int v=0; v<GeometryInfo<structdim>::vertices_per_cell; ++v)
-    p += vertex(v);
-  return p/GeometryInfo<structdim>::vertices_per_cell;
-}
-
-
 template <int structdim, int dim, int spacedim>
 bool
 TriaAccessor<structdim, dim, spacedim>::
index 2835521035e3585123762e22d48f89f0f84f3793..28b156c3d71c22ba61572c9561989286dce0e65e 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
-namespace Manifolds
-{
-
-  Quadrature<3>
-  get_default_quadrature(const TriaIterator<CellAccessor<3, 3> > &obj)
-  {
-    std::vector<Point<3> > sp;
-    std::vector<double> wp;
-
-    const int dim = 3;
-
-    const unsigned int np =
-      GeometryInfo<dim>::vertices_per_cell+
-      GeometryInfo<dim>::lines_per_cell+
-      GeometryInfo<dim>::faces_per_cell;
-
-    sp.resize(np);
-    wp.resize(np);
-    unsigned int j=0;
-
-    // note that the exact weights are chosen such as to minimize the
-    // distortion of the eight new hexes from the optimal shape; their
-    // derivation and values is copied over from the
-    // @p{MappingQ::set_laplace_on_vector} function
-    for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i, ++j)
-      {
-        sp[j] = obj->vertex(i);
-        wp[j] = 1.0/128.0;
-      }
-    for (unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i, ++j)
-      {
-        sp[j] = (obj->line(i)->has_children() ? obj->line(i)->child(0)->vertex(1) :
-                 obj->line(i)->get_manifold().get_new_point_on_line(obj->line(i)));
-        wp[j] = 7.0/192.0;
-      }
-    for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i, ++j)
-      {
-        sp[j] = (obj->face(i)->has_children() ? obj->face(i)->isotropic_child(0)->vertex(3) :
-                 obj->face(i)->get_manifold().get_new_point_on_face(obj->face(i)));
-        wp[j] = 1.0/12.0;
-      }
-    return Quadrature<3>(sp,wp);
-  }
-
-}
-
 using namespace Manifolds;
 
 /* -------------------------- Manifold --------------------- */
@@ -121,9 +75,7 @@ Point<spacedim>
 Manifold<dim, spacedim>::
 get_new_point_on_line (const typename Triangulation<dim, spacedim>::line_iterator &line) const
 {
-  return get_new_point
-         (get_default_quadrature<const typename Triangulation<dim, spacedim>::line_iterator,
-          spacedim>(line, false));
+  return get_new_point (get_default_quadrature(line));
 }
 
 
@@ -133,9 +85,7 @@ Point<spacedim>
 Manifold<dim, spacedim>::
 get_new_point_on_quad (const typename Triangulation<dim, spacedim>::quad_iterator &quad) const
 {
-  return get_new_point
-         (get_default_quadrature<const typename Triangulation<dim, spacedim>::quad_iterator,
-          spacedim>(quad, false));
+  return get_new_point (get_default_quadrature(quad));
 }
 
 
@@ -251,7 +201,7 @@ Point<3>
 Manifold<3,3>::
 get_new_point_on_hex (const Triangulation<3, 3>::hex_iterator &hex) const
 {
-  return get_new_point(get_default_quadrature(hex));
+  return get_new_point(get_default_quadrature(hex, true));
 }
 
 
index e9a0e9c2d1bcd185f43cd73b4205bed8fb4ecebe..a10024d89c81156f09ed90f90a0e0fdb21cf25f4 100644 (file)
@@ -28,9 +28,7 @@ template <int dim, int spacedim>
 SphericalManifold<dim,spacedim>::SphericalManifold(const Point<spacedim> center):
   ChartManifold<dim,spacedim,spacedim>(SphericalManifold<dim,spacedim>::get_periodicity()),
   center(center)
-{
-  Assert(spacedim != 1, ExcImpossibleInDim(1));
-}
+{}
 
 
 
index 00a6d200c0cc7de9c4e11e7366ecc7dcb4ea9e9d..7825aa6764f55e4537db397a7a91c26ac7ba093c 100644 (file)
@@ -17,7 +17,7 @@
 
 for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
   {
-#if deal_II_space_dimension > 1 && deal_II_dimension <= deal_II_space_dimension
+#if deal_II_dimension <= deal_II_space_dimension
     template class SphericalManifold<deal_II_dimension, deal_II_space_dimension>;
 #endif
 #if deal_II_dimension <= deal_II_space_dimension
index 4e197f3bf1d6b8b88a2438eb1083e86a0bbddef8..2e86b69cd45308848c13afb52be7d1dd00eb8882 100644 (file)
 #include <deal.II/grid/tria_accessor.templates.h>
 #include <deal.II/grid/tria_iterator.templates.h>
 #include <deal.II/base/geometry_info.h>
+#include <deal.II/base/quadrature.h>
 #include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/manifold.h>
 #include <deal.II/fe/mapping_q1.h>
+#include <deal.II/fe/fe_q.h>
 
 #include <cmath>
 
 DEAL_II_NAMESPACE_OPEN
 
-
 // anonymous namespace for helper functions
 namespace
 {
@@ -1019,6 +1021,27 @@ namespace
     Assert (false, ExcNotImplemented());
     return std::numeric_limits<double>::quiet_NaN();
   }
+
+  template <int dim, int spacedim>
+  Point<spacedim> get_new_point_on_obj(const TriaAccessor<1, dim, spacedim> &obj)
+  {
+    TriaIterator<TriaAccessor<1,dim,spacedim> > it(obj);
+    return obj.get_manifold().get_new_point_on_line(it);
+  }
+
+  template <int dim, int spacedim>
+  Point<spacedim> get_new_point_on_obj(const TriaAccessor<2, dim, spacedim> &obj)
+  {
+    TriaIterator<TriaAccessor<2,dim,spacedim> > it(obj);
+    return obj.get_manifold().get_new_point_on_quad(it);
+  }
+
+  template <int dim, int spacedim>
+  Point<spacedim> get_new_point_on_obj(const TriaAccessor<3, dim, spacedim> &obj)
+  {
+    TriaIterator<TriaAccessor<3,dim,spacedim> > it(obj);
+    return obj.get_manifold().get_new_point_on_hex(it);
+  }
 }
 
 
@@ -1161,6 +1184,51 @@ set_all_manifold_ids (const types::manifold_id manifold_ind) const
 }
 
 
+template <int structdim, int dim, int spacedim>
+Point<spacedim>
+TriaAccessor<structdim, dim, spacedim>::point (const Point<structdim> &coordinates) const
+{
+  // We use an FE_Q<structdim>(1) to extract the "weights" of each
+  // vertex, used to get a point from the manifold.
+  static FE_Q<structdim> fe(1);
+
+  // Surrounding points and weights.
+  std::vector<Point<spacedim> > p(GeometryInfo<structdim>::vertices_per_cell);
+  std::vector<double>   w(GeometryInfo<structdim>::vertices_per_cell);
+
+  for (unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
+    {
+      p[i] = this->vertex(i);
+      w[i] = fe.shape_value(i, coordinates);
+    }
+
+  Quadrature<spacedim> quadrature(p, w);
+  return this->get_manifold().get_new_point(quadrature);
+}
+
+
+template <int structdim, int dim, int spacedim>
+Point<spacedim>
+TriaAccessor<structdim, dim, spacedim>::center (const bool respect_manifold,
+                                                const bool use_laplace) const
+{
+  if (respect_manifold == false)
+    {
+      Assert(use_laplace == false, ExcNotImplemented());
+      Point<spacedim> p;
+      for (unsigned int v=0; v<GeometryInfo<structdim>::vertices_per_cell; ++v)
+        p += vertex(v);
+      return p/GeometryInfo<structdim>::vertices_per_cell;
+    }
+  else
+    {
+      TriaRawIterator<TriaAccessor<structdim, dim, spacedim> > it(*this);
+      Quadrature<spacedim> quadrature = Manifolds::get_default_quadrature(it, use_laplace);
+      return this->get_manifold().get_new_point(quadrature);
+    }
+}
+
+
 /*------------------------ Functions: CellAccessor<1> -----------------------*/
 
 
diff --git a/tests/manifold/tria_accessor_point_01.cc b/tests/manifold/tria_accessor_point_01.cc
new file mode 100644 (file)
index 0000000..0742e07
--- /dev/null
@@ -0,0 +1,68 @@
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+//    Copyright (C) 2011, 2013, 2014 by the mathLab team.
+//
+//    This file is subject to LGPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+
+
+// Test spherical manifold on hyper shells.
+
+#include "../tests.h"
+#include <fstream>
+#include <base/logstream.h>
+
+
+// all include files you need here
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/manifold_lib.h>
+#include <deal.II/grid/grid_out.h>
+
+// Helper function
+template <int dim, int spacedim>
+void test(unsigned int ref=1)
+{
+  deallog << "Testing dim " << dim 
+         << ", spacedim " << spacedim << std::endl;
+
+  SphericalManifold<dim,spacedim> manifold;
+  
+  Triangulation<dim,spacedim> tria;
+  GridGenerator::hyper_shell (tria, Point<spacedim>(), .3, .6, 12);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    cell->set_all_manifold_ids(1);
+  }
+  
+  tria.set_manifold(1, manifold);
+  tria.refine_global(1);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    for(unsigned int f=0;f<GeometryInfo<dim>::faces_per_cell; ++f)
+      if(cell->face(f)->at_boundary())
+       deallog << "Center: " << cell->face(f)->center(true) 
+               << ", Norm: " << cell->face(f)->center(true).norm() << std::endl;
+  }
+  
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+  
+  test<2,2>();
+  test<3,3>();
+  
+  return 0;
+}
+
diff --git a/tests/manifold/tria_accessor_point_01.output b/tests/manifold/tria_accessor_point_01.output
new file mode 100644 (file)
index 0000000..627da2c
--- /dev/null
@@ -0,0 +1,147 @@
+
+DEAL::Testing dim 2, spacedim 2
+DEAL::Center: 0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: 0.554328 0.229610, Norm: 0.600000
+DEAL::Center: 0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: 0.277164 0.114805, Norm: 0.300000
+DEAL::Center: 0.476012 0.365257, Norm: 0.600000
+DEAL::Center: 0.365257 0.476012, Norm: 0.600000
+DEAL::Center: 0.238006 0.182628, Norm: 0.300000
+DEAL::Center: 0.182628 0.238006, Norm: 0.300000
+DEAL::Center: 0.229610 0.554328, Norm: 0.600000
+DEAL::Center: 0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: 0.114805 0.277164, Norm: 0.300000
+DEAL::Center: 0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: -0.229610 0.554328, Norm: 0.600000
+DEAL::Center: -0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.114805 0.277164, Norm: 0.300000
+DEAL::Center: -0.365257 0.476012, Norm: 0.600000
+DEAL::Center: -0.476012 0.365257, Norm: 0.600000
+DEAL::Center: -0.182628 0.238006, Norm: 0.300000
+DEAL::Center: -0.238006 0.182628, Norm: 0.300000
+DEAL::Center: -0.554328 0.229610, Norm: 0.600000
+DEAL::Center: -0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: -0.277164 0.114805, Norm: 0.300000
+DEAL::Center: -0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: -0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: -0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: -0.297433 -0.0391579, Norm: 0.300000
+DEAL::Center: -0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: -0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: -0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: -0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: -0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: -0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: -0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: -0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: -0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: 0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: 0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: 0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: 0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: 0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: 0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: 0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: 0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: 0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: 0.297433 -0.0391579, Norm: 0.300000
+DEAL::Testing dim 3, spacedim 3
+DEAL::Center: -0.109034 -0.279484 0.00000, Norm: 0.300000
+DEAL::Center: -0.200401 -0.200401 0.0983833, Norm: 0.300000
+DEAL::Center: -0.200401 -0.200401 -0.0983833, Norm: 0.300000
+DEAL::Center: -0.279484 -0.109034 0.00000, Norm: 0.300000
+DEAL::Center: -0.218068 -0.558969 0.00000, Norm: 0.600000
+DEAL::Center: -0.400801 -0.400801 0.196767, Norm: 0.600000
+DEAL::Center: -0.400801 -0.400801 -0.196767, Norm: 0.600000
+DEAL::Center: -0.558969 -0.218068 0.00000, Norm: 0.600000
+DEAL::Center: -0.200401 -0.0983833 0.200401, Norm: 0.300000
+DEAL::Center: -0.109034 0.00000 0.279484, Norm: 0.300000
+DEAL::Center: -0.279484 0.00000 0.109034, Norm: 0.300000
+DEAL::Center: -0.200401 0.0983833 0.200401, Norm: 0.300000
+DEAL::Center: -0.400801 -0.196767 0.400801, Norm: 0.600000
+DEAL::Center: -0.218068 0.00000 0.558969, Norm: 0.600000
+DEAL::Center: -0.558969 0.00000 0.218068, Norm: 0.600000
+DEAL::Center: -0.400801 0.196767 0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 -0.279484 0.109034, Norm: 0.300000
+DEAL::Center: 0.0983833 -0.200401 0.200401, Norm: 0.300000
+DEAL::Center: -0.0983833 -0.200401 0.200401, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109034 0.279484, Norm: 0.300000
+DEAL::Center: 0.00000 -0.558969 0.218068, Norm: 0.600000
+DEAL::Center: 0.196767 -0.400801 0.400801, Norm: 0.600000
+DEAL::Center: -0.196767 -0.400801 0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218068 0.558969, Norm: 0.600000
+DEAL::Center: 0.200401 -0.200401 -0.0983833, Norm: 0.300000
+DEAL::Center: 0.279484 -0.109034 0.00000, Norm: 0.300000
+DEAL::Center: 0.109034 -0.279484 0.00000, Norm: 0.300000
+DEAL::Center: 0.200401 -0.200401 0.0983833, Norm: 0.300000
+DEAL::Center: 0.400801 -0.400801 -0.196767, Norm: 0.600000
+DEAL::Center: 0.558969 -0.218068 0.00000, Norm: 0.600000
+DEAL::Center: 0.218068 -0.558969 0.00000, Norm: 0.600000
+DEAL::Center: 0.400801 -0.400801 0.196767, Norm: 0.600000
+DEAL::Center: 0.279484 0.00000 0.109034, Norm: 0.300000
+DEAL::Center: 0.200401 0.0983833 0.200401, Norm: 0.300000
+DEAL::Center: 0.200401 -0.0983833 0.200401, Norm: 0.300000
+DEAL::Center: 0.109034 0.00000 0.279484, Norm: 0.300000
+DEAL::Center: 0.558969 0.00000 0.218068, Norm: 0.600000
+DEAL::Center: 0.400801 0.196767 0.400801, Norm: 0.600000
+DEAL::Center: 0.400801 -0.196767 0.400801, Norm: 0.600000
+DEAL::Center: 0.218068 0.00000 0.558969, Norm: 0.600000
+DEAL::Center: 0.0983833 0.200401 0.200401, Norm: 0.300000
+DEAL::Center: 0.00000 0.279484 0.109034, Norm: 0.300000
+DEAL::Center: 0.00000 0.109034 0.279484, Norm: 0.300000
+DEAL::Center: -0.0983833 0.200401 0.200401, Norm: 0.300000
+DEAL::Center: 0.196767 0.400801 0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 0.558969 0.218068, Norm: 0.600000
+DEAL::Center: 0.00000 0.218068 0.558969, Norm: 0.600000
+DEAL::Center: -0.196767 0.400801 0.400801, Norm: 0.600000
+DEAL::Center: 0.279484 0.109034 0.00000, Norm: 0.300000
+DEAL::Center: 0.200401 0.200401 -0.0983833, Norm: 0.300000
+DEAL::Center: 0.200401 0.200401 0.0983833, Norm: 0.300000
+DEAL::Center: 0.109034 0.279484 0.00000, Norm: 0.300000
+DEAL::Center: 0.558969 0.218068 0.00000, Norm: 0.600000
+DEAL::Center: 0.400801 0.400801 -0.196767, Norm: 0.600000
+DEAL::Center: 0.400801 0.400801 0.196767, Norm: 0.600000
+DEAL::Center: 0.218068 0.558969 0.00000, Norm: 0.600000
+DEAL::Center: 0.200401 -0.0983833 -0.200401, Norm: 0.300000
+DEAL::Center: 0.109034 0.00000 -0.279484, Norm: 0.300000
+DEAL::Center: 0.279484 0.00000 -0.109034, Norm: 0.300000
+DEAL::Center: 0.200401 0.0983833 -0.200401, Norm: 0.300000
+DEAL::Center: 0.400801 -0.196767 -0.400801, Norm: 0.600000
+DEAL::Center: 0.218068 0.00000 -0.558969, Norm: 0.600000
+DEAL::Center: 0.558969 0.00000 -0.218068, Norm: 0.600000
+DEAL::Center: 0.400801 0.196767 -0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 0.109034 -0.279484, Norm: 0.300000
+DEAL::Center: -0.0983833 0.200401 -0.200401, Norm: 0.300000
+DEAL::Center: 0.0983833 0.200401 -0.200401, Norm: 0.300000
+DEAL::Center: 0.00000 0.279484 -0.109034, Norm: 0.300000
+DEAL::Center: 0.00000 0.218068 -0.558969, Norm: 0.600000
+DEAL::Center: -0.196767 0.400801 -0.400801, Norm: 0.600000
+DEAL::Center: 0.196767 0.400801 -0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 0.558969 -0.218068, Norm: 0.600000
+DEAL::Center: -0.200401 0.200401 -0.0983833, Norm: 0.300000
+DEAL::Center: -0.279484 0.109034 0.00000, Norm: 0.300000
+DEAL::Center: -0.109034 0.279484 0.00000, Norm: 0.300000
+DEAL::Center: -0.200401 0.200401 0.0983833, Norm: 0.300000
+DEAL::Center: -0.400801 0.400801 -0.196767, Norm: 0.600000
+DEAL::Center: -0.558969 0.218068 0.00000, Norm: 0.600000
+DEAL::Center: -0.218068 0.558969 0.00000, Norm: 0.600000
+DEAL::Center: -0.400801 0.400801 0.196767, Norm: 0.600000
+DEAL::Center: -0.109034 0.00000 -0.279484, Norm: 0.300000
+DEAL::Center: -0.200401 -0.0983833 -0.200401, Norm: 0.300000
+DEAL::Center: -0.200401 0.0983833 -0.200401, Norm: 0.300000
+DEAL::Center: -0.279484 0.00000 -0.109034, Norm: 0.300000
+DEAL::Center: -0.218068 0.00000 -0.558969, Norm: 0.600000
+DEAL::Center: -0.400801 -0.196767 -0.400801, Norm: 0.600000
+DEAL::Center: -0.400801 0.196767 -0.400801, Norm: 0.600000
+DEAL::Center: -0.558969 0.00000 -0.218068, Norm: 0.600000
+DEAL::Center: 0.0983833 -0.200401 -0.200401, Norm: 0.300000
+DEAL::Center: 0.00000 -0.279484 -0.109034, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109034 -0.279484, Norm: 0.300000
+DEAL::Center: -0.0983833 -0.200401 -0.200401, Norm: 0.300000
+DEAL::Center: 0.196767 -0.400801 -0.400801, Norm: 0.600000
+DEAL::Center: 0.00000 -0.558969 -0.218068, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218068 -0.558969, Norm: 0.600000
+DEAL::Center: -0.196767 -0.400801 -0.400801, Norm: 0.600000
diff --git a/tests/manifold/tria_accessor_point_02.cc b/tests/manifold/tria_accessor_point_02.cc
new file mode 100644 (file)
index 0000000..73da3a5
--- /dev/null
@@ -0,0 +1,68 @@
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+//    Copyright (C) 2011, 2013, 2014 by the mathLab team.
+//
+//    This file is subject to LGPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+
+
+// Test spherical manifold on hyper shells.
+
+#include "../tests.h"
+#include <fstream>
+#include <base/logstream.h>
+
+
+// all include files you need here
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/manifold_lib.h>
+#include <deal.II/grid/grid_out.h>
+
+// Helper function
+template <int dim, int spacedim>
+void test(unsigned int ref=1)
+{
+  deallog << "Testing dim " << dim 
+         << ", spacedim " << spacedim << std::endl;
+
+  SphericalManifold<dim,spacedim> manifold;
+  
+  Triangulation<dim,spacedim> tria;
+  GridGenerator::hyper_shell (tria, Point<spacedim>(), .3, .6, 12);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    cell->set_all_manifold_ids(1);
+  }
+  
+  tria.set_manifold(1, manifold);
+  tria.refine_global(1);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    for(unsigned int f=0;f<GeometryInfo<dim>::faces_per_cell; ++f)
+      if(cell->face(f)->at_boundary())
+       deallog << "Center: " << cell->face(f)->center(true, true) 
+               << ", Norm: " << cell->face(f)->center(true, true).norm() << std::endl;
+  }
+  
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+  
+  test<2,2>();
+  test<3,3>();
+  
+  return 0;
+}
+
diff --git a/tests/manifold/tria_accessor_point_02.output b/tests/manifold/tria_accessor_point_02.output
new file mode 100644 (file)
index 0000000..4365298
--- /dev/null
@@ -0,0 +1,147 @@
+
+DEAL::Testing dim 2, spacedim 2
+DEAL::Center: 0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: 0.554328 0.229610, Norm: 0.600000
+DEAL::Center: 0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: 0.277164 0.114805, Norm: 0.300000
+DEAL::Center: 0.476012 0.365257, Norm: 0.600000
+DEAL::Center: 0.365257 0.476012, Norm: 0.600000
+DEAL::Center: 0.238006 0.182628, Norm: 0.300000
+DEAL::Center: 0.182628 0.238006, Norm: 0.300000
+DEAL::Center: 0.229610 0.554328, Norm: 0.600000
+DEAL::Center: 0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: 0.114805 0.277164, Norm: 0.300000
+DEAL::Center: 0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: -0.229610 0.554328, Norm: 0.600000
+DEAL::Center: -0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.114805 0.277164, Norm: 0.300000
+DEAL::Center: -0.365257 0.476012, Norm: 0.600000
+DEAL::Center: -0.476012 0.365257, Norm: 0.600000
+DEAL::Center: -0.182628 0.238006, Norm: 0.300000
+DEAL::Center: -0.238006 0.182628, Norm: 0.300000
+DEAL::Center: -0.554328 0.229610, Norm: 0.600000
+DEAL::Center: -0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: -0.277164 0.114805, Norm: 0.300000
+DEAL::Center: -0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: -0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: -0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: -0.297433 -0.0391579, Norm: 0.300000
+DEAL::Center: -0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: -0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: -0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: -0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: -0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: -0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: -0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: -0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: -0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: 0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: 0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: 0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: 0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: 0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: 0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: 0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: 0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: 0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: 0.297433 -0.0391579, Norm: 0.300000
+DEAL::Testing dim 3, spacedim 3
+DEAL::Center: -0.109091 -0.279462 0.00000, Norm: 0.300000
+DEAL::Center: -0.200412 -0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: -0.200412 -0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: -0.279462 -0.109091 0.00000, Norm: 0.300000
+DEAL::Center: -0.218183 -0.558924 0.00000, Norm: 0.600000
+DEAL::Center: -0.400824 -0.400824 0.196675, Norm: 0.600000
+DEAL::Center: -0.400824 -0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: -0.558924 -0.218183 0.00000, Norm: 0.600000
+DEAL::Center: -0.200412 -0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: -0.109091 0.00000 0.279462, Norm: 0.300000
+DEAL::Center: -0.279462 0.00000 0.109091, Norm: 0.300000
+DEAL::Center: -0.200412 0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: -0.400824 -0.196675 0.400824, Norm: 0.600000
+DEAL::Center: -0.218183 0.00000 0.558924, Norm: 0.600000
+DEAL::Center: -0.558924 0.00000 0.218183, Norm: 0.600000
+DEAL::Center: -0.400824 0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.279462 0.109091, Norm: 0.300000
+DEAL::Center: 0.0983373 -0.200412 0.200412, Norm: 0.300000
+DEAL::Center: -0.0983373 -0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109091 0.279462, Norm: 0.300000
+DEAL::Center: 0.00000 -0.558924 0.218183, Norm: 0.600000
+DEAL::Center: 0.196675 -0.400824 0.400824, Norm: 0.600000
+DEAL::Center: -0.196675 -0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218183 0.558924, Norm: 0.600000
+DEAL::Center: 0.200412 -0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: 0.279462 -0.109091 0.00000, Norm: 0.300000
+DEAL::Center: 0.109091 -0.279462 0.00000, Norm: 0.300000
+DEAL::Center: 0.200412 -0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: 0.400824 -0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: 0.558924 -0.218183 0.00000, Norm: 0.600000
+DEAL::Center: 0.218183 -0.558924 0.00000, Norm: 0.600000
+DEAL::Center: 0.400824 -0.400824 0.196675, Norm: 0.600000
+DEAL::Center: 0.279462 0.00000 0.109091, Norm: 0.300000
+DEAL::Center: 0.200412 0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: 0.200412 -0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: 0.109091 0.00000 0.279462, Norm: 0.300000
+DEAL::Center: 0.558924 0.00000 0.218183, Norm: 0.600000
+DEAL::Center: 0.400824 0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.400824 -0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.218183 0.00000 0.558924, Norm: 0.600000
+DEAL::Center: 0.0983373 0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 0.279462 0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 0.109091 0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.196675 0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.558924 0.218183, Norm: 0.600000
+DEAL::Center: 0.00000 0.218183 0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.279462 0.109091 0.00000, Norm: 0.300000
+DEAL::Center: 0.200412 0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: 0.200412 0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: 0.109091 0.279462 0.00000, Norm: 0.300000
+DEAL::Center: 0.558924 0.218183 0.00000, Norm: 0.600000
+DEAL::Center: 0.400824 0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: 0.400824 0.400824 0.196675, Norm: 0.600000
+DEAL::Center: 0.218183 0.558924 0.00000, Norm: 0.600000
+DEAL::Center: 0.200412 -0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: 0.109091 0.00000 -0.279462, Norm: 0.300000
+DEAL::Center: 0.279462 0.00000 -0.109091, Norm: 0.300000
+DEAL::Center: 0.200412 0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: 0.400824 -0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: 0.218183 0.00000 -0.558924, Norm: 0.600000
+DEAL::Center: 0.558924 0.00000 -0.218183, Norm: 0.600000
+DEAL::Center: 0.400824 0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.109091 -0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.0983373 0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 0.279462 -0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 0.218183 -0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.196675 0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.558924 -0.218183, Norm: 0.600000
+DEAL::Center: -0.200412 0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: -0.279462 0.109091 0.00000, Norm: 0.300000
+DEAL::Center: -0.109091 0.279462 0.00000, Norm: 0.300000
+DEAL::Center: -0.200412 0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: -0.400824 0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: -0.558924 0.218183 0.00000, Norm: 0.600000
+DEAL::Center: -0.218183 0.558924 0.00000, Norm: 0.600000
+DEAL::Center: -0.400824 0.400824 0.196675, Norm: 0.600000
+DEAL::Center: -0.109091 0.00000 -0.279462, Norm: 0.300000
+DEAL::Center: -0.200412 -0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: -0.200412 0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: -0.279462 0.00000 -0.109091, Norm: 0.300000
+DEAL::Center: -0.218183 0.00000 -0.558924, Norm: 0.600000
+DEAL::Center: -0.400824 -0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: -0.400824 0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: -0.558924 0.00000 -0.218183, Norm: 0.600000
+DEAL::Center: 0.0983373 -0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 -0.279462 -0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109091 -0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 -0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.196675 -0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.558924 -0.218183, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218183 -0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 -0.400824 -0.400824, Norm: 0.600000
diff --git a/tests/manifold/tria_accessor_point_03.cc b/tests/manifold/tria_accessor_point_03.cc
new file mode 100644 (file)
index 0000000..73da3a5
--- /dev/null
@@ -0,0 +1,68 @@
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+//    Copyright (C) 2011, 2013, 2014 by the mathLab team.
+//
+//    This file is subject to LGPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  spherical_manifold_01.cc  ---------------------------
+
+
+// Test spherical manifold on hyper shells.
+
+#include "../tests.h"
+#include <fstream>
+#include <base/logstream.h>
+
+
+// all include files you need here
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/manifold_lib.h>
+#include <deal.II/grid/grid_out.h>
+
+// Helper function
+template <int dim, int spacedim>
+void test(unsigned int ref=1)
+{
+  deallog << "Testing dim " << dim 
+         << ", spacedim " << spacedim << std::endl;
+
+  SphericalManifold<dim,spacedim> manifold;
+  
+  Triangulation<dim,spacedim> tria;
+  GridGenerator::hyper_shell (tria, Point<spacedim>(), .3, .6, 12);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    cell->set_all_manifold_ids(1);
+  }
+  
+  tria.set_manifold(1, manifold);
+  tria.refine_global(1);
+
+  for(typename Triangulation<dim,spacedim>::active_cell_iterator cell = tria.begin_active(); cell != tria.end(); ++cell) {
+    for(unsigned int f=0;f<GeometryInfo<dim>::faces_per_cell; ++f)
+      if(cell->face(f)->at_boundary())
+       deallog << "Center: " << cell->face(f)->center(true, true) 
+               << ", Norm: " << cell->face(f)->center(true, true).norm() << std::endl;
+  }
+  
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+  
+  test<2,2>();
+  test<3,3>();
+  
+  return 0;
+}
+
diff --git a/tests/manifold/tria_accessor_point_03.output b/tests/manifold/tria_accessor_point_03.output
new file mode 100644 (file)
index 0000000..4365298
--- /dev/null
@@ -0,0 +1,147 @@
+
+DEAL::Testing dim 2, spacedim 2
+DEAL::Center: 0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: 0.554328 0.229610, Norm: 0.600000
+DEAL::Center: 0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: 0.277164 0.114805, Norm: 0.300000
+DEAL::Center: 0.476012 0.365257, Norm: 0.600000
+DEAL::Center: 0.365257 0.476012, Norm: 0.600000
+DEAL::Center: 0.238006 0.182628, Norm: 0.300000
+DEAL::Center: 0.182628 0.238006, Norm: 0.300000
+DEAL::Center: 0.229610 0.554328, Norm: 0.600000
+DEAL::Center: 0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: 0.114805 0.277164, Norm: 0.300000
+DEAL::Center: 0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.0783157 0.594867, Norm: 0.600000
+DEAL::Center: -0.229610 0.554328, Norm: 0.600000
+DEAL::Center: -0.0391579 0.297433, Norm: 0.300000
+DEAL::Center: -0.114805 0.277164, Norm: 0.300000
+DEAL::Center: -0.365257 0.476012, Norm: 0.600000
+DEAL::Center: -0.476012 0.365257, Norm: 0.600000
+DEAL::Center: -0.182628 0.238006, Norm: 0.300000
+DEAL::Center: -0.238006 0.182628, Norm: 0.300000
+DEAL::Center: -0.554328 0.229610, Norm: 0.600000
+DEAL::Center: -0.594867 0.0783157, Norm: 0.600000
+DEAL::Center: -0.277164 0.114805, Norm: 0.300000
+DEAL::Center: -0.297433 0.0391579, Norm: 0.300000
+DEAL::Center: -0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: -0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: -0.297433 -0.0391579, Norm: 0.300000
+DEAL::Center: -0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: -0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: -0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: -0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: -0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: -0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: -0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: -0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: -0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.0783157 -0.594867, Norm: 0.600000
+DEAL::Center: 0.229610 -0.554328, Norm: 0.600000
+DEAL::Center: 0.0391579 -0.297433, Norm: 0.300000
+DEAL::Center: 0.114805 -0.277164, Norm: 0.300000
+DEAL::Center: 0.365257 -0.476012, Norm: 0.600000
+DEAL::Center: 0.476012 -0.365257, Norm: 0.600000
+DEAL::Center: 0.182628 -0.238006, Norm: 0.300000
+DEAL::Center: 0.238006 -0.182628, Norm: 0.300000
+DEAL::Center: 0.554328 -0.229610, Norm: 0.600000
+DEAL::Center: 0.594867 -0.0783157, Norm: 0.600000
+DEAL::Center: 0.277164 -0.114805, Norm: 0.300000
+DEAL::Center: 0.297433 -0.0391579, Norm: 0.300000
+DEAL::Testing dim 3, spacedim 3
+DEAL::Center: -0.109091 -0.279462 0.00000, Norm: 0.300000
+DEAL::Center: -0.200412 -0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: -0.200412 -0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: -0.279462 -0.109091 0.00000, Norm: 0.300000
+DEAL::Center: -0.218183 -0.558924 0.00000, Norm: 0.600000
+DEAL::Center: -0.400824 -0.400824 0.196675, Norm: 0.600000
+DEAL::Center: -0.400824 -0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: -0.558924 -0.218183 0.00000, Norm: 0.600000
+DEAL::Center: -0.200412 -0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: -0.109091 0.00000 0.279462, Norm: 0.300000
+DEAL::Center: -0.279462 0.00000 0.109091, Norm: 0.300000
+DEAL::Center: -0.200412 0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: -0.400824 -0.196675 0.400824, Norm: 0.600000
+DEAL::Center: -0.218183 0.00000 0.558924, Norm: 0.600000
+DEAL::Center: -0.558924 0.00000 0.218183, Norm: 0.600000
+DEAL::Center: -0.400824 0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.279462 0.109091, Norm: 0.300000
+DEAL::Center: 0.0983373 -0.200412 0.200412, Norm: 0.300000
+DEAL::Center: -0.0983373 -0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109091 0.279462, Norm: 0.300000
+DEAL::Center: 0.00000 -0.558924 0.218183, Norm: 0.600000
+DEAL::Center: 0.196675 -0.400824 0.400824, Norm: 0.600000
+DEAL::Center: -0.196675 -0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218183 0.558924, Norm: 0.600000
+DEAL::Center: 0.200412 -0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: 0.279462 -0.109091 0.00000, Norm: 0.300000
+DEAL::Center: 0.109091 -0.279462 0.00000, Norm: 0.300000
+DEAL::Center: 0.200412 -0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: 0.400824 -0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: 0.558924 -0.218183 0.00000, Norm: 0.600000
+DEAL::Center: 0.218183 -0.558924 0.00000, Norm: 0.600000
+DEAL::Center: 0.400824 -0.400824 0.196675, Norm: 0.600000
+DEAL::Center: 0.279462 0.00000 0.109091, Norm: 0.300000
+DEAL::Center: 0.200412 0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: 0.200412 -0.0983373 0.200412, Norm: 0.300000
+DEAL::Center: 0.109091 0.00000 0.279462, Norm: 0.300000
+DEAL::Center: 0.558924 0.00000 0.218183, Norm: 0.600000
+DEAL::Center: 0.400824 0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.400824 -0.196675 0.400824, Norm: 0.600000
+DEAL::Center: 0.218183 0.00000 0.558924, Norm: 0.600000
+DEAL::Center: 0.0983373 0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 0.279462 0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 0.109091 0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 0.200412 0.200412, Norm: 0.300000
+DEAL::Center: 0.196675 0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.558924 0.218183, Norm: 0.600000
+DEAL::Center: 0.00000 0.218183 0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 0.400824 0.400824, Norm: 0.600000
+DEAL::Center: 0.279462 0.109091 0.00000, Norm: 0.300000
+DEAL::Center: 0.200412 0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: 0.200412 0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: 0.109091 0.279462 0.00000, Norm: 0.300000
+DEAL::Center: 0.558924 0.218183 0.00000, Norm: 0.600000
+DEAL::Center: 0.400824 0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: 0.400824 0.400824 0.196675, Norm: 0.600000
+DEAL::Center: 0.218183 0.558924 0.00000, Norm: 0.600000
+DEAL::Center: 0.200412 -0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: 0.109091 0.00000 -0.279462, Norm: 0.300000
+DEAL::Center: 0.279462 0.00000 -0.109091, Norm: 0.300000
+DEAL::Center: 0.200412 0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: 0.400824 -0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: 0.218183 0.00000 -0.558924, Norm: 0.600000
+DEAL::Center: 0.558924 0.00000 -0.218183, Norm: 0.600000
+DEAL::Center: 0.400824 0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.109091 -0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.0983373 0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 0.279462 -0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 0.218183 -0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.196675 0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 0.558924 -0.218183, Norm: 0.600000
+DEAL::Center: -0.200412 0.200412 -0.0983373, Norm: 0.300000
+DEAL::Center: -0.279462 0.109091 0.00000, Norm: 0.300000
+DEAL::Center: -0.109091 0.279462 0.00000, Norm: 0.300000
+DEAL::Center: -0.200412 0.200412 0.0983373, Norm: 0.300000
+DEAL::Center: -0.400824 0.400824 -0.196675, Norm: 0.600000
+DEAL::Center: -0.558924 0.218183 0.00000, Norm: 0.600000
+DEAL::Center: -0.218183 0.558924 0.00000, Norm: 0.600000
+DEAL::Center: -0.400824 0.400824 0.196675, Norm: 0.600000
+DEAL::Center: -0.109091 0.00000 -0.279462, Norm: 0.300000
+DEAL::Center: -0.200412 -0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: -0.200412 0.0983373 -0.200412, Norm: 0.300000
+DEAL::Center: -0.279462 0.00000 -0.109091, Norm: 0.300000
+DEAL::Center: -0.218183 0.00000 -0.558924, Norm: 0.600000
+DEAL::Center: -0.400824 -0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: -0.400824 0.196675 -0.400824, Norm: 0.600000
+DEAL::Center: -0.558924 0.00000 -0.218183, Norm: 0.600000
+DEAL::Center: 0.0983373 -0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.00000 -0.279462 -0.109091, Norm: 0.300000
+DEAL::Center: 0.00000 -0.109091 -0.279462, Norm: 0.300000
+DEAL::Center: -0.0983373 -0.200412 -0.200412, Norm: 0.300000
+DEAL::Center: 0.196675 -0.400824 -0.400824, Norm: 0.600000
+DEAL::Center: 0.00000 -0.558924 -0.218183, Norm: 0.600000
+DEAL::Center: 0.00000 -0.218183 -0.558924, Norm: 0.600000
+DEAL::Center: -0.196675 -0.400824 -0.400824, Norm: 0.600000

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.