]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added clone members to (deprecated) Boundary objects.
authorLuca Heltai <luca.heltai@sissa.it>
Sun, 8 Apr 2018 17:27:24 +0000 (19:27 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Mon, 9 Apr 2018 12:32:13 +0000 (14:32 +0200)
include/deal.II/grid/tria_boundary.h
include/deal.II/grid/tria_boundary_lib.h
source/grid/tria_boundary.cc
source/grid/tria_boundary_lib.cc
tests/grid/oned_vertex_manifold.cc
tests/manifold/chart_manifold_03_embedded.cc

index 8aa25aba8751185ae543ff3aeaf2d9fdd8a52ad1..d50da52db7b9a742db8a09b1bea72f74c5600585 100644 (file)
@@ -91,6 +91,10 @@ public:
    */
   virtual ~Boundary () = default;
 
+  /**
+   * Clone this Boundary.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
 
   /**
    * Return intermediate points on a line spaced according to the interior
@@ -249,6 +253,11 @@ public:
    */
   StraightBoundary ();
 
+  /**
+   * Clone this Boundary.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Let the new point be the arithmetic mean of the two vertices of the line.
    *
@@ -256,7 +265,7 @@ public:
    * the base class for more information.
    */
   virtual Point<spacedim>
-  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const override;
 
   /**
    * Let the new point be the arithmetic mean of the four vertices of this
@@ -268,7 +277,7 @@ public:
    */
   virtual
   Point<spacedim>
-  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const override;
 
   /**
    * Return <tt>n=points.size()</tt> points that split the StraightBoundary
@@ -280,7 +289,7 @@ public:
   virtual
   void
   get_intermediate_points_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Return <tt>n=points.size()=m*m</tt> points that splits the
@@ -292,7 +301,7 @@ public:
   virtual
   void
   get_intermediate_points_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Implementation of the function declared in the base class.
@@ -303,7 +312,7 @@ public:
   virtual
   Tensor<1,spacedim>
   normal_vector (const typename Triangulation<dim,spacedim>::face_iterator &face,
-                 const Point<spacedim> &p) const;
+                 const Point<spacedim> &p) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -314,7 +323,7 @@ public:
   virtual
   void
   get_normals_at_vertices (const typename Triangulation<dim,spacedim>::face_iterator &face,
-                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const override;
 
   /**
    * Given a candidate point and a line segment characterized by the iterator,
@@ -333,7 +342,7 @@ public:
   virtual
   Point<spacedim>
   project_to_surface (const typename Triangulation<dim,spacedim>::line_iterator &line,
-                      const Point<spacedim> &candidate) const;
+                      const Point<spacedim> &candidate) const override;
 
   /**
    * Same function as above but for a point that is to be projected onto the
@@ -349,7 +358,7 @@ public:
   virtual
   Point<spacedim>
   project_to_surface (const typename Triangulation<dim,spacedim>::quad_iterator &quad,
-                      const Point<spacedim> &candidate) const;
+                      const Point<spacedim> &candidate) const override;
 
   /**
    * Same function as above but for a point that is to be projected onto the
@@ -366,7 +375,7 @@ public:
   virtual
   Point<spacedim>
   project_to_surface (const typename Triangulation<dim,spacedim>::hex_iterator &hex,
-                      const Point<spacedim> &candidate) const;
+                      const Point<spacedim> &candidate) const override;
 };
 
 
index 321f6d0a85a301669dbe70e761dab301a7961315..ff256667a25814362c5b6629a77da009cc722c27 100644 (file)
@@ -72,19 +72,24 @@ public:
                     const Point<spacedim> &direction,
                     const Point<spacedim> &point_on_axis);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Refer to the general documentation of this class and the documentation of
    * the base class.
    */
   virtual Point<spacedim>
-  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
    * the base class.
    */
   virtual Point<spacedim>
-  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -94,7 +99,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -104,7 +109,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -114,7 +119,7 @@ public:
    */
   virtual void
   get_normals_at_vertices (const typename Triangulation<dim,spacedim>::face_iterator &face,
-                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const override;
 
   /**
    * Return the radius of the cylinder.
@@ -216,6 +221,11 @@ public:
                 const Point<dim> x_0,
                 const Point<dim> x_1);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,dim> > clone() const override;
+
   /**
    * Return the radius of the (truncated) cone at given point <tt>x</tt> on
    * the axis.
@@ -228,7 +238,7 @@ public:
    */
   virtual
   Point<dim>
-  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -236,7 +246,7 @@ public:
    */
   virtual
   Point<dim>
-  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -247,7 +257,7 @@ public:
   virtual
   void
   get_intermediate_points_on_line (const typename Triangulation<dim>::line_iterator &line,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -258,7 +268,7 @@ public:
   virtual
   void
   get_intermediate_points_on_quad (const typename Triangulation<dim>::quad_iterator &quad,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -269,7 +279,7 @@ public:
   virtual
   void
   get_normals_at_vertices (const typename Triangulation<dim>::face_iterator &face,
-                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const override;
 
   /**
    * Implementation of the function declared in the base class.
@@ -280,7 +290,7 @@ public:
   virtual
   Tensor<1,dim>
   normal_vector (const typename Triangulation<dim>::face_iterator &face,
-                 const Point<dim> &p) const;
+                 const Point<dim> &p) const override;
 
 
 protected:
@@ -351,13 +361,18 @@ public:
   HyperBallBoundary (const Point<spacedim> p      = Point<spacedim>(),
                      const double     radius = 1.0);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Refer to the general documentation of this class and the documentation of
    * the base class.
    */
   virtual
   Point<spacedim>
-  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -365,7 +380,7 @@ public:
    */
   virtual
   Point<spacedim>
-  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -376,7 +391,7 @@ public:
   virtual
   void
   get_intermediate_points_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -387,7 +402,7 @@ public:
   virtual
   void
   get_intermediate_points_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad,
-                                   std::vector<Point<spacedim> > &points) const;
+                                   std::vector<Point<spacedim> > &points) const override;
 
   /**
    * Implementation of the function declared in the base class.
@@ -398,7 +413,7 @@ public:
   virtual
   Tensor<1,spacedim>
   normal_vector (const typename Triangulation<dim,spacedim>::face_iterator &face,
-                 const Point<spacedim> &p) const;
+                 const Point<spacedim> &p) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -409,7 +424,7 @@ public:
   virtual
   void
   get_normals_at_vertices (const typename Triangulation<dim,spacedim>::face_iterator &face,
-                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const override;
 
   /**
    * Return the center of the ball.
@@ -493,17 +508,22 @@ public:
   HalfHyperBallBoundary (const Point<dim> p      = Point<dim>(),
                          const double     radius = 1.0);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,dim> > clone() const override;
+
   /**
    * Check if on the line <tt>x==0</tt>, otherwise pass to the base class.
    */
   virtual Point<dim>
-  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const override;
 
   /**
    * Check if on the line <tt>x==0</tt>, otherwise pass to the base class.
    */
   virtual Point<dim>
-  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -513,7 +533,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_line (const typename Triangulation<dim>::line_iterator &line,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -523,7 +543,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_quad (const typename Triangulation<dim>::quad_iterator &quad,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -533,7 +553,7 @@ public:
    */
   virtual void
   get_normals_at_vertices (const typename Triangulation<dim>::face_iterator &face,
-                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const override;
 };
 
 
@@ -562,6 +582,11 @@ public:
    * radius as argument. This radius will be ignored
    */
   HyperShellBoundary (const Point<dim> &center = Point<dim>());
+
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,dim> > clone() const override;
 };
 
 
@@ -597,17 +622,22 @@ public:
                           const double inner_radius = -1,
                           const double outer_radius = -1);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,dim> > clone() const override;
+
   /**
    * Construct a new point on a line.
    */
   virtual Point<dim>
-  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim>::line_iterator &line) const override;
 
   /**
    * Construct a new point on a quad.
    */
   virtual Point<dim>
-  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim>::quad_iterator &quad) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -617,7 +647,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_line (const typename Triangulation<dim>::line_iterator &line,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Refer to the general documentation of this class and the documentation of
@@ -627,7 +657,7 @@ public:
    */
   virtual void
   get_intermediate_points_on_quad (const typename Triangulation<dim>::quad_iterator &quad,
-                                   std::vector<Point<dim> > &points) const;
+                                   std::vector<Point<dim> > &points) const override;
 
   /**
    * Compute the normals to the boundary at the vertices of the given face.
@@ -637,7 +667,7 @@ public:
    */
   virtual void
   get_normals_at_vertices (const typename Triangulation<dim>::face_iterator &face,
-                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const;
+                           typename Boundary<dim>::FaceVertexNormals &face_vertex_normals) const override;
 
 private:
   /**
@@ -668,32 +698,37 @@ public:
    */
   TorusBoundary (const double R, const double r);
 
+  /**
+   * Clone this Boundary object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
 //Boundary Refinement Functions
   /**
    * Construct a new point on a line.
    */
   virtual Point<spacedim>
-  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const;
+  get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const override;
 
   /**
    * Construct a new point on a quad.
    */
   virtual Point<spacedim>
-  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const;
+  get_new_point_on_quad (const typename Triangulation<dim,spacedim>::quad_iterator &quad) const override;
 
   /**
    * Construct a new points on a line.
    */
   virtual void   get_intermediate_points_on_line (
     const typename Triangulation< dim, spacedim >::line_iterator   &line,
-    std::vector< Point< spacedim > >         &points) const;
+    std::vector< Point< spacedim > >         &points) const override;
 
   /**
    * Construct a new points on a quad.
    */
   virtual void  get_intermediate_points_on_quad (
     const typename Triangulation< dim, spacedim >::quad_iterator &quad,
-    std::vector< Point< spacedim > >         &points ) const;
+    std::vector< Point< spacedim > >         &points ) const override;
 
   /**
    * Get the normal from cartesian coordinates. This normal does not have unit
@@ -701,7 +736,7 @@ public:
    */
   virtual void get_normals_at_vertices (
     const typename Triangulation< dim, spacedim >::face_iterator &face,
-    typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const;
+    typename Boundary<dim,spacedim>::FaceVertexNormals &face_vertex_normals) const override;
 
 private:
   //Handy functions
