Skip to content

Commit

Permalink
Address some comments
Browse files Browse the repository at this point in the history
  • Loading branch information
fdrmrc committed Apr 29, 2022
1 parent 4864299 commit fe0a87c
Show file tree
Hide file tree
Showing 5 changed files with 26 additions and 19 deletions.
4 changes: 0 additions & 4 deletions doc/news/changes/minor/20220428MarcoFederLucaHeltai

This file was deleted.

4 changes: 2 additions & 2 deletions include/deal.II/cgal/surface_mesh.h
Original file line number Diff line number Diff line change
Expand Up @@ -44,7 +44,7 @@ namespace CGALWrappers
* More information on this class is available at
* https://doc.cgal.org/latest/Surface_mesh/index.html
*
* The function will throw an exception in dimension one. In dimension two it
* The function will throw an exception in dimension one. In dimension two, it
* generates a surface mesh of the quadrilateral cell or of the triangle cell,
* while in dimension three it will generate the surface mesh of the cell,
* i.e., a polyhedral mesh containing the faces of the input cell.
Expand All @@ -62,7 +62,7 @@ namespace CGALWrappers
*/
template <typename CGALPointType, int dim, int spacedim>
void
to_cgal_mesh(
convert_to_cgal_surface_mesh(
const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
const dealii::Mapping<dim, spacedim> & mapping,
CGAL::Surface_mesh<CGALPointType> & mesh);
Expand Down
13 changes: 9 additions & 4 deletions source/cgal/surface_mesh.cc
Original file line number Diff line number Diff line change
Expand Up @@ -35,7 +35,6 @@ namespace
const unsigned nv = face->n_vertices();
std::vector<typename CGAL_Mesh::Vertex_index> indices;


switch (nv)
{
case 2:
Expand Down Expand Up @@ -80,11 +79,16 @@ namespace CGALWrappers
{
template <typename CGALPointType, int dim, int spacedim>
void
to_cgal_mesh(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
const Mapping<dim, spacedim> & mapping,
CGAL::Surface_mesh<CGALPointType> &mesh)
convert_to_cgal_surface_mesh(
const typename Triangulation<dim, spacedim>::cell_iterator &cell,
const Mapping<dim, spacedim> & mapping,
CGAL::Surface_mesh<CGALPointType> & mesh)
{
Assert(dim > 1, ExcImpossibleInDim(dim));
Assert(
mesh.is_empty(),
ExcMessage(
"The CGAL::Surface_mesh object must be empty upon calling this function."));
using Mesh = CGAL::Surface_mesh<CGALPointType>;
using Vertex_index = typename Mesh::Vertex_index;
const auto &vertices = mapping.get_vertices(cell);
Expand All @@ -108,6 +112,7 @@ namespace CGALWrappers
mesh,
(f % 2 == 0 || cell->n_vertices() != 8));
}

// explicit instantiations
# include "surface_mesh.inst"

Expand Down
3 changes: 2 additions & 1 deletion source/cgal/surface_mesh.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,8 @@
for (dim : DIMENSIONS; spacedim : SPACE_DIMENSIONS; cgal_kernel : CGAL_KERNELS)
{
#if dim <= spacedim
template void to_cgal_mesh<typename cgal_kernel::Point_3, dim, spacedim>(
template void
convert_to_cgal_surface_mesh<typename cgal_kernel::Point_3, dim, spacedim>(
const typename Triangulation<dim, spacedim>::cell_iterator &cell,
const Mapping<dim, spacedim> & mapping,
CGAL::Surface_mesh<typename cgal_kernel::Point_3> & mesh);
Expand Down
21 changes: 13 additions & 8 deletions tests/cgal/cgal_surface_mesh_01.cc
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,7 @@
//
// ---------------------------------------------------------------------

// Convert a deal.II cell to a cgal Surface_mesh
// Convert a deal.II cell to a cgal Surface_mesh.

#include <deal.II/base/config.h>

Expand All @@ -38,21 +38,26 @@ void
test()
{
deallog << "dim= " << dim << ",\t spacedim= " << spacedim << std::endl;
std::vector<std::vector<unsigned int>> d2t = {{}, {2}, {3, 4}, {4, 5, 6, 8}};
for (const auto nv : d2t[dim])
using namespace ReferenceCells;
std::vector<std::vector<ReferenceCell>> ref_cells = {
{},
{Line},
{Triangle, Quadrilateral},
{Tetrahedron, Pyramid, Wedge, Hexahedron}};
for (const auto r_cell : ref_cells[dim])
{
Triangulation<dim, spacedim> tria;
CGAL::Surface_mesh<CGALPoint> mesh;

const auto ref = ReferenceCell::n_vertices_to_type(dim, nv);
const auto mapping = ref.template get_default_mapping<dim, spacedim>(1);
GridGenerator::reference_cell(tria, ref);
const auto mapping =
r_cell.template get_default_mapping<dim, spacedim>(1);
GridGenerator::reference_cell(tria, r_cell);

const auto cell = tria.begin_active();
to_cgal_mesh(cell, *mapping, mesh);
convert_to_cgal_surface_mesh(cell, *mapping, mesh);

Assert(mesh.is_valid(), dealii::ExcMessage("The CGAL mesh is not valid"));
deallog << "deal vertices: " << nv << ", cgal vertices "
deallog << "deal vertices: " << cell->n_vertices() << ", cgal vertices "
<< mesh.num_vertices() << std::endl;
deallog << "deal faces: " << cell->n_faces() << ", cgal faces "
<< mesh.num_faces() << std::endl;
Expand Down

0 comments on commit fe0a87c

Please sign in to comment.