]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Data members of TensorProductManifold converted to std::unique_ptr<...> 7361/head
authorStefano Dominici <sfndmn@gmail.com>
Sun, 21 Oct 2018 08:10:33 +0000 (10:10 +0200)
committerStefano Dominici <sfndmn@gmail.com>
Mon, 22 Oct 2018 06:51:51 +0000 (08:51 +0200)
Mimor modifications to:
 - Remove unnecessary parts.

include/deal.II/grid/manifold_lib.h
include/deal.II/grid/tensor_product_manifold.h
tests/manifold/tensor_product_manifold_03.cc
tests/manifold/tensor_product_manifold_03.output

index a4498acaf39f6449c9e17146925f2a5a3d9783de..f6da4832029ac3b6fd3ef12da157e74c5b952f40 100644 (file)
@@ -73,11 +73,6 @@ public:
    */
   PolarManifold(const Point<spacedim> center = Point<spacedim>());
 
-  /**
-   * Virtual destructor
-   */
-  virtual ~PolarManifold() = default;
-
   /**
    * Make a clone of this Manifold object.
    */
@@ -404,10 +399,6 @@ public:
                       const Point<spacedim> &    point_on_axis,
                       const double               tolerance = 1e-10);
 
-  /**
-   * Virtual destructor.
-   */
-  virtual ~CylindricalManifold() = default;
   /**
    * Make a clone of this Manifold object.
    */
@@ -507,7 +498,6 @@ public:
     const Tensor<1, chartdim> &periodicity = Tensor<1, chartdim>(),
     const double               tolerance   = 1e-10);
 
-
   /**
    * Expressions constructor. Takes the expressions of the push_forward
    * function of spacedim components, and of the pull_back function of @p
@@ -677,11 +667,6 @@ public:
    */
   TorusManifold(const double R, const double r);
 
-  /**
-   * Virtual destructor
-   */
-  virtual ~TorusManifold() = default;
-
   /**
    * Make a clone of this Manifold object.
    */
index 73c65611325ec3a46b7027b6019de13121ec27fc..5f371c3437f041d7fd93138ba442acdd2f3ac0c7 100644 (file)
@@ -21,6 +21,7 @@
 #include <deal.II/base/point.h>
 #include <deal.II/base/std_cxx14/memory.h>
 #include <deal.II/base/subscriptor.h>
+#include <deal.II/base/utilities.h>
 
 #include <deal.II/grid/manifold.h>
 
@@ -89,11 +90,6 @@ public:
     const ChartManifold<dim_A, spacedim_A, chartdim_A> &manifold_A,
     const ChartManifold<dim_B, spacedim_B, chartdim_B> &manifold_B);
 
-  /**
-   * Virtual destructor
-   */
-  virtual ~TensorProductManifold() = default;
-
   /**
    * Clone this manifold.
    */
@@ -119,10 +115,10 @@ public:
   push_forward_gradient(const Point<chartdim> &chart_point) const override;
 
 private:
-  std::shared_ptr<const ChartManifold<dim_A, spacedim_A, chartdim_A>>
+  std::unique_ptr<const ChartManifold<dim_A, spacedim_A, chartdim_A>>
     manifold_A;
 
-  std::shared_ptr<const ChartManifold<dim_B, spacedim_B, chartdim_B>>
+  std::unique_ptr<const ChartManifold<dim_B, spacedim_B, chartdim_B>>
     manifold_B;
 };
 
@@ -196,18 +192,13 @@ TensorProductManifold<dim,
       internal::TensorProductManifoldImplementation::concat(
         manifold_A.get_periodicity(),
         manifold_B.get_periodicity()))
-{
-  std::shared_ptr<Manifold<dim_A, spacedim_A>> tmpA(
-    std::move(manifold_A.clone()));
-  std::shared_ptr<Manifold<dim_B, spacedim_B>> tmpB(
-    std::move(manifold_B.clone()));
-  this->manifold_A =
-    std::static_pointer_cast<ChartManifold<dim_A, spacedim_A, chartdim_A>>(
-      tmpA);
-  this->manifold_B =
-    std::static_pointer_cast<ChartManifold<dim_B, spacedim_B, chartdim_B>>(
-      tmpB);
-}
+  , manifold_A(Utilities::dynamic_unique_cast<
+               ChartManifold<dim_A, spacedim_A, chartdim_A>,
+               Manifold<dim_A, spacedim_A>>(manifold_A.clone()))
+  , manifold_B(Utilities::dynamic_unique_cast<
+               ChartManifold<dim_B, spacedim_B, chartdim_B>,
+               Manifold<dim_B, spacedim_B>>(manifold_B.clone()))
+{}
 
 template <int dim,
           int dim_A,
index ecfea9cccbf8d7dbe54110ca39d0358f7a3b8296..6d96b0b4af73f5ce888db4998110e98571a04483 100644 (file)
 #include "../tests.h"
 
 
-void
+bool
 test()
 {
   Triangulation<3, 3> tria;
   GridGenerator::hyper_cube(tria);
   {
-    FunctionManifold<1, 1>                     F("x", "x");
-    PolarManifold<2, 2>                        G(Point<2>(0.5, 0.5));
-    TensorProductManifold<3, 2, 2, 2, 1, 1, 1> manifold(G, F);
-    tria.set_all_manifold_ids(0);
-    tria.set_manifold(0, manifold);
+    FunctionManifold<1, 1> F("x", "x");
+    PolarManifold<2, 2>    G(Point<2>(0.5, 0.5));
+    try
+      {
+        TensorProductManifold<3, 2, 2, 2, 1, 1, 1> manifold(G, F);
+        tria.set_all_manifold_ids(0);
+        tria.set_manifold(0, manifold);
+      }
+    catch (std::bad_cast &e)
+      {
+        (void)e;
+        return false;
+      }
   }
   tria.refine_global(1);
+  return true;
 }
 
 
@@ -46,7 +55,16 @@ test()
 int
 main()
 {
-  test();
+  initlog();
+
+  if (test())
+    {
+      deallog << "OK" << std::endl;
+    }
+  else
+    {
+      deallog << "FAILED" << std::endl;
+    }
 
   return 0;
 }
index e69de29bb2d1d6434b8b29ae775ad8c2e48c5391..0fd8fc12f0b442283edd8868867114c242b04d11 100644 (file)
@@ -0,0 +1,2 @@
+
+DEAL::OK

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.