double
measure() const;
+ /**
+ * Test whether the cell has children (this is the criterion for activity of
+ * a cell).
+ */
+ bool
+ active() const;
+
+ /**
+ * Return the level within the mesh hierarchy at which this cell is
+ * located.
+ */
+ int
+ level() const;
+
+ /**
+ * Return the index of the element presently pointed to on the present
+ * level.
+ */
+ int
+ index() const;
+
+ /**
+ * Return whether the neighbor is coarser then the present cell.
+ */
+ bool
+ neighbor_is_coarser(const unsigned int neighbor) const;
+
+ /**
+ * Return the how-many'th neighbor this cell is of cell->neighbor(neighbor),
+ * i.e. return the face_no such that
+ * cell->neighbor(neighbor)->neighbor(face_no)==cell.
+ */
+ unsigned int
+ neighbor_of_neighbor(const unsigned int neighbor) const;
+
+ /**
+ * Return the global index of i-th vertex of a cell.
+ */
+ unsigned int
+ vertex_index(const unsigned int i) const;
+
/**
* Exception.
*/
PointWrapper
transform_real_to_unit_cell(CellAccessorWrapper &cell, PointWrapper &point);
+ /**
+ * Transform the point on the real cell to the corresponding point on the
+ * unit cell, and then projects it to a dim-1 point on the face with the
+ * given face number 'face_no'. Ideally the point 'point' is near the
+ * face 'face_no', but any point in the cell can technically be projected.
+ * The returned point is of dimension dim with dim-1 coodinate value
+ * explicitly set to zero.
+ */
+ PointWrapper
+ project_real_point_to_unit_point_on_face(CellAccessorWrapper &cell,
+ const unsigned int face_no,
+ PointWrapper & point);
+
/**
* Get the underlying mapping.
*/
double
norm_square() const;
+ /**
+ * Convert point's coordinates to a python list with [x,y] or [x,y,z]
+ * for 2-D and 3-D, respectively.
+ */
+ boost::python::list
+ to_list() const;
+
/**
* Assignment operator. The dimension of the point is changed if it is
* different than the one of @p other.
PointWrapper & p,
MappingQGenericWrapper mapping = MappingQGenericWrapper());
+ /**
+ * Find and return a list of active cells that surround a given vertex with
+ * index vertex_index.
+ */
+ boost::python::list
+ find_cells_adjacent_to_vertex(const unsigned int vertex_index);
+
/**
* Assign a manifold object to a certain part of the triangulation.
* If an object with manifold number is refined, this object
boost::python::list
active_cells();
+ /**
+ * Return the diameter of the smallest active cell of a triangulation
+ */
+ double
+ minimal_cell_diameter() const;
+
+ /**
+ * Return the diameter of the largest active cell of a triangulation
+ */
+ double
+ maximal_cell_diameter() const;
+
/**
* Write mesh to the output file @filename according to the given data
* format.
- template <int dim, int spacedim>
- bool
- at_boundary(const void *cell_accessor)
- {
- const CellAccessor<dim, spacedim> *cell =
- static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor);
- return cell->at_boundary();
- }
-
-
-
- template <int dim, int spacedim>
- bool
- has_boundary_lines(const void *cell_accessor)
- {
- const CellAccessor<dim, spacedim> *cell =
- static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor);
- return cell->has_boundary_lines();
- }
-
-
-
template <int dim, int spacedim>
const CellAccessor<dim, spacedim> *
neighbor(const int i, const void *cell_accessor)
template <int dim, int spacedim>
- double
- measure(const void *cell_accessor)
+ const CellAccessor<dim, spacedim> *
+ cell_cast(const void *cell_accessor)
{
- const CellAccessor<dim, spacedim> *cell =
- static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor);
-
- return cell->measure();
+ return static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor);
}
} // namespace internal
CellAccessorWrapper::at_boundary() const
{
if ((dim == 2) && (spacedim == 2))
- return internal::at_boundary<2, 2>(cell_accessor);
+ return internal::cell_cast<2, 2>(cell_accessor)->at_boundary();
else if ((dim == 2) && (spacedim == 3))
- return internal::at_boundary<2, 3>(cell_accessor);
+ return internal::cell_cast<2, 3>(cell_accessor)->at_boundary();
else
- return internal::at_boundary<3, 3>(cell_accessor);
+ return internal::cell_cast<3, 3>(cell_accessor)->at_boundary();
}
CellAccessorWrapper::has_boundary_lines() const
{
if ((dim == 2) && (spacedim == 2))
- return internal::has_boundary_lines<2, 2>(cell_accessor);
+ return internal::cell_cast<2, 2>(cell_accessor)->has_boundary_lines();
else if ((dim == 2) && (spacedim == 3))
- return internal::has_boundary_lines<2, 3>(cell_accessor);
+ return internal::cell_cast<2, 3>(cell_accessor)->has_boundary_lines();
else
- return internal::has_boundary_lines<3, 3>(cell_accessor);
+ return internal::cell_cast<3, 3>(cell_accessor)->has_boundary_lines();
}
CellAccessorWrapper::measure() const
{
if ((dim == 2) && (spacedim == 2))
- return internal::measure<2, 2>(cell_accessor);
+ return internal::cell_cast<2, 2>(cell_accessor)->measure();
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)->measure();
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)->measure();
+ }
+
+
+
+ bool
+ CellAccessorWrapper::active() const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)->active();
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)->active();
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)->active();
+ }
+
+
+
+ int
+ CellAccessorWrapper::level() const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)->level();
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)->level();
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)->level();
+ }
+
+
+
+ int
+ CellAccessorWrapper::index() const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)->index();
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)->index();
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)->index();
+ }
+
+
+
+ bool
+ CellAccessorWrapper::neighbor_is_coarser(const unsigned int neighbor) const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)
+ ->neighbor_is_coarser(neighbor);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)
+ ->neighbor_is_coarser(neighbor);
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)
+ ->neighbor_is_coarser(neighbor);
+ }
+
+
+
+ unsigned int
+ CellAccessorWrapper::neighbor_of_neighbor(const unsigned int neighbor) const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)
+ ->neighbor_of_neighbor(neighbor);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::cell_cast<2, 3>(cell_accessor)
+ ->neighbor_of_neighbor(neighbor);
+ else
+ return internal::cell_cast<3, 3>(cell_accessor)
+ ->neighbor_of_neighbor(neighbor);
+ }
+
+
+
+ unsigned int
+ CellAccessorWrapper::vertex_index(const unsigned int i) const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::cell_cast<2, 2>(cell_accessor)->vertex_index(i);
else if ((dim == 2) && (spacedim == 3))
- return internal::measure<2, 3>(cell_accessor);
+ return internal::cell_cast<2, 3>(cell_accessor)->vertex_index(i);
else
- return internal::measure<3, 3>(cell_accessor);
+ return internal::cell_cast<3, 3>(cell_accessor)->vertex_index(i);
}
} // namespace python
const char coarsen_flag_docstring[] =
- "Get/Set the coarsen_flag of the cell \n";
+ "Get/Set the coarsen_flag of the cell. \n";
const char material_id_docstring[] =
- "Get/Set the material_id of the cell \n";
+ "Get/Set the material_id of the cell. \n";
const char manifold_id_docstring[] =
- "Get/Set the manifold_id of the cell \n";
+ "Get/Set the manifold_id of the cell. \n";
const char barycenter_docstring[] =
- "Return the barycenter of the current cell \n";
+ "Return the barycenter of the current cell. \n";
const char set_vertex_docstring[] =
- " Set the ith vertex of the cell to point_wrapper \n";
+ " Set the ith vertex of the cell to point_wrapper. \n";
const char get_vertex_docstring[] =
- " Get the ith vertex of the cell \n";
+ " Get the ith vertex of the cell. \n";
const char at_boundary_docstring[] =
- " Return whether the cell is at the boundary \n";
+ " Return whether the cell is at the boundary. \n";
+ const char vertex_index_docstring[] =
+ " Return the global index of i-th vertex of a cell. \n";
+
+
+
+ const char neighbor_of_neighbor_docstring[] =
+ " Return the how-many'th neighbor this cell is of \n"
+ " cell.neighbor(neighbor), i.e. return the face_no such that \n"
+ " cell.neighbor(neighbor).neighbor(face_no)==cell. \n";
+
+
+
+ const char neighbor_is_coarser_docstring[] =
+ " Return whether the neighbor is coarser then the present cell. \n";
+
+
+
+ const char index_docstring[] =
+ " Return the index of the element presently pointed to \n"
+ " on the present. \n";
+
+
+
+ const char level_docstring[] =
+ " Return the level within the mesh hierarchy at which this cell \n"
+ " is located. \n";
+
+
+
+ const char active_docstring[] =
+ " Test whether the cell has children (this is the criterion \n"
+ " for activity of a cell). \n";
+
+
void
export_cell_accessor()
{
.def("measure",
&CellAccessorWrapper::measure,
measure_docstring,
- boost::python::args("self"));
+ boost::python::args("self"))
+ .def("active",
+ &CellAccessorWrapper::active,
+ active_docstring,
+ boost::python::args("self"))
+ .def("level",
+ &CellAccessorWrapper::level,
+ level_docstring,
+ boost::python::args("self"))
+ .def("index",
+ &CellAccessorWrapper::index,
+ index_docstring,
+ boost::python::args("self"))
+ .def("neighbor_is_coarser",
+ &CellAccessorWrapper::neighbor_is_coarser,
+ neighbor_is_coarser_docstring,
+ boost::python::args("self", "neighbor"))
+ .def("neighbor_of_neighbor",
+ &CellAccessorWrapper::neighbor_of_neighbor,
+ neighbor_of_neighbor_docstring,
+ boost::python::args("self", "neighbor"))
+ .def("vertex_index",
+ &CellAccessorWrapper::vertex_index,
+ vertex_index_docstring,
+ boost::python::args("self", "vertex"));
}
} // namespace python
" Map the point p on the unit cell to the corresponding point \n"
" on the real cell. \n";
+
const char transform_real_to_unit_cell_docstring[] =
" Map the point p on the real cell to the corresponding point \n"
" on the unit cell. \n";
+ const char project_real_point_to_unit_point_on_face_docstring[] =
+ " Transform the point on the real cell to the corresponding point \n"
+ " on the unit cell, and then projects it to a dim-1 point on the \n"
+ " face with the given face number face_no. Ideally the point is \n"
+ " near the face face_no, but any point in the cell can technically \n"
+ " be projected. The returned point is of dimension dim with \n"
+ " dim-1 coodinate value explicitly set to zero. \n";
+
void
export_mapping()
.def("transform_real_to_unit_cell",
&MappingQGenericWrapper::transform_real_to_unit_cell,
transform_real_to_unit_cell_docstring,
- boost::python::args("self", "cell_wrapper", "point_wrapper"))
+ boost::python::args("self", "cell", "point"))
.def("transform_unit_to_real_cell",
&MappingQGenericWrapper::transform_unit_to_real_cell,
transform_unit_to_real_cell_docstring,
- boost::python::args("self", "cell_wrapper", "point_wrapper"));
+ boost::python::args("self", "cell", "point"))
+ .def("project_real_point_to_unit_point_on_face",
+ &MappingQGenericWrapper::project_real_point_to_unit_point_on_face,
+ project_real_point_to_unit_point_on_face_docstring,
+ boost::python::args("self", "cell", "face_no", "point"));
}
} // namespace python
namespace python
{
+ const char to_list_docstring[] =
+ "Convert point's coordinates to a python list with [x,y] or [x,y,z] \n"
+ "for 2-D or 3-D, respectively. \n";
+
+
+
const char distance_docstring[] =
"Return the Euclidean distance of this point to the point p \n";
{
boost::python::class_<PointWrapper>(
"Point", boost::python::init<boost::python::list>())
+ .def("to_list",
+ &PointWrapper::to_list,
+ to_list_docstring,
+ boost::python::args("self"))
.def("distance",
&PointWrapper::distance,
distance_docstring,
2)
const char n_active_cells_docstring[] =
- "Return the number of active cells \n";
+ "Return the number of active cells. \n";
+
+
+
+ const char dim_docstring[] =
+ "Return the dimension of the Triangulation \n";
+
+
+
+ const char spacedim_docstring[] =
+ "Return the space dimension of the Triangulation. \n";
const char generate_hyper_cube_docstring[] =
"Generate a hyper_cube (square in 2D and cube in 3D) \n"
- "with exactly one cell \n";
+ "with exactly one cell. \n";
const char generate_simplex_docstring[] =
- "Generate a simplex with (dim+1) vertices and mesh cells \n";
+ "Generate a simplex with (dim+1) vertices and mesh cells. \n";
const char generate_hyper_rectangle_docstring[] =
"Generate a coordinate-parallel brick from the two diagonally opposite \n"
- "corners points p1 and p2 \n";
+ "corners points p1 and p2. \n";
const char generate_subdivided_hyper_rectangle_docstring[] =
"Generate a coordinate-parallel brick from the two diagonally opposite \n"
"corners point p1 and p2. In direction i, repetitions[i] cells are \n"
- "created \n";
+ "created. \n";
const char shift_docstring[] =
- "Shift every vertex of the Triangulation by the given shift vector \n";
+ "Shift every vertex of the Triangulation by the given shift vector. \n";
const char merge_docstring[] =
"Given two triangulations, create the triangulation that contains \n"
- "the cells of both triangulations \n";
+ "the cells of both triangulations. \n";
const char refine_global_docstring[] =
- "Refine all the cells times time \n";
+ "Refine all the cells times time. \n";
const char execute_coarsening_and_refinement_docstring[] =
- "Execute both refinement and coarsening of the Triangulation \n";
+ "Execute both refinement and coarsening of the Triangulation. \n";
const char active_cells_docstring[] =
- "Return the list of active cell accessors of the Triangulation \n";
+ "Return the list of active cell accessors of the Triangulation. \n";
const char save_docstring[] =
- "Write the Triangulation to a file \n";
+ "Write the Triangulation to a file. \n";
const char load_docstring[] =
- "Load the Triangulation from a file \n";
+ "Load the Triangulation from a file. \n";
+ const char find_cells_adjacent_to_vertex_docstring[] =
+ "Find and return a list of active cells that surround a given \n"
+ "vertex with index vertex_index. \n";
+
+
+
+ const char minimal_cell_diameter_docstring[] =
+ "Return the diameter of the smallest active cell of a triangulation. \n";
+
+
+
+ const char maximal_cell_diameter_docstring[] =
+ "Return the diameter of the largest active cell of a triangulation. \n";
+
+
void
export_triangulation()
{
&TriangulationWrapper::n_active_cells,
n_active_cells_docstring,
boost::python::args("self"))
+ .def("dim",
+ &TriangulationWrapper::get_dim,
+ dim_docstring,
+ boost::python::args("self"))
+ .def("spacedim",
+ &TriangulationWrapper::get_spacedim,
+ spacedim_docstring,
+ boost::python::args("self"))
+ .def("n_active_cells",
+ &TriangulationWrapper::n_active_cells,
+ n_active_cells_docstring,
+ boost::python::args("self"))
+ .def("minimal_cell_diameter",
+ &TriangulationWrapper::minimal_cell_diameter,
+ minimal_cell_diameter_docstring,
+ boost::python::args("self"))
+ .def("maximal_cell_diameter",
+ &TriangulationWrapper::maximal_cell_diameter,
+ maximal_cell_diameter_docstring,
+ boost::python::args("self"))
.def("generate_hyper_cube",
&TriangulationWrapper::generate_hyper_cube,
generate_hyper_cube_overloads(
find_active_cell_around_point_overloads(
boost::python::args("self", "point", "mapping"),
find_active_cell_around_point_docstring))
+ .def("find_cells_adjacent_to_vertex",
+ &TriangulationWrapper::find_cells_adjacent_to_vertex,
+ find_cells_adjacent_to_vertex_docstring,
+ boost::python::args("self", "vertex_index"))
.def("refine_global",
&TriangulationWrapper::refine_global,
refine_global_docstring,
return PointWrapper(coord_list);
}
+
+
+
+ template <int dim, int spacedim>
+ PointWrapper
+ project_real_point_to_unit_point_on_face(void *mapping_ptr,
+ void *cell_accessor_ptr,
+ const unsigned int face_no,
+ void * point_ptr)
+ {
+ const MappingQGeneric<dim, spacedim> *mapping =
+ static_cast<const MappingQGeneric<dim, spacedim> *>(mapping_ptr);
+
+ const CellAccessor<dim, spacedim> *cell_accessor =
+ static_cast<const CellAccessor<dim, spacedim> *>(cell_accessor_ptr);
+
+ const Point<spacedim> *point =
+ static_cast<const Point<spacedim> *>(point_ptr);
+
+ typename Triangulation<dim, spacedim>::active_cell_iterator cell(
+ &cell_accessor->get_triangulation(),
+ cell_accessor->level(),
+ cell_accessor->index());
+
+ Point<dim - 1> p_face =
+ mapping->project_real_point_to_unit_point_on_face(cell,
+ face_no,
+ *point);
+
+ boost::python::list coord_list;
+ for (int i = 0; i < dim - 1; ++i)
+ coord_list.append(p_face[i]);
+ coord_list.append(0.);
+
+ return PointWrapper(coord_list);
+ }
+
} // namespace internal
+ PointWrapper
+ MappingQGenericWrapper::project_real_point_to_unit_point_on_face(
+ CellAccessorWrapper &cell,
+ const unsigned int face_no,
+ PointWrapper & p)
+ {
+ AssertThrow(
+ spacedim == p.get_dim(),
+ ExcMessage(
+ "Dimension of the point is not equal to the space dimension of the mapping."));
+
+ if ((dim == 2) && (spacedim == 2))
+ return internal::project_real_point_to_unit_point_on_face<2, 2>(
+ mapping_ptr, cell.cell_accessor, face_no, p.point);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::project_real_point_to_unit_point_on_face<2, 3>(
+ mapping_ptr, cell.cell_accessor, face_no, p.point);
+ else
+ return internal::project_real_point_to_unit_point_on_face<3, 3>(
+ mapping_ptr, cell.cell_accessor, face_no, p.point);
+ }
+
+
+
void *
MappingQGenericWrapper::get_mapping()
{
{
namespace internal
{
+ template <int dim>
+ boost::python::list
+ to_list(const void *point)
+ {
+ const Point<dim> &p = *static_cast<const Point<dim> *>(point);
+
+ boost::python::list p_list;
+ for (int d = 0; d < dim; ++d)
+ p_list.append(p[d]);
+
+ return p_list;
+ }
+
template <int dim>
double
distance(const Point<dim> &p1, const Point<dim> &p2)
+ boost::python::list
+ PointWrapper::to_list() const
+ {
+ if (dim == 2)
+ return internal::to_list<2>(point);
+ else
+ return internal::to_list<3>(point);
+ }
+
+
+
double
PointWrapper::distance(const PointWrapper &p) const
{
}
+ template <int dim, int spacedim>
+ boost::python::list
+ find_cells_adjacent_to_vertex(const unsigned int vertex_index,
+ TriangulationWrapper &triangulation_wrapper)
+ {
+ Triangulation<dim, spacedim> *tria =
+ static_cast<Triangulation<dim, spacedim> *>(
+ triangulation_wrapper.get_triangulation());
+
+ std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
+ adjacent_cells =
+ GridTools::find_cells_adjacent_to_vertex(*tria, vertex_index);
+
+ boost::python::list cells;
+ for (auto &cell : adjacent_cells)
+ cells.append(CellAccessorWrapper(triangulation_wrapper,
+ cell->level(),
+ cell->index()));
+
+ return cells;
+ }
+
+
template <int dim, int spacedim>
boost::python::list
+ template <int dim, int spacedim>
+ double
+ maximal_cell_diameter(const void *triangulation)
+ {
+ const Triangulation<dim, spacedim> *tria =
+ static_cast<const Triangulation<dim, spacedim> *>(triangulation);
+ return GridTools::maximal_cell_diameter(*tria);
+ }
+
+
+
+ template <int dim, int spacedim>
+ double
+ minimal_cell_diameter(const void *triangulation)
+ {
+ const Triangulation<dim, spacedim> *tria =
+ static_cast<const Triangulation<dim, spacedim> *>(triangulation);
+ return GridTools::minimal_cell_diameter(*tria);
+ }
+
+
+
template <int dim, int spacedim>
void
write(const std::string &filename,
+ boost::python::list
+ TriangulationWrapper::find_cells_adjacent_to_vertex(
+ const unsigned int vertex_index)
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::find_cells_adjacent_to_vertex<2, 2>(vertex_index, *this);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::find_cells_adjacent_to_vertex<2, 3>(vertex_index, *this);
+ else
+ return internal::find_cells_adjacent_to_vertex<3, 3>(vertex_index, *this);
+ }
+
+
+
void
TriangulationWrapper::refine_global(const unsigned int n)
{
+ double
+ TriangulationWrapper::minimal_cell_diameter() const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::minimal_cell_diameter<2, 2>(triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::minimal_cell_diameter<2, 3>(triangulation);
+ else
+ return internal::minimal_cell_diameter<3, 3>(triangulation);
+ }
+
+
+
+ double
+ TriangulationWrapper::maximal_cell_diameter() const
+ {
+ if ((dim == 2) && (spacedim == 2))
+ return internal::maximal_cell_diameter<2, 2>(triangulation);
+ else if ((dim == 2) && (spacedim == 3))
+ return internal::maximal_cell_diameter<2, 3>(triangulation);
+ else
+ return internal::maximal_cell_diameter<3, 3>(triangulation);
+ }
+
+
+
void
TriangulationWrapper::write(const std::string &filename,
const std::string format) const
void
TriangulationWrapper::save(const std::string &filename) const
{
- std::ofstream ofs(filename);
- boost::archive::text_oarchive oa(ofs);
+ std::ofstream ofs(filename);
+ boost::archive::binary_oarchive oa(ofs);
if ((dim == 2) && (spacedim == 2))
{
Triangulation<2, 2> *tria =
static_cast<Triangulation<2, 2> *>(triangulation);
- oa << *tria;
+ tria->save(oa, 0);
}
else if ((dim == 2) && (spacedim == 3))
{
Triangulation<2, 3> *tria =
static_cast<Triangulation<2, 3> *>(triangulation);
- oa << *tria;
+ tria->save(oa, 0);
}
}
else
Triangulation<3, 3> *tria =
static_cast<Triangulation<3, 3> *>(triangulation);
- oa << *tria;
+ tria->save(oa, 0);
}
}
void
TriangulationWrapper::load(const std::string &filename)
{
- std::ifstream ifs(filename);
- boost::archive::text_iarchive ia(ifs);
+ std::ifstream ifs(filename);
+ boost::archive::binary_iarchive ia(ifs);
if ((dim == 2) && (spacedim == 2))
{
Triangulation<2, 2> *tria =
static_cast<Triangulation<2, 2> *>(triangulation);
- ia >> *tria;
+ tria->load(ia, 0);
}
else if ((dim == 2) && (spacedim == 3))
{
Triangulation<2, 3> *tria =
static_cast<Triangulation<2, 3> *>(triangulation);
- ia >> *tria;
+ tria->load(ia, 0);
}
else
{
Triangulation<3> *tria =
static_cast<Triangulation<3, 3> *>(triangulation);
- ia >> *tria;
+ tria->load(ia, 0);
}
}
--- /dev/null
+Added: Python wrappers for the CellAccessor's methods: active, level,
+index, vertex_index, neighbor_is_coarser, neighbor_of_neighbor and
+GridTools' methods: find_cells_adjacent_to_vertex, minimal_cell_diameter,
+maximal_cell_diameter, and Mapping::project_real_point_to_unit_point_on_face
+<br>
+(Alexander Grayver, 2020/01/06)