]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Simplified interface of Manifold. Removed single argument versions of normal_vector...
authorheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 29 Jan 2014 10:44:06 +0000 (10:44 +0000)
committerheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 29 Jan 2014 10:44:06 +0000 (10:44 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_manifold_id@32341 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/grid/manifold.h
deal.II/include/deal.II/grid/tria_boundary.h
deal.II/include/deal.II/grid/tria_boundary_lib.h
deal.II/source/grid/manifold.cc
deal.II/source/grid/tria_boundary_lib.cc

index 5743e42d6c188cd218c32e4f7e68994be101fd46..1b0f147d5c991682bf4dd7c26411d7a564ebe14a 100644 (file)
@@ -92,32 +92,12 @@ public:
    *
    * Derived classes do not need to implement this function unless
    * mesh smoothing algorithms are used with a particular boundary
-   * object. The default implementation of this function calls the
-   * method with the same name, but without surrounding_points.
-   *
-   * This version of the projection function can be overwritten when
-   * only an approximation of the manifold is available, and this
-   * approximation requires surrounding points to be defined. 
+   * object.
    */
   virtual
   Point<spacedim> project_to_manifold (const std::vector<Point<spacedim> > &surrounding_points,
                                       const Point<spacedim> &candidate) const;
 
-  /**
-   * Given a point which lies close to the given manifold, it modifies
-   * it and projects it to manifold itself. This function is used in
-   * some mesh smoothing algorithms that try to move around points in
-   * order to improve the mesh quality but need to ensure that points
-   * that were on the manfold remain on the manifold.
-   *
-   * Derived classes do not need to implement this function unless
-   * mesh smoothing algorithms are used with a particular boundary
-   * object. The default implementation of this function throws an
-   * exception of type ExcPureFunctionCalled.
-   */
-  virtual
-  Point<spacedim> project_to_manifold (const Point<spacedim> &candidate) const;
-
   /**
     * Given a vector of points, return the normals to the Manifold in
     * those points. Notice that the Manifold may or may not be a
@@ -130,20 +110,6 @@ public:
     void get_normals_at_points(const std::vector<Point<spacedim> > &points,
                               std::vector<Point<spacedim> > &normals) const;
 
-    /**
-     * Compute the normal to the manifold at the given point. The
-     * default implementation of this function just calls
-     * normal_vector. The difference from this function and that one,
-     * is that one can use this function when only an approximation of
-     * the normal is required, for example when computing normals at
-     * quadrature points. Note that this function only makes sense
-     * when the manifold is a codimension one manifold.
-     */
-    virtual
-    Point<spacedim>
-    normal_vector(const std::vector<Point<spacedim> > &vertices,
-                 const Point<spacedim> &point) const;
-
     /**
      * Compute the normal to the manifold at the given point. The
      * default implementation of this function throws an
@@ -151,7 +117,9 @@ public:
      * manifold is a codimension one manifold.
      */
     virtual
-    Point<spacedim> normal_vector(const Point<spacedim> &point) const;
+    Point<spacedim>
+    normal_vector(const std::vector<Point<spacedim> > &vertices,
+                 const Point<spacedim> &point) const;
 };
 
 
@@ -205,8 +173,6 @@ public:
                  const std::vector<double> &weights) const;
   
   
-  using Manifold<spacedim>::project_to_manifold;
-  
   /**
    *  Project to FlatManifold: do nothing. Note however that this
    *  function can be overloaded by derived classes, which will then
@@ -215,9 +181,24 @@ public:
    *  this class.
    */
   virtual
-  Point<spacedim> project_to_manifold (const Point<spacedim> &candidate) const;
+  Point<spacedim> project_to_manifold (const std::vector<Point<spacedim> > &points, 
+                                      const Point<spacedim> &candidate) const;
 
 private:
+  /**
+   * The periodicity of this Manifold. Periodicity affects the way a
+   * middle point is computed. It is assumed that if two points are
+   * more than half period distant, then the distance should be
+   * computed by crossing the periodicity boundary, i.e., the average
+   * is computed by adding a full period to the sum of the two. For
+   * example, if along direction 0 we have 2*pi periodicity, then the
+   * average of (2*pi-eps) and (eps) is not pi, but 2*pi (or zero),
+   * since, on a periodic manifold, these two points are at distance
+   * 2*eps and not (2*pi-eps).
+   *
+   * A periodicity 0 along one direction means no periodicity. This is
+   * the default value for all directions.
+   */
   const Point<spacedim> periodicity;
 };
 
index 72b31d4db457cdcfd47747bb8384a221e4093aa1..882304fa536a66ed31cb7abecb5f9760d9f6336d 100644 (file)
@@ -155,14 +155,15 @@ public:
     virtual
     void get_normals_at_points(const std::vector<Point<spacedim> > &points,
                               std::vector<Point<spacedim> > &normals) const;
