void
generate_cheese(boost::python::list &holes);
+ /*! @copydoc GridGenerator::plate_with_a_hole
+ */
+ void
+ generate_plate_with_a_hole(const double inner_radius = 0.4,
+ const double outer_radius = 1.,
+ const double pad_bottom = 2.,
+ const double pad_top = 2.,
+ const double pad_left = 1.,
+ const double pad_right = 1.,
+ const PointWrapper ¢er = PointWrapper(),
+ const int polar_manifold_id = 0,
+ const int tfi_manifold_id = 1,
+ const double L = 1.,
+ const unsigned int n_slices = 2,
+ const bool colorize = false);
+
/*! @copydoc GridGenerator::general_cell
*/
void
generate_hyper_cube_with_cylindrical_hole,
0,
5)
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(generate_plate_with_a_hole_overloads,
+ generate_plate_with_a_hole,
+ 0,
+ 12)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(generate_hyper_ball_overloads,
generate_hyper_ball,
1,
+ const char generate_plate_with_a_hole_docstring[] =
+ "Generate a rectangular plate with an (offset) cylindrical hole. \n"
+ "The geometry consists of 2 regions: The first is a square region \n"
+ "with length outer_radius and a hole of radius inner_radius. \n"
+ "The second region describes the remainder of the bulk material. \n";
+
+
+
const char shift_docstring[] =
"Shift every vertex of the Triangulation by the given shift vector. \n";
&TriangulationWrapper::generate_cheese,
generate_cheese_docstring,
boost::python::args("self", "holes"))
+ .def("generate_plate_with_a_hole",
+ &TriangulationWrapper::generate_plate_with_a_hole,
+ generate_plate_with_a_hole_overloads(
+ boost::python::args("self",
+ "inner_radius",
+ "outer_radius",
+ "pad_bottom",
+ "pad_top",
+ "pad_left",
+ "pad_right",
+ "center",
+ "polar_manifold_id",
+ "tfi_manifold_id",
+ "L",
+ "n_slices",
+ "colorize"),
+ generate_plate_with_a_hole_docstring))
.def("generate_general_cell",
&TriangulationWrapper::generate_general_cell,
generate_general_cell_overloads(
{
dim = other.dim;
- AssertThrow(other.point != nullptr,
- ExcMessage("Underlying point does not exist."));
-
if (dim == 2)
{
Point<2> *other_point = static_cast<Point<2> *>(other.point);
new Point<3>((*other_point)[0], (*other_point)[1], (*other_point)[2]);
}
else
- AssertThrow(false,
- ExcMessage("The dimension of the point should be 2 or 3."));
+ {
+ point = nullptr;
+ }
}
} // namespace python
+ template <int dim>
+ void
+ generate_plate_with_a_hole(const double inner_radius,
+ const double outer_radius,
+ const double pad_bottom,
+ const double pad_top,
+ const double pad_left,
+ const double pad_right,
+ const PointWrapper ¢er,
+ const int polar_manifold_id,
+ const int tfi_manifold_id,
+ const double L,
+ const unsigned int n_slices,
+ const bool colorize,
+ void *triangulation)
+ {
+ // Cast the PointWrapper object to Point<dim>
+ Point<dim> point =
+ center.get_point() ?
+ *(static_cast<const Point<dim> *>(center.get_point())) :
+ Point<dim>();
+
+ Triangulation<dim, dim> *tria =
+ static_cast<Triangulation<dim, dim> *>(triangulation);
+ tria->clear();
+ GridGenerator::plate_with_a_hole(*tria,
+ inner_radius,
+ outer_radius,
+ pad_bottom,
+ pad_top,
+ pad_left,
+ pad_right,
+ point,
+ polar_manifold_id,
+ tfi_manifold_id,
+ L,
+ n_slices,
+ colorize);
+ }
+
+
+
template <int dim>
void
generate_general_cell(std::vector<PointWrapper> &wrapped_points,
+ void
+ TriangulationWrapper::generate_plate_with_a_hole(const double inner_radius,
+ const double outer_radius,
+ const double pad_bottom,
+ const double pad_top,
+ const double pad_left,
+ const double pad_right,
+ const PointWrapper ¢er,
+ const int polar_manifold_id,
+ const int tfi_manifold_id,
+ const double L,
+ const unsigned int n_slices,
+ const bool colorize)
+ {
+ if (dim == 2)
+ internal::generate_plate_with_a_hole<2>(inner_radius,
+ outer_radius,
+ pad_bottom,
+ pad_top,
+ pad_left,
+ pad_right,
+ center,
+ polar_manifold_id,
+ tfi_manifold_id,
+ L,
+ n_slices,
+ colorize,
+ triangulation);
+ else
+ internal::generate_plate_with_a_hole<3>(inner_radius,
+ outer_radius,
+ pad_bottom,
+ pad_top,
+ pad_left,
+ pad_right,
+ center,
+ polar_manifold_id,
+ tfi_manifold_id,
+ L,
+ n_slices,
+ colorize,
+ triangulation);
+ }
+
+
+
void
TriangulationWrapper::generate_general_cell(boost::python::list &vertices,
const bool colorize)
n_cells = triangulation.n_active_cells()
self.assertEqual(n_cells, 175)
+ def test_generate_plate_with_hole(self):
+ for dim in self.dim:
+ triangulation = Triangulation(dim[0])
+ triangulation.generate_plate_with_a_hole()
+ n_cells = triangulation.n_active_cells()
+ self.assertEqual(n_cells, 28)
+
def test_generate_general_cell(self):
for dim in self.restricted_dim:
triangulation = Triangulation(dim[0], dim[1])