--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_manifold_wrapper_h
+#define dealii_manifold_wrapper_h
+
+#include <deal.II/base/config.h>
+
+#include <boost/python.hpp>
+
+#include <point_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ class ManifoldWrapper
+ {
+ public:
+ /**
+ * Copy constructor.
+ */
+ ManifoldWrapper(const ManifoldWrapper &other);
+
+ /**
+ * Constructor.
+ */
+ ManifoldWrapper(const int dim, const int spacedim);
+
+ /**
+ * Destructor.
+ */
+ ~ManifoldWrapper();
+
+ /**
+ * Create SphericalManifold
+ */
+ void
+ create_spherical(const PointWrapper center);
+
+ /**
+ * Create PolarManifold
+ */
+ void
+ create_polar(const PointWrapper center);
+
+ /**
+ * Create CylindricalManifold
+ */
+ void
+ create_cylindrical(const int axis, const double tolerance);
+
+ /**
+ * Return pointer to an underlying manifold object
+ */
+ void *
+ get_manifold() const;
+
+
+ private:
+ /**
+ * Dimension of the underlying Manifold object.
+ */
+ int dim;
+
+ /**
+ * Space dimension of the underlying Manifold object.
+ */
+ int spacedim;
+
+ /**
+ * Pointer that can be casted to the underlying Mapping object.
+ */
+ void *manifold_ptr;
+ };
+
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
#include <boost/python.hpp>
+#include <manifold_wrapper.h>
#include <point_wrapper.h>
#include <string>
void
flatten_triangulation(TriangulationWrapper &tria_out);
+ /**
+ *
+ */
+ void
+ set_manifold(const int number, ManifoldWrapper &manifold);
+
/**
* Refine all the cells @p n times.
*/
export_point.cc
export_triangulation.cc
export_mapping.cc
+ export_manifold.cc
cell_accessor_wrapper.cc
tria_accessor_wrapper.cc
point_wrapper.cc
triangulation_wrapper.cc
mapping_wrapper.cc
+ manifold_wrapper.cc
)
FOREACH(_build ${DEAL_II_BUILD_TYPES})
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <boost/python.hpp>
+
+#include <manifold_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ const char create_spherical_docstring[] =
+ " Create spherical manifold with a given center point \n";
+
+
+ const char create_polar_docstring[] =
+ " Create polar manifold with a given center point \n";
+
+
+ const char create_cylindrical_docstring[] =
+ " Create cylindrical manifold along a given axis \n"
+ " (0 - x, 1 - y, 2 - z) \n";
+
+
+
+ void
+ export_manifold()
+ {
+ boost::python::class_<ManifoldWrapper>(
+ "ManifoldWrapper",
+ boost::python::init<const int, const int>(
+ boost::python::args("dim", "spacedim")))
+ .def("create_spherical",
+ &ManifoldWrapper::create_spherical,
+ create_spherical_docstring,
+ boost::python::args("self", "center"))
+ .def("create_polar",
+ &ManifoldWrapper::create_polar,
+ create_polar_docstring,
+ boost::python::args("self", "center"))
+ .def("create_cylindrical",
+ &ManifoldWrapper::create_cylindrical,
+ create_cylindrical_docstring,
+ boost::python::args("self", "axis", "tolerance"));
+ }
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
+ const char set_manifold_docstring[] =
+ "Assign a manifold object to a certain part of the triangulation. \n"
+ "The manifold_object is not copied and MUST persist until the \n"
+ "triangulation is destroyed. \n";
+
+
+
void
export_triangulation()
{
.def("load",
&TriangulationWrapper::load,
load_docstring,
- boost::python::args("self", "filename"));
+ boost::python::args("self", "filename"))
+ .def("set_manifold",
+ &TriangulationWrapper::set_manifold,
+ set_manifold_docstring,
+ boost::python::args("self", "number", "manifold"));
}
} // namespace python
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/grid/manifold_lib.h>
+
+#include <manifold_wrapper.h>
+#include <point_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ namespace internal
+ {
+ template <int dim, int spacedim>
+ Manifold<dim, spacedim> *
+ create_spherical_manifold(PointWrapper center)
+ {
+ const Point<spacedim> *point =
+ static_cast<const Point<spacedim> *>(center.get_point());
+ return new SphericalManifold<dim, spacedim>(*point);
+ }
+
+
+
+ template <int dim, int spacedim>
+ Manifold<dim, spacedim> *
+ create_polar_manifold(PointWrapper center)
+ {
+ const Point<spacedim> *point =
+ static_cast<const Point<spacedim> *>(center.get_point());
+ return new PolarManifold<dim, spacedim>(*point);
+ }
+
+
+
+ template <int dim, int spacedim>
+ Manifold<dim, spacedim> *
+ create_cylindrical_manifold(const int axis, const double tolerance)
+ {
+ return new CylindricalManifold<dim, spacedim>(axis, tolerance);
+ }
+
+
+
+ template <int dim, int spacedim>
+ Manifold<dim, spacedim> *
+ clone(void *manifold_ptr)
+ {
+ const Manifold<dim, spacedim> *manifold =
+ static_cast<const Manifold<dim, spacedim> *>(manifold_ptr);
+
+ if (const SphericalManifold<dim, spacedim> *d =
+ dynamic_cast<const SphericalManifold<dim, spacedim> *>(manifold))
+ {
+ return new SphericalManifold<dim, spacedim>(*d);
+ }
+ else if (const PolarManifold<dim, spacedim> *d =
+ dynamic_cast<const PolarManifold<dim, spacedim> *>(manifold))
+ {
+ return new PolarManifold<dim, spacedim>(*d);
+ }
+ else
+ ExcMessage("Unsupported manifold type in clone.");
+
+ return nullptr;
+ }
+ } // namespace internal
+
+
+
+ ManifoldWrapper::ManifoldWrapper(const int dim, const int spacedim)
+ : dim(dim)
+ , spacedim(spacedim)
+ {
+ AssertThrow(((dim == 2) && (spacedim == 2)) ||
+ ((dim == 2) && (spacedim == 3)) ||
+ ((dim == 3) && (spacedim == 3)),
+ ExcMessage("Wrong dim-spacedim combination."));
+ }
+
+
+
+ ManifoldWrapper::ManifoldWrapper(const ManifoldWrapper &other)
+ {
+ dim = other.dim;
+ spacedim = other.spacedim;
+
+ AssertThrow(other.manifold_ptr != nullptr,
+ ExcMessage("Underlying manifold does not exist."));
+
+ if ((dim == 2) && (spacedim == 2))
+ {
+ manifold_ptr = internal::clone<2, 2>(other.manifold_ptr);
+ }
+ else if ((dim == 2) && (spacedim == 3))
+ {
+ manifold_ptr = internal::clone<2, 3>(other.manifold_ptr);
+ }
+ else if ((dim == 3) && (spacedim == 3))
+ {
+ manifold_ptr = internal::clone<3, 3>(other.manifold_ptr);
+ }
+ else
+ AssertThrow(false, ExcMessage("Wrong dim-spacedim combination."));
+ }
+
+
+
+ ManifoldWrapper::~ManifoldWrapper()
+ {
+ if (dim != -1)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ {
+ // We cannot call delete on a void pointer so cast the void pointer
+ // back first.
+ Manifold<2, 2> *tmp = static_cast<Manifold<2, 2> *>(manifold_ptr);
+ delete tmp;
+ }
+ else if ((dim == 2) && (spacedim == 3))
+ {
+ Manifold<2, 3> *tmp = static_cast<Manifold<2, 3> *>(manifold_ptr);
+ delete tmp;
+ }
+ else
+ {
+ Manifold<3, 3> *tmp = static_cast<Manifold<3, 3> *>(manifold_ptr);
+ delete tmp;
+ }
+
+ dim = -1;
+ spacedim = -1;
+ manifold_ptr = nullptr;
+ }
+ }
+
+
+
+ void
+ ManifoldWrapper::create_spherical(const PointWrapper center)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ {
+ manifold_ptr = internal::create_spherical_manifold<2, 2>(center);
+ }
+ else if ((dim == 2) && (spacedim == 3))
+ {
+ manifold_ptr = internal::create_spherical_manifold<2, 3>(center);
+ }
+ else if ((dim == 3) && (spacedim == 3))
+ {
+ manifold_ptr = internal::create_spherical_manifold<3, 3>(center);
+ }
+ else
+ AssertThrow(false, ExcMessage("Wrong dim-spacedim combination."));
+ }
+
+
+
+ void
+ ManifoldWrapper::create_polar(const PointWrapper center)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ {
+ manifold_ptr = internal::create_polar_manifold<2, 2>(center);
+ }
+ else if ((dim == 2) && (spacedim == 3))
+ {
+ manifold_ptr = internal::create_polar_manifold<2, 3>(center);
+ }
+ else if ((dim == 3) && (spacedim == 3))
+ {
+ manifold_ptr = internal::create_polar_manifold<3, 3>(center);
+ }
+ else
+ AssertThrow(false, ExcMessage("Wrong dim-spacedim combination."));
+ }
+
+
+
+ void
+ ManifoldWrapper::create_cylindrical(const int axis, const double tolerance)
+ {
+ if ((dim == 2) && (spacedim == 3))
+ {
+ manifold_ptr =
+ internal::create_cylindrical_manifold<2, 3>(axis, tolerance);
+ }
+ else if ((dim == 3) && (spacedim == 3))
+ {
+ manifold_ptr =
+ internal::create_cylindrical_manifold<3, 3>(axis, tolerance);
+ }
+ else
+ AssertThrow(false, ExcMessage("Wrong dim-spacedim combination."));
+ }
+
+
+ void *
+ ManifoldWrapper::get_manifold() const
+ {
+ return manifold_ptr;
+ }
+
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
}
+
+ void
+ TriangulationWrapper::set_manifold(const int number,
+ ManifoldWrapper &manifold)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ {
+ Triangulation<2, 2> *tria =
+ static_cast<Triangulation<2, 2> *>(triangulation);
+ Manifold<2, 2> *m =
+ static_cast<Manifold<2, 2> *>(manifold.get_manifold());
+ tria->set_manifold(number, *m);
+ }
+ else if ((dim == 2) && (spacedim == 3))
+ {
+ Triangulation<2, 3> *tria =
+ static_cast<Triangulation<2, 3> *>(triangulation);
+ Manifold<2, 3> *m =
+ static_cast<Manifold<2, 3> *>(manifold.get_manifold());
+ tria->set_manifold(number, *m);
+ }
+ else
+ {
+ Triangulation<3, 3> *tria =
+ static_cast<Triangulation<3, 3> *>(triangulation);
+ Manifold<3, 3> *m =
+ static_cast<Manifold<3, 3> *>(manifold.get_manifold());
+ tria->set_manifold(number, *m);
+ }
+ }
+
+
+
void
TriangulationWrapper::setup(const std::string &dimension,
const std::string &spacedimension)
export_triangulation();
void
export_mapping();
+ void
+ export_manifold();
} // namespace python
DEAL_II_NAMESPACE_CLOSE
dealii::python::export_point();
dealii::python::export_triangulation();
dealii::python::export_mapping();
+ dealii::python::export_manifold();
}
#else
dealii::python::export_point();
dealii::python::export_triangulation();
dealii::python::export_mapping();
+ dealii::python::export_manifold();
}
#endif
--- /dev/null
+Added: New python wrappers for the Manifold class
+with a possibility to assign manifold object to a
+python triangulation object.
+<br>
+(Alexander Grayver, 2019/11/28)