const unsigned int n_slices = 2,
const bool colorize = false);
+ /*! @copydoc GridGenerator::generate_channel_with_cylinder
+ */
+ void
+ generate_channel_with_cylinder(const double shell_region_width = 0.03,
+ const unsigned int n_shells = 2,
+ const double skewness = 2.,
+ const bool colorize = false);
+
/*! @copydoc GridGenerator::general_cell
*/
void
generate_plate_with_a_hole,
0,
12)
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(
+ generate_channel_with_cylinder_overloads,
+ generate_channel_with_cylinder,
+ 0,
+ 4)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(generate_hyper_ball_overloads,
generate_hyper_ball,
1,
+ const char generate_channel_with_cylinder_docstring[] =
+ "Generate a grid consisting of a channel with a cylinder. \n"
+ "The channel has three distinct regions: \n"
+ " 1. If n_shells is greater than zero, then there are that many \n"
+ " shells centered around the cylinder, \n"
+ " 2. a blending region between the shells and the rest of the \n"
+ " triangulation, and \n"
+ " 3. a bulk region consisting of Cartesian cells. \n";
+
+
+
const char shift_docstring[] =
"Shift every vertex of the Triangulation by the given shift vector. \n";
"n_slices",
"colorize"),
generate_plate_with_a_hole_docstring))
+ .def("generate_channel_with_cylinder",
+ &TriangulationWrapper::generate_channel_with_cylinder,
+ generate_channel_with_cylinder_overloads(
+ boost::python::args(
+ "shell_region_width", "n_shells", "skewness", "colorize"),
+ generate_channel_with_cylinder_docstring))
.def("generate_general_cell",
&TriangulationWrapper::generate_general_cell,
generate_general_cell_overloads(
+ template <int dim>
+ void
+ generate_channel_with_cylinder(const double shell_region_width,
+ const unsigned int n_shells,
+ const double skewness,
+ const bool colorize,
+ void *triangulation)
+ {
+ Triangulation<dim, dim> *tria =
+ static_cast<Triangulation<dim, dim> *>(triangulation);
+ tria->clear();
+ GridGenerator::channel_with_cylinder(
+ *tria, shell_region_width, n_shells, skewness, colorize);
+ }
+
+
+
template <int dim>
void
generate_general_cell(std::vector<PointWrapper> &wrapped_points,
+ void
+ TriangulationWrapper::generate_channel_with_cylinder(
+ const double shell_region_width,
+ const unsigned int n_shells,
+ const double skewness,
+ const bool colorize)
+ {
+ if (dim == 2)
+ internal::generate_channel_with_cylinder<2>(
+ shell_region_width, n_shells, skewness, colorize, triangulation);
+ else
+ internal::generate_channel_with_cylinder<3>(
+ shell_region_width, n_shells, skewness, colorize, triangulation);
+ }
+
+
+
void
TriangulationWrapper::generate_general_cell(boost::python::list &vertices,
const bool colorize)
n_cells = triangulation.n_active_cells()
self.assertEqual(n_cells, 28)
+ def test_generate_channel_with_cylinder(self):
+ for dim in self.dim:
+ triangulation = Triangulation(dim[0])
+ triangulation.generate_channel_with_cylinder()
+ n_cells = triangulation.n_active_cells()
+ if (dim[0] == '2D'):
+ self.assertEqual(n_cells, 108)
+ else:
+ self.assertEqual(n_cells, 432)
+
def test_generate_general_cell(self):
for dim in self.restricted_dim:
triangulation = Triangulation(dim[0], dim[1])