From: Alexander Grayver Date: Fri, 6 Dec 2019 15:04:59 +0000 (+0100) Subject: Add a series of new python wrappers X-Git-Tag: v9.2.0-rc1~781^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=8fb990a2be431504ed84b9677de40c5c6c38042a;p=dealii.git Add a series of new python wrappers --- diff --git a/contrib/python-bindings/include/cell_accessor_wrapper.h b/contrib/python-bindings/include/cell_accessor_wrapper.h index fd2d2bb284..0c7f471a06 100644 --- a/contrib/python-bindings/include/cell_accessor_wrapper.h +++ b/contrib/python-bindings/include/cell_accessor_wrapper.h @@ -89,6 +89,13 @@ namespace python 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. */ diff --git a/contrib/python-bindings/include/mapping_wrapper.h b/contrib/python-bindings/include/mapping_wrapper.h index b80752ce14..a51d437b64 100644 --- a/contrib/python-bindings/include/mapping_wrapper.h +++ b/contrib/python-bindings/include/mapping_wrapper.h @@ -61,6 +61,12 @@ namespace python PointWrapper transform_real_to_unit_cell(CellAccessorWrapper &cell, PointWrapper &point); + /** + * Get the underlying mapping + */ + void * + get_mapping(); + private: /** * Dimension of the underlying Mapping object. diff --git a/contrib/python-bindings/include/triangulation_wrapper.h b/contrib/python-bindings/include/triangulation_wrapper.h index 6910143fc7..63f28e0d7f 100644 --- a/contrib/python-bindings/include/triangulation_wrapper.h +++ b/contrib/python-bindings/include/triangulation_wrapper.h @@ -21,6 +21,7 @@ #include #include +#include #include #include @@ -316,6 +317,44 @@ namespace python 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 diff --git a/contrib/python-bindings/source/cell_accessor_wrapper.cc b/contrib/python-bindings/source/cell_accessor_wrapper.cc index 73991421da..d8ede1bafb 100644 --- a/contrib/python-bindings/source/cell_accessor_wrapper.cc +++ b/contrib/python-bindings/source/cell_accessor_wrapper.cc @@ -173,6 +173,22 @@ namespace python return PointWrapper(barycenter_list); } + template + PointWrapper + get_center(const bool respect_manifold, + const bool interpolate_from_surrounding, + const void *cell_accessor) + { + const CellAccessor *cell = + static_cast *>(cell_accessor); + Point 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 @@ -496,6 +512,26 @@ namespace python + 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) { diff --git a/contrib/python-bindings/source/export_cell_accessor.cc b/contrib/python-bindings/source/export_cell_accessor.cc index 877da21b98..6e088447db 100644 --- a/contrib/python-bindings/source/export_cell_accessor.cc +++ b/contrib/python-bindings/source/export_cell_accessor.cc @@ -22,6 +22,8 @@ DEAL_II_NAMESPACE_OPEN 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" @@ -68,6 +70,12 @@ namespace python + 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"; @@ -138,6 +146,13 @@ namespace python &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, diff --git a/contrib/python-bindings/source/export_triangulation.cc b/contrib/python-bindings/source/export_triangulation.cc index 5919562899..65835a57bb 100644 --- a/contrib/python-bindings/source/export_triangulation.cc +++ b/contrib/python-bindings/source/export_triangulation.cc @@ -97,6 +97,10 @@ namespace python 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"; @@ -283,6 +287,15 @@ namespace python + 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" @@ -303,6 +316,15 @@ namespace python + 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"; @@ -365,6 +387,17 @@ namespace python + 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() { @@ -484,10 +517,27 @@ namespace python &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, diff --git a/contrib/python-bindings/source/mapping_wrapper.cc b/contrib/python-bindings/source/mapping_wrapper.cc index 0c7ad75c67..97c317fe24 100644 --- a/contrib/python-bindings/source/mapping_wrapper.cc +++ b/contrib/python-bindings/source/mapping_wrapper.cc @@ -222,6 +222,14 @@ namespace python p.point); } + + + void * + MappingQGenericWrapper::get_mapping() + { + return mapping_ptr; + } + } // namespace python DEAL_II_NAMESPACE_CLOSE diff --git a/contrib/python-bindings/source/triangulation_wrapper.cc b/contrib/python-bindings/source/triangulation_wrapper.cc index 1123f7e42b..a9adb51806 100644 --- a/contrib/python-bindings/source/triangulation_wrapper.cc +++ b/contrib/python-bindings/source/triangulation_wrapper.cc @@ -484,6 +484,73 @@ namespace python + template + void + distort_random(const double factor, + const bool keep_boundary, + void * triangulation) + { + Triangulation *tria = + static_cast *>(triangulation); + + GridTools::distort_random(factor, *tria, keep_boundary); + } + + + + template + void + transform(boost::python::object &transformation, void *triangulation) + { + Triangulation *tria = + static_cast *>(triangulation); + + auto func = [&transformation](const Point &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( + transformation(p_list_in)); + + Point p_out; + + for (int d = 0; d < spacedim; ++d) + p_out[d] = boost::python::extract(p_list_out[d]); + + return p_out; + }; + + GridTools::transform(func, *tria); + } + + + + template + std::pair + find_active_cell_around_point(PointWrapper & p, + MappingQGenericWrapper &mapping_wrapper, + void * triangulation) + { + Triangulation *tria = + static_cast *>(triangulation); + + Point point = *(static_cast *>(p.get_point())); + + const MappingQGeneric *mapping = + static_cast *>( + 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 boost::python::list active_cells(TriangulationWrapper &triangulation_wrapper) @@ -1199,6 +1266,91 @@ namespace python + 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); + Triangulation<3, 3> *tria_out = static_cast *>( + 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 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) { diff --git a/doc/news/changes/minor/20191206Grayver b/doc/news/changes/minor/20191206Grayver new file mode 100644 index 0000000000..4b1de45787 --- /dev/null +++ b/doc/news/changes/minor/20191206Grayver @@ -0,0 +1,5 @@ +Added: Python wrappers for the CellAccessor::center, GridTools::transform, +GridTools::distort_random, GridGenerator::extrude_triangulation, +GridTools::find_active_cell_around_point +
+(Alexander Grayver, 2019/12/06)