+  
   /**
-    * Given a point and a vector of vertices, compute an approximation
-    * of the normal to the manifold.
-    */
-    virtual
-    Point<spacedim> 
-    normal_vector(const std::vector<Point<spacedim> > &vertices,
-                 const Point<spacedim> &point) const;
+   * Given a point and a vector of vertices, compute an approximation
+   * of the normal to the manifold.
+   */
+  virtual
+  Point<spacedim> 
+  normal_vector(const std::vector<Point<spacedim> > &vertices,
+               const Point<spacedim> &point) const;
 };
 
 
index a83b3d89a6b7398bc8928e910001cf38a5ad23eb..8daf006889b8a4e651342e3fd7723a089461f979 100644 (file)
@@ -71,21 +71,24 @@ public:
   CylinderBoundary (const double           radius,
                     const Point<spacedim> &direction,
                     const Point<spacedim> &point_on_axis);
-    
+
   /**
    * Since this class is derived from FlatManifold, we need to
    * overload only the project_to_manifold function in order for the
    * boundary to function properly.
    */
     virtual Point<spacedim>
-    project_to_manifold (const Point<spacedim> &candidate) const;
+    project_to_manifold (const std::vector<Point<spacedim> > & vertices, 
+                        const Point<spacedim> &candidate) const;
+
 
   /**
    * Compute the normal to the surface. 
    */
-    virtual Point<spacedim>
-    normal_vector (const Point<spacedim> &point) const;
-    
+  virtual Point<spacedim>
+  normal_vector (const std::vector<Point<spacedim> > & vertices, 
+                const Point<spacedim> &point) const;
+  
   /**
    * Return the radius of the cylinder.
    */
@@ -186,7 +189,8 @@ public:
    * boundary to function properly.
    */
     virtual Point<dim>
-    project_to_manifold (const Point<dim> &candidate) const;
+    project_to_manifold (const std::vector<Point<dim> > & vertices, 
+                        const Point<dim> &candidate) const;
 
 
   /**
@@ -197,7 +201,8 @@ public:
    */
   virtual
   Point<dim> 
-  normal_vector (const Point<dim> &point) const;
+  normal_vector (const std::vector<Point<dim> > & vertices, 
+                const Point<dim> &point) const;
 
 protected:
   /**
@@ -266,7 +271,8 @@ public:
    */
   virtual
   Point<spacedim>
-  normal_vector (const Point<spacedim> &p) const;
+  normal_vector (const std::vector<Point<spacedim> > & vertices, 
+                const Point<spacedim> &p) const;
 
 
   /**
@@ -370,7 +376,8 @@ public:
    */
   virtual
   Point<dim>
-  normal_vector (const Point<dim> &p) const;
+  normal_vector (const std::vector<Point<dim> > & vertices, 
+                const Point<dim> &p) const;
 };
 
 
@@ -442,7 +449,8 @@ public:
    * Compute the normal to the surface. 
    */
     virtual Point<dim>
-    normal_vector (const Point<dim> &point) const;
+    normal_vector (const std::vector<Point<dim> > & vertices, 
+                  const Point<dim> &point) const;
     
 private:
   /**
@@ -495,7 +503,8 @@ public:
    * Get the normal from cartesian coordinates. This normal does not
    * have unit length.
    */
-  Point<spacedim>  normal_vector(const Point<spacedim> &p) const;
+  Point<spacedim>  normal_vector(const std::vector<Point<spacedim> > & vertices, 
+                                const Point<spacedim> &p) const;
 
 private:
   double           get_correct_angle(const double angle,const double x,const double y) const;
index dc59e03bf18347110b8671605e1ab65ce4fc5dcb..fe6243d883008ef789fdf5dabbca3ffa8a347d10 100644 (file)
@@ -30,22 +30,13 @@ Manifold<spacedim>::~Manifold ()
 {}
 
 
-
-
-template <int spacedim>
-Point<spacedim>
-Manifold<spacedim>::project_to_manifold (const Point<spacedim> &candidate) const 
-{
-  Assert (false, ExcPureFunctionCalled());
-  return candidate;
-}
-
 template <int spacedim>
 Point<spacedim>
 Manifold<spacedim>::project_to_manifold (const std::vector<Point<spacedim> > &,
                                         const Point<spacedim> &candidate) const 