index 82a64a2829208b2a41da885633ae74bd16bd066b..6cbab1f677029584729d12446b0bf6ff9fa8a7ae 100644 (file)
@@ -29,6 +29,16 @@ DEAL_II_NAMESPACE_OPEN
 /* -------------------------- Boundary --------------------- */
 
 
+template <int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+Boundary<dim,spacedim>::clone() const
+{
+  AssertThrow(false, ExcPureFunctionCalled());
+  return std::unique_ptr<Manifold<dim,spacedim> >();
+}
+
+
+
 template <int dim, int spacedim>
 void
 Boundary<dim, spacedim>::
@@ -194,6 +204,14 @@ get_new_point_on_line (const typename Triangulation<dim, spacedim>::line_iterato
 }
 
 
+template <int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+StraightBoundary<dim,spacedim>::clone() const
+{
+  return std_cxx14::make_unique<StraightBoundary<dim,spacedim> >();
+}
+
+
 namespace
 {
   // compute the new midpoint of a quad --
index 62b44b7783bdf763ab1c902118977ca5dd8e5522..f0864857cf463a91413b8961616fcf85a28e10ca 100644 (file)
@@ -18,6 +18,7 @@
 #include <deal.II/grid/tria_iterator.h>
 #include <deal.II/grid/tria_accessor.h>
 #include <deal.II/base/tensor.h>
+#include <deal.II/base/std_cxx14/memory.h>
 #include <cmath>
 
 
@@ -47,6 +48,16 @@ CylinderBoundary<dim,spacedim>::CylinderBoundary (const double           radius,
 {}
 
 
+
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+CylinderBoundary<dim,spacedim>::clone() const
+{
+  return std_cxx14::make_unique<CylinderBoundary<dim, spacedim> >(radius, direction, point_on_axis);
+}
+
+
+
 template <int dim, int spacedim>
 Point<spacedim>
 CylinderBoundary<dim,spacedim>::get_axis_vector (const unsigned int axis)
@@ -285,6 +296,14 @@ ConeBoundary<dim>::ConeBoundary (const double radius_0,
 
 
 
+template<int dim>
+std::unique_ptr<Manifold<dim,dim> > ConeBoundary<dim>::clone() const
+{
+  return std_cxx14::make_unique<ConeBoundary<dim> >(radius_0,radius_1,x_0,x_1);
+}
+
+
+
 template <int dim>
 double ConeBoundary<dim>::get_radius (Point<dim> x) const
 {
@@ -515,6 +534,15 @@ HyperBallBoundary<dim,spacedim>::HyperBallBoundary (const Point<spacedim> p,
 
 
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+HyperBallBoundary<dim,spacedim>::clone() const
+{
+  return std_cxx14::make_unique<HyperBallBoundary<dim,spacedim> >(center, radius);
+}
+
+
+
 template <int dim, int spacedim>
 Point<spacedim>
 HyperBallBoundary<dim,spacedim>::get_new_point_on_line (const typename Triangulation<dim,spacedim>::line_iterator &line) const
@@ -820,6 +848,15 @@ HalfHyperBallBoundary<dim>::HalfHyperBallBoundary (const Point<dim> center,
 
 
 
+template<int dim>
+std::unique_ptr<Manifold<dim,dim> > HalfHyperBallBoundary<dim>::clone() const
+{
+  return std_cxx14::make_unique<HalfHyperBallBoundary<dim> >(this->get_center(),
+                                                             this->get_radius());
+}
+
+
+
 template <int dim>
 Point<dim>
 HalfHyperBallBoundary<dim>::
@@ -960,6 +997,16 @@ HyperShellBoundary<dim>::HyperShellBoundary (const Point<dim> &center)
 }
 
 
+
+template<int dim>
+std::unique_ptr<Manifold<dim,dim> >
+HyperShellBoundary<dim>::clone() const
+{
+  return std_cxx14::make_unique<HyperShellBoundary<dim> >(this->get_center());
+}
+
+
+
 /* ---------------------------------------------------------------------- */
 
 
@@ -983,6 +1030,17 @@ HalfHyperShellBoundary<dim>::HalfHyperShellBoundary (const Point<dim> &center,
 
 
 
+template<int dim>
+std::unique_ptr<Manifold<dim,dim> >
+HalfHyperShellBoundary<dim>::clone() const
+{
+  return std_cxx14::make_unique<HalfHyperShellBoundary<dim> >(this->get_center(),
+                                                              inner_radius,
+                                                              outer_radius);
+}
+
+
+
 template <int dim>
 Point<dim>
 HalfHyperShellBoundary<dim>::
@@ -1240,6 +1298,13 @@ TorusBoundary<2,3>::TorusBoundary (const double R_,
   Assert (R>r, ExcMessage("Outer radius must be greater than inner radius."));
 }
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+TorusBoundary<dim,spacedim>::clone() const
+{
+  return std_cxx14::make_unique<TorusBoundary<dim,spacedim> >(R,r);
+}
+
 
 
 template <int dim, int spacedim>
index 79838a2fcceb2c12949f55e2e44699628dfde582..c14d2c521df52b8ff0a416dbb3e2b8ab39c58d50 100644 (file)
@@ -51,15 +51,6 @@ int main()
           deallog << "current vertex manifold id: "
                   << cell->face(vertex_n)->manifold_id()
                   << std::endl;
-
-          if (cell->face(vertex_n)->manifold_id() == 42)
-            {
-              // since the reference is const, the pointer we compare it to
-              // must also be a to a const object
-              const auto &spherical_reference = spherical_manifold;
-              Assert(&spherical_reference == &cell->face(vertex_n)->get_manifold(),
-                     ExcMessage("manifolds should be the same."));
-            }
         }
     }
 }
index f73486a43140dcc425d46da5ca8073c96ebbdf7d..c2cff33f6981c737126236653038efa739e68e04 100644 (file)
@@ -36,10 +36,15 @@ public:
     ChartManifold<dim,spacedim,spacedim+1> (periodicity)
   {}
 
+  virtual std::unique_ptr<Manifold<dim,spacedim> >
+  clone () const override
+  {
+    return std::unique_ptr<Manifold<dim,spacedim> >(new MyFlatManifold(this->get_periodicity()));
+  }
 
   virtual
   Point<spacedim+1>
-  pull_back(const Point<spacedim> &space_point) const
+  pull_back(const Point<spacedim> &space_point) const override
   {
     Point<spacedim+1> p;
     for (unsigned int d=0; d<spacedim; ++d)
@@ -50,7 +55,7 @@ public:
 
   virtual
   Point<spacedim>
-  push_forward(const Point<spacedim+1> &chart_point) const
+  push_forward(const Point<spacedim+1> &chart_point) const override
   {
     Point<spacedim> p;
     for (unsigned int d=0; d<spacedim; ++d)
@@ -60,7 +65,7 @@ public:
 
   virtual
   DerivativeForm<1,spacedim+1,spacedim>
-  push_forward_gradient(const Point<spacedim+1> &chart_point) const
+  push_forward_gradient(const Point<spacedim+1> &chart_point) const override
   {
     DerivativeForm<1,spacedim+1,spacedim> x;
     for (unsigned int d=0; d<spacedim; ++d)

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.