PointWrapper
get_barycenter() const;
+ /**
+ * Get the barycenter of the cell.
+ */
+ PointWrapper
+ get_center(const bool respect_manifold = false,
+ const bool interpolate_from_surrounding = false) const;
+
/**
* Set the material id.
*/
PointWrapper
transform_real_to_unit_cell(CellAccessorWrapper &cell, PointWrapper &point);
+ /**
+ * Get the underlying mapping
+ */
+ void *
+ get_mapping();
+
private:
/**
* Dimension of the underlying Mapping object.
#include <boost/python.hpp>
#include <manifold_wrapper.h>
+#include <mapping_wrapper.h>
#include <point_wrapper.h>
#include <string>
void
flatten_triangulation(TriangulationWrapper &tria_out);
+ /**
+ * Take a 2d Triangulation that is being extruded in z direction by
+ * the total height of height using n_slices slices (minimum is 2).
+ * The boundary indicators of the faces of input are going to be
+ * assigned to the corresponding side walls in z direction. The
+ * bottom and top get the next two free boundary indicators.
+ */
+ void
+ extrude_triangulation(const unsigned int n_slices,
+ const double height,
+ TriangulationWrapper &tria_out);
+
+ /**
+ * Distort the given triangulation by randomly moving around all the
+ * vertices of the grid. The direction of movement of each vertex is
+ * random, while the length of the shift vector has a value of factor
+ * times the minimal length of the active edges adjacent to this vertex.
+ * Note that factor should obviously be well below 0.5.
+ */
+ void
+ distort_random(const double factor, const bool keep_boundary = true);
+
+ /**
+ * Transform the vertices of the given triangulation by applying the
+ * function object provided as first argument to all its vertices.
+ */
+ void
+ transform(boost::python::object &transformation);
+
+ /**
+ * Find and return an active cell that surrounds a given point p.
+ * The mapping used to determine whether the given point is inside a given
+ * cell.
+ */
+ CellAccessorWrapper
+ find_active_cell_around_point(PointWrapper & p,
+ MappingQGenericWrapper &mapping);
+
/**
* Assign a manifold object to a certain part of the triangulation.
* If an object with manifold number is refined, this object
return PointWrapper(barycenter_list);
}
+ template <int dim, int spacedim>
+ PointWrapper
+ get_center(const bool respect_manifold,
+ const bool interpolate_from_surrounding,
+ const void *cell_accessor)
+ {
+ const CellAccessor<dim, spacedim> *cell =
+ static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor);
+ Point<spacedim> center =
+ cell->center(respect_manifold, interpolate_from_surrounding);
+ boost::python::list center_list;
+ for (int i = 0; i < spacedim; ++i)
+ center_list.append(center[i]);
+
+ return PointWrapper(center_list);
+ }
template <int dim, int spacedim>
+ PointWrapper
+ CellAccessorWrapper::get_center(const bool respect_manifold,
+ const bool interpolate_from_surrounding) const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::get_center<2, 2>(respect_manifold,
+ interpolate_from_surrounding,
+ cell_accessor);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::get_center<2, 3>(respect_manifold,
+ interpolate_from_surrounding,
+ cell_accessor);
+ else
+ return internal::get_center<3, 3>(respect_manifold,
+ interpolate_from_surrounding,
+ cell_accessor);
+ }
+
+
+
void
CellAccessorWrapper::set_material_id(const int material_id)
{
namespace python
{
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(get_center_overloads, get_center, 0, 2)
+
const char refine_flag_docstring[] =
"Get/Set the refine_flag of the cell. In 2D, the possibilities are: \n"
" - isotropic \n"
+ const char center_docstring[] =
+ "Center of the object. This method can take into account manifold \n"
+ "attached to a cell. \n";
+
+
+
const char set_vertex_docstring[] =
" Set the ith vertex of the cell to point_wrapper \n";
&CellAccessorWrapper::get_barycenter,
barycenter_docstring,
boost::python::args("self"))
+ .def("center",
+ &CellAccessorWrapper::get_center,
+ get_center_overloads(
+ boost::python::args("self",
+ "respect_manifold",
+ "interpolate_from_surrounding"),
+ center_docstring))
.def("set_vertex",
&CellAccessorWrapper::set_vertex,
set_vertex_docstring,
generate_hyper_shell,
3,
5)
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(distort_random_overloads,
+ distort_random,
+ 1,
+ 2)
const char n_active_cells_docstring[] =
"Return the number of active cells \n";
+ const char extrude_docstring[] =
+ "Take a 2d Triangulation that is being extruded in z direction by \n"
+ "the total height of height using n_slices slices (minimum is 2). \n"
+ "The boundary indicators of the faces of input are going to be \n"
+ "assigned to the corresponding side walls in z direction. The \n"
+ "bottom and top get the next two free boundary indicators. \n";
+
+
+
const char flatten_triangulation_docstring[] =
"Create a new flat triangulation out_tria which contains a single \n"
"level with all active cells of the input triangulation. If the spacedim\n"
+ const char distort_random_docstring[] =
+ "Distort the given triangulation by randomly moving around all the \n"
+ "vertices of the grid. The direction of movement of each vertex is \n"
+ "random, while the length of the shift vector has a value of factor \n"
+ "times the minimal length of the active edges adjacent to this vertex. \n"
+ "Note that factor should obviously be well below 0.5. \n";
+
+
+
const char refine_global_docstring[] =
"Refine all the cells times time \n";
+ const char transform_docstring[] =
+ "Transform the vertices of the given triangulation by applying the \n"
+ "function object provided as first argument to all its vertices. \n";
+
+
+
+ const char find_active_cell_around_point_docstring[] =
+ "Find and return an active cell that surrounds a given point p. \n";
+
+
+
void
export_triangulation()
{
&TriangulationWrapper::merge_triangulations,
merge_docstring,
boost::python::args("self", "triangulation_1", "triangulation_2"))
+ .def("extrude_triangulation",
+ &TriangulationWrapper::extrude_triangulation,
+ extrude_docstring,
+ boost::python::args("self", "n_slices", "depth", "tria_out"))
.def("flatten_triangulation",
&TriangulationWrapper::flatten_triangulation,
flatten_triangulation_docstring,
boost::python::args("self", "tria_out"))
+ .def("distort_random",
+ &TriangulationWrapper::distort_random,
+ distort_random_overloads(
+ boost::python::args("self", "factor", "keep_boundary"),
+ distort_random_docstring))
+ .def("transform",
+ &TriangulationWrapper::transform,
+ transform_docstring,
+ boost::python::args("self", "transformation"))
+ .def("find_active_cell_around_point",
+ &TriangulationWrapper::find_active_cell_around_point,
+ find_active_cell_around_point_docstring,
+ boost::python::args("self", "point", "mapping"))
.def("refine_global",
&TriangulationWrapper::refine_global,
refine_global_docstring,
p.point);
}
+
+
+ void *
+ MappingQGenericWrapper::get_mapping()
+ {
+ return mapping_ptr;
+ }
+
} // namespace python
DEAL_II_NAMESPACE_CLOSE
+ template <int dim, int spacedim>
+ void
+ distort_random(const double factor,
+ const bool keep_boundary,
+ void * triangulation)
+ {
+ Triangulation<dim, spacedim> *tria =
+ static_cast<Triangulation<dim, spacedim> *>(triangulation);
+
+ GridTools::distort_random(factor, *tria, keep_boundary);
+ }
+
+
+
+ template <int dim, int spacedim>
+ void
+ transform(boost::python::object &transformation, void *triangulation)
+ {
+ Triangulation<dim, spacedim> *tria =
+ static_cast<Triangulation<dim, spacedim> *>(triangulation);
+
+ auto func = [&transformation](const Point<spacedim> &p_in) {
+ boost::python::list p_list_in;
+
+ for (int d = 0; d < spacedim; ++d)
+ p_list_in.append(p_in[d]);
+
+ boost::python::list p_list_out =
+ boost::python::extract<boost::python::list>(
+ transformation(p_list_in));
+
+ Point<spacedim> p_out;
+
+ for (int d = 0; d < spacedim; ++d)
+ p_out[d] = boost::python::extract<double>(p_list_out[d]);
+
+ return p_out;
+ };
+
+ GridTools::transform(func, *tria);
+ }
+
+
+
+ template <int dim, int spacedim>
+ std::pair<int, int>
+ find_active_cell_around_point(PointWrapper & p,
+ MappingQGenericWrapper &mapping_wrapper,
+ void * triangulation)
+ {
+ Triangulation<dim, spacedim> *tria =
+ static_cast<Triangulation<dim, spacedim> *>(triangulation);
+
+ Point<spacedim> point = *(static_cast<Point<spacedim> *>(p.get_point()));
+
+ const MappingQGeneric<dim, spacedim> *mapping =
+ static_cast<const MappingQGeneric<dim, spacedim> *>(
+ mapping_wrapper.get_mapping());
+
+ auto cell_pair =
+ GridTools::find_active_cell_around_point(*mapping, *tria, point);
+
+ return std::make_pair(cell_pair.first->level(), cell_pair.first->index());
+ }
+
+
+
template <int dim, int spacedim>
boost::python::list
active_cells(TriangulationWrapper &triangulation_wrapper)
+ void
+ TriangulationWrapper::extrude_triangulation(
+ const unsigned int n_slices,
+ const double height,
+ TriangulationWrapper &triangulation_out)
+ {
+ AssertThrow((dim == 2) && (spacedim == 2),
+ ExcMessage(
+ "Extrude can only be applied to the dim and spacedim two."));
+ AssertThrow((triangulation_out.get_dim() == 3) &&
+ (triangulation_out.get_spacedim() == 3),
+ ExcMessage(
+ "The output Triangulation must be of dimension three"));
+
+ if ((dim == 2) && (spacedim == 2))
+ {
+ Triangulation<2, 2> *tria =
+ static_cast<Triangulation<2, 2> *>(triangulation);
+ Triangulation<3, 3> *tria_out = static_cast<Triangulation<3, 3> *>(
+ triangulation_out.get_triangulation());
+ GridGenerator::extrude_triangulation(*tria,
+ n_slices,
+ height,
+ *tria_out);
+ }
+ }
+
+
+
+ void
+ TriangulationWrapper::distort_random(const double factor,
+ const bool keep_boundary)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ internal::distort_random<2, 2>(factor, keep_boundary, triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ internal::distort_random<2, 3>(factor, keep_boundary, triangulation);
+ else
+ internal::distort_random<3, 3>(factor, keep_boundary, triangulation);
+ }
+
+
+
+ void
+ TriangulationWrapper::transform(boost::python::object &transformation)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ internal::transform<2, 2>(transformation, triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ internal::transform<2, 3>(transformation, triangulation);
+ else
+ internal::transform<3, 3>(transformation, triangulation);
+ }
+
+
+
+ CellAccessorWrapper
+ TriangulationWrapper::find_active_cell_around_point(
+ PointWrapper & p,
+ MappingQGenericWrapper &mapping)
+ {
+ std::pair<int, int> level_index_pair;
+ if ((dim == 2) && (spacedim == 2))
+ level_index_pair =
+ internal::find_active_cell_around_point<2, 2>(p,
+ mapping,
+ triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ level_index_pair =
+ internal::find_active_cell_around_point<2, 3>(p,
+ mapping,
+ triangulation);
+ else
+ level_index_pair =
+ internal::find_active_cell_around_point<3, 3>(p,
+ mapping,
+ triangulation);
+
+ return CellAccessorWrapper(*this,
+ level_index_pair.first,
+ level_index_pair.second);
+ }
+
+
+
void
TriangulationWrapper::refine_global(const unsigned int n)
{
--- /dev/null
+Added: Python wrappers for the CellAccessor::center, GridTools::transform,
+GridTools::distort_random, GridGenerator::extrude_triangulation,
+GridTools::find_active_cell_around_point
+<br>
+(Alexander Grayver, 2019/12/06)