-{
-  return project_to_manifold(candidate);
+{  
+  Assert (false, ExcPureFunctionCalled());
+  return candidate;
 }
 
 
@@ -56,20 +47,13 @@ Manifold<spacedim>::get_normals_at_points(const std::vector<Point<spacedim> > &p
 {
   AssertDimension(normals.size(), points.size());
   for(unsigned int i=0; i<normals.size(); ++i)
-    normals[i] = normal_vector(points[i]);
+    normals[i] = normal_vector(points, points[i]);
 }
 
 template <int spacedim>
 Point<spacedim>
 Manifold<spacedim>::normal_vector(const std::vector<Point<spacedim> > &,
                                  const Point<spacedim> &point) const
-{
-  return normal_vector(point);
-}
-
-template <int spacedim>
-Point<spacedim>
-Manifold<spacedim>::normal_vector(const Point<spacedim> &point) const
 {
   Assert (false, ExcPureFunctionCalled());
   return point;
@@ -122,7 +106,8 @@ get_new_point (const std::vector<Point<spacedim> > &surrounding_points,
 
 template <int spacedim>
 Point<spacedim> 
-FlatManifold<spacedim>::project_to_manifold (const Point<spacedim> &candidate) const 
+FlatManifold<spacedim>::project_to_manifold (const std::vector<Point<spacedim> > & vertices, 
+                                            const Point<spacedim> &candidate) const 
 {
   return candidate;
 }
index 831775a029a3c09bdc739f92d8e1e0f5183721d1..77044ede7888f0984bd06b6ba2a462930be0dee7 100644 (file)
@@ -61,7 +61,8 @@ CylinderBoundary<dim,spacedim>::get_axis_vector (const unsigned int axis)
 template <int dim, int spacedim>
 Point<spacedim>
 CylinderBoundary<dim,spacedim>::
-project_to_manifold (const Point<spacedim> &middle) const
+project_to_manifold (const std::vector<Point<spacedim> > & vertices, 
+                    const Point<spacedim> &middle) const
 {
   // The temptative point is computed from the get_new_point function
   // of the FlatManifold class, which in turns calls this class
@@ -85,7 +86,8 @@ project_to_manifold (const Point<spacedim> &middle) const
 template <int dim, int spacedim>
 Point<spacedim>
 CylinderBoundary<dim,spacedim>::
-normal_vector (const Point<spacedim> &vertex) const
+normal_vector (const std::vector<Point<spacedim> > & vertices, 
+              const Point<spacedim> &vertex) const
 {
   const Point<spacedim> vector_from_axis = (vertex-point_on_axis) -
                                           ((vertex-point_on_axis) * direction) * direction;
@@ -134,7 +136,8 @@ double ConeBoundary<dim>::get_radius (Point<dim> x) const
 template<int dim>
 Point<dim>
 ConeBoundary<dim>::
-project_to_manifold (const Point<dim> &middle) const
+project_to_manifold (const std::vector<Point<dim> > & vertices, 
+                    const Point<dim> &middle) const
 {
   const Point<dim> axis = x_1 - x_0;
   // To project it on the boundary of the cone we first compute the orthogonal
@@ -149,7 +152,8 @@ project_to_manifold (const Point<dim> &middle) const
 template<int dim>
 Point<dim>
 ConeBoundary<dim>::
-normal_vector (const Point<dim> &vertex) const 
+normal_vector (const std::vector<Point<dim> > & vertices, 
+              const Point<dim> &vertex) const 
 {
   const Point<dim> axis = x_1 - x_0;
 
@@ -270,7 +274,8 @@ HyperBallBoundary<dim,spacedim>::get_intermediate_point(const Point<spacedim> &p
 template <int dim, int spacedim>
 Point<spacedim>
 HyperBallBoundary<dim,spacedim>::
-normal_vector (const Point<spacedim> &vertex) const
+normal_vector (const std::vector<Point<spacedim> > & vertices, 
+              const Point<spacedim> &vertex) const
 {
   return vertex-center;
 }
@@ -327,7 +332,8 @@ HalfHyperBallBoundary<dim>::get_new_point(const std::vector<Point<dim> > &surrou
 template <int dim>
 Point<dim>
 HalfHyperBallBoundary<dim>::
-normal_vector (const Point<dim> &vertex) const
+normal_vector (const std::vector<Point<dim> > & vertices, 
+              const Point<dim> &vertex) const
 {
   // Check whether the point is inside the flat part, or outside. If
   // it is both at distance R from the center and with zero first
@@ -397,7 +403,8 @@ HalfHyperShellBoundary<dim>::HalfHyperShellBoundary (const Point<dim> &center,
 template <int dim>
 Point<dim>
 HalfHyperShellBoundary<dim>::
-normal_vector (const Point<dim> &vertex) const
+normal_vector (const std::vector<Point<dim> > & vertices, 
+              const Point<dim> &vertex) const
 {
   if (vertex(0) == this->center(0)) 
     {
@@ -406,7 +413,7 @@ normal_vector (const Point<dim> &vertex) const
       return n;
     }
   else
-    return HyperShellBoundary<dim>::normal_vector (vertex);
+    return HyperShellBoundary<dim>::normal_vector (vertices, vertex);
 }
 
 
@@ -523,7 +530,8 @@ TorusBoundary<2,3>:: get_surf_norm_from_sp(const Point<2> &surfP) const
 //Normal field without unit length
 template <>
 Point<3>
-TorusBoundary<2,3>::normal_vector(const Point<3> &p) const
+TorusBoundary<2,3>::normal_vector(const std::vector<Point<3> > & vertices, 
+                                 const Point<3> &p) const
 {
 
   Point<2> surfP=pull_back(p);

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.