]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement clone for each standard Manifold class.
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 6 Apr 2018 13:04:08 +0000 (15:04 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Mon, 9 Apr 2018 12:32:12 +0000 (14:32 +0200)
include/deal.II/grid/manifold_lib.h
source/grid/manifold_lib.cc

index 9bfc7b3f6612d6712a89ef4be4ce2cee69a05b20..4fd7cf0b6109cd87d569d0a45ca3f70491c54dae 100644 (file)
@@ -71,13 +71,18 @@ public:
    */
   PolarManifold(const Point<spacedim> center = Point<spacedim>());
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Pull back the given point from the Euclidean space. Will return the polar
    * coordinates associated with the point @p space_point. Only used when
    * spacedim = 2.
    */
   virtual Point<spacedim>
-  pull_back(const Point<spacedim> &space_point) const;
+  pull_back(const Point<spacedim> &space_point) const override;
 
   /**
    * Given a point in the spherical coordinate system, this method returns the
@@ -85,7 +90,7 @@ public:
    * Only used when spacedim = 3.
    */
   virtual Point<spacedim>
-  push_forward(const Point<spacedim> &chart_point) const;
+  push_forward(const Point<spacedim> &chart_point) const override;
 
   /**
    * Given a point in the spacedim dimensional Euclidean space, this
@@ -101,7 +106,7 @@ public:
    */
   virtual
   DerivativeForm<1,spacedim,spacedim>
-  push_forward_gradient(const Point<spacedim> &chart_point) const;
+  push_forward_gradient(const Point<spacedim> &chart_point) const override;
 
   /**
    * The center of the spherical coordinate system.
@@ -220,6 +225,11 @@ public:
    */
   SphericalManifold(const Point<spacedim> center = Point<spacedim>());
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Given any two points in space, first project them on the surface
    * of a sphere with unit radius, then connect them with a geodesic
@@ -388,6 +398,11 @@ public:
                        const Point<spacedim> &point_on_axis,
                        const double tolerance = 1e-10);
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Compute the Cartesian coordinates for a point given in cylindrical
    * coordinates.
@@ -509,13 +524,18 @@ public:
    */
   ~FunctionManifold();
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Given a point in the @p chartdim coordinate system, uses the
    * push_forward_function to compute the push_forward of points in @p
    * chartdim space dimensions to @p spacedim space dimensions.
    */
   virtual Point<spacedim>
-  push_forward(const Point<chartdim> &chart_point) const;
+  push_forward(const Point<chartdim> &chart_point) const override;
 
   /**
    * Given a point in the chartdim dimensional Euclidean space, this
@@ -539,7 +559,7 @@ public:
    */
   virtual
   DerivativeForm<1,chartdim,spacedim>
-  push_forward_gradient(const Point<chartdim> &chart_point) const;
+  push_forward_gradient(const Point<chartdim> &chart_point) const override;
 
   /**
    * Given a point in the spacedim coordinate system, uses the
@@ -547,7 +567,7 @@ public:
    * space dimensions to @p chartdim space dimensions.
    */
   virtual Point<chartdim>
-  pull_back(const Point<spacedim> &space_point) const;
+  pull_back(const Point<spacedim> &space_point) const override;
 
 private:
   /**
@@ -582,10 +602,36 @@ private:
    * pointers.
    */
   const bool owns_pointers;
+
+  /**
+   * The expresssion used to construct the push_forward function.
+   */
+  const std::string push_forward_expression;
+
+  /**
+   * The expresssion used to construct the pull_back function.
+   */
+  const std::string pull_back_expression;
+
+  /**
+   * Variable names in the chart domain.
+   */
+  const std::string chart_vars;
+
+  /**
+   * Variable names in the space domain.
+   */
+  const std::string space_vars;
+
+  /**
+   * The finite difference step to use internally.
+   */
+  const double finite_difference_step;
 };
 
 
 
+
 /**
  * Manifold description for the surface of a Torus in three dimensions. The
  * Torus is assumed to be in the x-z plane. The reference coordinate system
@@ -614,24 +660,29 @@ public:
    */
   TorusManifold (const double R, const double r);
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim, 3> > clone() const override;
+
   /**
    * Pull back operation.
    */
   virtual Point<3>
-  pull_back(const Point<3> &p) const;
+  pull_back(const Point<3> &p) const override;
 
   /**
    * Push forward operation.
    */
   virtual Point<3>
-  push_forward(const Point<3> &chart_point) const;
+  push_forward(const Point<3> &chart_point) const override;
 
   /**
    * Gradient.
    */
   virtual
   DerivativeForm<1,3,3>
-  push_forward_gradient(const Point<3> &chart_point) const;
+  push_forward_gradient(const Point<3> &chart_point) const override;
 
 private:
   double r, R;
@@ -781,6 +832,11 @@ public:
    */
   TransfiniteInterpolationManifold();
 
