* Create CylindricalManifold
*/
void
- create_cylindrical(const int axis, const double tolerance);
+ create_cylindrical(const int axis = 0, const double tolerance = 1e-10);
/**
* Return pointer to an underlying manifold object
void
generate_half_hyper_ball(PointWrapper ¢er, const double radius = 1.);
+ /**
+ * Produce a hyper-shell, the region between two spheres around center,
+ * with given inner_radius and outer_radius. The number n_cells indicates
+ * the number of cells of the resulting triangulation, i.e., how many
+ * cells form the ring (in 2d) or the shell (in 3d).
+ * The appropriate manifold class is SphericalManifold.
+ */
+ void
+ generate_hyper_shell(PointWrapper & center,
+ const double inner_radius,
+ const double outer_radius,
+ const unsigned n_cells = 0);
+
/**
* Shift each vertex of the Triangulation by the given @p shift_list.
*/
void
write(const std::string &filename, const std::string format) const;
+ /**
+ * Read mesh from the file @filename using the given data
+ * format.
+ */
+ void
+ read(const std::string &filename, const std::string format) const;
+
/**
* Write the Triangulation in file.
*/
namespace python
{
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(create_cylindrical_overloads,
+ create_cylindrical,
+ 0,
+ 2)
+
+
const char create_spherical_docstring[] =
" Create spherical manifold with a given center point \n";
export_manifold()
{
boost::python::class_<ManifoldWrapper>(
- "ManifoldWrapper",
+ "Manifold",
boost::python::init<const int, const int>(
boost::python::args("dim", "spacedim")))
.def("create_spherical",
boost::python::args("self", "center"))
.def("create_cylindrical",
&ManifoldWrapper::create_cylindrical,
- create_cylindrical_docstring,
- boost::python::args("self", "axis", "tolerance"));
+ create_cylindrical_overloads(
+ boost::python::args("self", "axis", "tolerance"),
+ create_cylindrical_docstring));
}
} // namespace python
1,
2)
-
+ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(generate_hyper_shell_overloads,
+ generate_hyper_shell,
+ 3,
+ 4)
const char n_active_cells_docstring[] =
"Return the number of active cells \n";
+ const char generate_hyper_shell_docstring[] =
+ "Produce a hyper-shell, the region between two spheres around center, \n"
+ "with given inner_radius and outer_radius. The number n_cells indicates \n"
+ "the number of cells of the resulting triangulation, i.e., how many \n"
+ "cells form the ring (in 2d) or the shell (in 3d). \n"
+ "The appropriate manifold class is SphericalManifold. \n";
+
+
+
const char generate_half_hyper_ball_docstring[] =
"Generate a half hyper-ball around center, which contains four \n"
"elements in 2d and 6 in 3d. The cut plane is perpendicular to the \n"
+ const char read_docstring[] =
+ "Read a mesh from the file according to the given data format. \n"
+ "The possible formats are: \n"
+ " - msh \n"
+ " - vtk \n";
+
+
+
const char save_docstring[] =
"Write the Triangulation to a file \n";
generate_half_hyper_ball_overloads(
boost::python::args("self", "center", "radius"),
generate_half_hyper_ball_docstring))
+ .def("generate_hyper_shell",
+ &TriangulationWrapper::generate_hyper_shell,
+ generate_hyper_shell_overloads(boost::python::args("self",
+ "center",
+ "inner_radius",
+ "outer_radius",
+ "n_cells"),
+ generate_hyper_shell_docstring))
.def("shift",
&TriangulationWrapper::shift,
shift_docstring,
&TriangulationWrapper::write,
write_docstring,
boost::python::args("self", "filename", "format"))
+ .def("read",
+ &TriangulationWrapper::read,
+ read_docstring,
+ boost::python::args("self", "filename", "format"))
.def("save",
&TriangulationWrapper::save,
save_docstring,
{
template <int dim, int spacedim>
Manifold<dim, spacedim> *
- create_spherical_manifold(PointWrapper center)
+ create_spherical_manifold(const PointWrapper ¢er)
{
const Point<spacedim> *point =
static_cast<const Point<spacedim> *>(center.get_point());
template <int dim, int spacedim>
Manifold<dim, spacedim> *
- create_polar_manifold(PointWrapper center)
+ create_polar_manifold(const PointWrapper ¢er)
{
const Point<spacedim> *point =
static_cast<const Point<spacedim> *>(center.get_point());
{
return new PolarManifold<dim, spacedim>(*d);
}
+ else if (const CylindricalManifold<dim, spacedim> *d =
+ dynamic_cast<const CylindricalManifold<dim, spacedim> *>(
+ manifold))
+ {
+ return new CylindricalManifold<dim, spacedim>(*d);
+ }
else
ExcMessage("Unsupported manifold type in clone.");
#include <deal.II/base/types.h>
#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_in.h>
#include <deal.II/grid/grid_out.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/grid/tria.h>
}
+ template <int dim>
+ void
+ generate_hyper_shell(PointWrapper & center,
+ const double inner_radius,
+ const double outer_radius,
+ const unsigned n_cells,
+ void * triangulation)
+ {
+ // Cast the PointWrapper object to Point<dim>
+ Point<dim> center_point =
+ *(static_cast<Point<dim> *>(center.get_point()));
+
+ Triangulation<dim> *tria =
+ static_cast<Triangulation<dim> *>(triangulation);
+ tria->clear();
+ GridGenerator::hyper_shell(
+ *tria, center_point, inner_radius, outer_radius, n_cells);
+ }
+
template <int dim>
void
mesh_writer.write(*tria, ofs, output_format);
ofs.close();
}
+
+
+
+ template <int dim, int spacedim>
+ void
+ read(const std::string &filename,
+ const std::string &format,
+ void * triangulation)
+ {
+ Triangulation<dim, spacedim> *tria =
+ static_cast<Triangulation<dim, spacedim> *>(triangulation);
+
+ tria->clear();
+
+ typename GridIn<dim, spacedim>::Format input_format;
+ if (format.compare("msh") == 0)
+ input_format = GridIn<dim, spacedim>::Format::msh;
+ else if (format.compare("vtk") == 0)
+ input_format = GridIn<dim, spacedim>::Format::vtk;
+ else
+ ExcMessage("Cannot read triangulation of the given format.");
+
+ GridIn<dim, spacedim> mesh_reader;
+ mesh_reader.attach_triangulation(*tria);
+ std::ifstream ifs(filename);
+ mesh_reader.read(ifs, input_format);
+ ifs.close();
+ }
} // namespace internal
}
+
+ void
+ TriangulationWrapper::generate_hyper_shell(PointWrapper & center,
+ const double inner_radius,
+ const double outer_radius,
+ const unsigned n_cells)
+ {
+ AssertThrow(
+ dim == spacedim,
+ ExcMessage(
+ "This function is only implemented for dim equal to spacedim."));
+ if (dim == 2)
+ internal::generate_hyper_shell<2>(
+ center, inner_radius, outer_radius, n_cells, triangulation);
+ else
+ internal::generate_hyper_shell<3>(
+ center, inner_radius, outer_radius, n_cells, triangulation);
+ }
+
+
+
void
TriangulationWrapper::generate_hyper_sphere(PointWrapper ¢er,
const double radius)
+ void
+ TriangulationWrapper::read(const std::string &filename,
+ const std::string format) const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ internal::read<2, 2>(filename, format, triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ internal::read<2, 3>(filename, format, triangulation);
+ else
+ internal::read<3, 3>(filename, format, triangulation);
+ }
+
+
+
void
TriangulationWrapper::save(const std::string &filename) const
{