*/
PolarManifold(const Point<spacedim> center = Point<spacedim>());
- /**
- * Virtual destructor
- */
- virtual ~PolarManifold() = default;
-
/**
* Make a clone of this Manifold object.
*/
const Point<spacedim> & point_on_axis,
const double tolerance = 1e-10);
- /**
- * Virtual destructor.
- */
- virtual ~CylindricalManifold() = default;
/**
* Make a clone of this Manifold object.
*/
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
*/
TorusManifold(const double R, const double r);
- /**
- * Virtual destructor
- */
- virtual ~TorusManifold() = default;
-
/**
* Make a clone of this Manifold object.
*/
#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>
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.
*/
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;
};
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,
#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;
}
int
main()
{
- test();
+ initlog();
+
+ if (test())
+ {
+ deallog << "OK" << std::endl;
+ }
+ else
+ {
+ deallog << "FAILED" << std::endl;
+ }
return 0;
}