+  /**
+   * Make a clone of this Manifold object.
+   */
+  virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
   /**
    * Initializes the manifold with a coarse mesh. The prerequisite for using
    * this class is that the input triangulation is uniformly refined and the
index 06547a4d5490dbbe611a0c9cf2b1a321127bc4e8..36dff043d27a8efb0397bf6ad327fb43c6c4ebbd 100644 (file)
@@ -122,6 +122,15 @@ PolarManifold<dim,spacedim>::PolarManifold(const Point<spacedim> center):
 
 
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+PolarManifold<dim,spacedim>::clone() const
+{
+  return std::unique_ptr<Manifold<dim,spacedim> >(new PolarManifold<dim,spacedim>(center));
+}
+
+
+
 template <int dim, int spacedim>
 Tensor<1,spacedim>
 PolarManifold<dim,spacedim>::get_periodicity()
@@ -268,6 +277,15 @@ SphericalManifold<dim,spacedim>::SphericalManifold(const Point<spacedim> center)
 
 
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+SphericalManifold<dim,spacedim>::clone() const
+{
+  return std::unique_ptr<Manifold<dim,spacedim> >(new SphericalManifold<dim,spacedim>(center));
+}
+
+
+
 template <int dim, int spacedim>
 Point<spacedim>
 SphericalManifold<dim,spacedim>::
@@ -948,6 +966,16 @@ CylindricalManifold<dim, spacedim>::CylindricalManifold(const Tensor<1, spacedim
 
 
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+CylindricalManifold<dim,spacedim>::clone() const
+{
+  return std::unique_ptr<Manifold<dim,spacedim> >
+         (new CylindricalManifold<dim,spacedim>(direction, point_on_axis, tolerance));
+}
+
+
+
 template <int dim, int spacedim>
 Point<spacedim>
 CylindricalManifold<dim,spacedim>::
@@ -1077,7 +1105,8 @@ FunctionManifold<dim,spacedim,chartdim>::FunctionManifold
   push_forward_function(&push_forward_function),
   pull_back_function(&pull_back_function),
   tolerance(tolerance),
-  owns_pointers(false)
+  owns_pointers(false),
+  finite_difference_step(0)
 {
   AssertDimension(push_forward_function.n_components, spacedim);
   AssertDimension(pull_back_function.n_components, chartdim);
@@ -1098,7 +1127,12 @@ FunctionManifold<dim,spacedim,chartdim>::FunctionManifold
   ChartManifold<dim,spacedim,chartdim>(periodicity),
   const_map(const_map),
   tolerance(tolerance),
-  owns_pointers(true)
+  owns_pointers(true),
+  push_forward_expression(push_forward_expression),
+  pull_back_expression(pull_back_expression),
+  chart_vars(chart_vars),
+  space_vars(space_vars),
+  finite_difference_step(h)
 {
   FunctionParser<chartdim> *pf = new FunctionParser<chartdim>(spacedim, 0.0, h);
   FunctionParser<spacedim> *pb = new FunctionParser<spacedim>(chartdim, 0.0, h);
@@ -1127,6 +1161,32 @@ FunctionManifold<dim,spacedim,chartdim>::~FunctionManifold()
 
 
 
+template<int dim, int spacedim, int chartdim>
+std::unique_ptr<Manifold<dim, spacedim> >
+FunctionManifold<dim,spacedim,chartdim>::clone() const
+{
+  if (owns_pointers == true)
+    {
+      return std::unique_ptr<Manifold<dim,spacedim> >
+             (new FunctionManifold<dim,spacedim,chartdim>(push_forward_expression,
+                                                          pull_back_expression,
+                                                          this->get_periodicity(),
+                                                          const_map,
+                                                          chart_vars,
+                                                          space_vars,
+                                                          tolerance,
+                                                          finite_difference_step));
+    }
+  else
+    return std::unique_ptr<Manifold<dim,spacedim> >
+           (new FunctionManifold<dim,spacedim,chartdim>(*push_forward_function,
+                                                        *pull_back_function,
+                                                        this->get_periodicity(),
+                                                        tolerance));
+}
+
+
+
 template <int dim, int spacedim, int chartdim>
 Point<spacedim>
 FunctionManifold<dim,spacedim,chartdim>::push_forward(const Point<chartdim> &chart_point) const
@@ -1229,6 +1289,15 @@ TorusManifold<dim>::TorusManifold (const double R, const double r)
 
 
 
+template<int dim>
+std::unique_ptr<Manifold<dim, 3> >
+TorusManifold<dim>::clone() const
+{
+  return std::unique_ptr<Manifold<dim,3> >(new TorusManifold<dim>(R,r));
+}
+
+
+
 template <int dim>
 DerivativeForm<1,3,3>
 TorusManifold<dim>::push_forward_gradient(const Point<3> &chart_point) const
@@ -1270,6 +1339,18 @@ TransfiniteInterpolationManifold<dim,spacedim>::TransfiniteInterpolationManifold
 
 
 
+template<int dim, int spacedim>
+std::unique_ptr<Manifold<dim, spacedim> >
+TransfiniteInterpolationManifold<dim,spacedim>::clone() const
+{
+  auto ptr = new TransfiniteInterpolationManifold<dim,spacedim>();
+  if (triangulation)
+    ptr->initialize(*triangulation);
+  return std::unique_ptr<Manifold<dim,spacedim> >(ptr);
+}
+
+
+
 template <int dim, int spacedim>
 void
 TransfiniteInterpolationManifold<dim,spacedim>

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.