CGAL 6.0 - Tetrahedral Remeshing
Loading...
Searching...
No Matches
User Manual

Authors
Jane Tournois, Noura Faraj, Jean-Marc Thiery, Tamy Boubekeur


Multi-Material Isotropic Tetrahedral Remeshing

This package implements an algorithm for quality tetrahedral remeshing, introduced by Faraj et al in [1]. This practical iterative remeshing algorithm is designed to remesh multi-material tetrahedral meshes, by iteratively performing a sequence of elementary operations such as edge splits, edge collapses, edge flips, and vertex relocations following a Laplacian smoothing. The algorithm results in high-quality uniform isotropic meshes, with the desired mesh density, while preserving the input geometric curve and surface features.

Specific remeshing rules have been designed to satisfy the following criteria. First, the algorithm preserves the geometric complex topology, including multi-material surface patches and polyline features. Polyline features can be defined as intersections between more than two subdomains, or listed by the user. Second, it has been made possible to remesh only a selection of cells, instead of remeshing the whole domain, while preserving or remeshing the interface surfaces between the preserved and the remeshed tetrahedra.

All the local atomic operations that are performed by the algorithm preserve the input topology of the geometric complex.

The tetrahedral remeshing algorithm improves the quality of dihedral angles, while targeting the user-defined uniform sizing field and preserving the topology of the feature complex, as highlighted by Figure Figure 61.1.

Experimental evidence shows that a higher number of remeshing iterations leads to a mesh with a improved fidelity to the sizing criterion, and higher quality dihedral angles.

Figure 61.1 Tetrahedral mesh, modified by our uniform tetrahedral remeshing method. (Left) Before remeshing, dihedral angles were in the interval [1.3; 177.8]. (Right) After remeshing and keeping the same density, dihedral angles are in the interval [12,7; 157.7].


API

The tetrahedral remeshing algorithm is implemented as a single free function CGAL::tetrahedral_isotropic_remeshing() that takes only two required parameters: the input triangulation, and the desired edge length, which drives the remeshing process.

Named Parameters are used to deal with optional parameters. The page Named Parameters describes their usage.

Examples

Tetrahedral Remeshing Example

The following example shows the simplest use of the tetrahedral remeshing function. The only required parameter is a given target edge length that drives the remeshing process towards a high-quality tetrahedral mesh with improved dihedral angles, and a more uniform mesh, with edge lengths getting closer to the input parameter value.

By default, the cells with a non-zero Subdomain_index are selected for remeshing. The cells with Subdomain_index equal to zero are ignored by the remeshing algorithm.


File Tetrahedral_remeshing/tetrahedral_remeshing_example.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Tetrahedral_remeshing/Remeshing_triangulation_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include "tetrahedral_remeshing_generate_input.h"
#include <iostream>
#include <fstream>
#include <string>
int main(int argc, char* argv[])
{
const double target_edge_length = (argc > 1) ? atof(argv[1]) : 0.1;
const std::size_t nbv = (argc > 2) ? atoi(argv[2]) : 1000;
Remeshing_triangulation tr;
CGAL::Tetrahedral_remeshing::insert_random_points_in_cube(nbv, tr);
for (auto cell : tr.finite_cell_handles())
cell->set_subdomain_index(1);
CGAL::tetrahedral_isotropic_remeshing(tr, target_edge_length);
return EXIT_SUCCESS;
}
The class Remeshing_triangulation_3 is a class template which provides a valid triangulation type tha...
Definition: Remeshing_triangulation_3.h:59
void tetrahedral_isotropic_remeshing(CGAL::Triangulation_3< Traits, TDS, SLDS > &tr, const double &target_edge_length, const NamedParameters &np=parameters::default_values())
remeshes a tetrahedral mesh.
Definition: tetrahedral_remeshing.h:174

Tetrahedral Remeshing of A Selection

Optional BGL named parameters offer more precise control on the remeshing process. In this example, a triangulation with two subdomains (defined by indices stored in cells) is given as input, but only one (defined by the Subdomain_index 2) of its subdomains is remeshed.

Only the cells with a non-zero Subdomain_index will be remeshed. The named parameter cell_is_selected_map can be used to change this behavior.


File Tetrahedral_remeshing/tetrahedral_remeshing_of_one_subdomain.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Tetrahedral_remeshing/Remeshing_triangulation_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include "tetrahedral_remeshing_generate_input.h"
using Subdomain_index = Remeshing_triangulation::Tds::Cell::Subdomain_index;
template<typename Tr>
struct Cells_of_subdomain_pmap
{
const int m_subdomain;
public:
using key_type = typename Tr::Cell_handle;
using value_type = bool;
using reference = bool;
using category = boost::read_write_property_map_tag;
friend value_type get(const Cells_of_subdomain_pmap& map,
const key_type& c)
{
return (map.m_subdomain == c->subdomain_index());
}
friend void put(Cells_of_subdomain_pmap&,
const key_type&,
const value_type)
{
; //nothing to do : subdomain indices are updated in remeshing
}
};
Subdomain_index subdomain_on_side_of_plane(const K::Point_3& p,
const K::Plane_3& plane)
{
if (plane.has_on_positive_side(p))
return 1;
else
return 2;
}
K::Point_3 centroid(const Remeshing_triangulation::Cell_handle c)
{
return CGAL::centroid(c->vertex(0)->point(),
c->vertex(1)->point(),
c->vertex(2)->point(),
c->vertex(3)->point());
}
int main(int argc, char* argv[])
{
const double target_edge_length = (argc > 1) ? atof(argv[1]) : 0.05;
const std::size_t nbv = (argc > 2) ? atoi(argv[2]) : 1000;
const std::size_t nbv_on_plane = (argc > 3) ? atoi(argv[3]) : 100;
Remeshing_triangulation tr;
CGAL::Tetrahedral_remeshing::insert_random_points_in_cube(nbv, tr);
const K::Plane_3 plane(0, 0, 1, 0);
CGAL::Tetrahedral_remeshing::insert_points_on_plane(plane, nbv_on_plane, tr);
for (auto cell : tr.finite_cell_handles())
{
const K::Point_3 cc = centroid(cell);
cell->set_subdomain_index(subdomain_on_side_of_plane(cc, plane));
}
target_edge_length,
CGAL::parameters::cell_is_selected_map(
Cells_of_subdomain_pmap<Remeshing_triangulation>{2}));
std::ofstream os("out.mesh");
return EXIT_SUCCESS;
}
void write_MEDIT(std::ostream &os, const T3 &t3, const NamedParameters &np=parameters::default_values())
CGAL::Point_2< Kernel > centroid(const CGAL::Point_2< Kernel > &p, const CGAL::Point_2< Kernel > &q, const CGAL::Point_2< Kernel > &r)

Tetrahedral Remeshing With Polyline Features

Optional BGL named parameters offer more precise control on the remeshing process. In this example, a triangulation with polyline features that should be preserved - though resampled - during the remeshing process, is given as input. Preserving all surfaces exactly could also be achieved by setting the named parameter remesh_boundaries to false.


File Tetrahedral_remeshing/tetrahedral_remeshing_with_features.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Tetrahedral_remeshing/Remeshing_triangulation_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include <CGAL/property_map.h>
#include <unordered_set>
#include <iostream>
#include <utility>
#include <cassert>
#include "tetrahedral_remeshing_generate_input.h"
using Point = Remeshing_triangulation::Point;
using Vertex_handle = Remeshing_triangulation::Vertex_handle;
using Vertex_pair = std::pair<Vertex_handle, Vertex_handle>;
using Constraints_set = std::unordered_set<Vertex_pair, boost::hash<Vertex_pair>>;
using Constraints_pmap = CGAL::Boolean_property_map<Constraints_set>;
int main(int argc, char* argv[])
{
const double target_edge_length = (argc > 1) ? atof(argv[1]) : 0.02;
const int nb_iter = (argc > 2) ? atoi(argv[2]) : 1;
const int nbv = (argc > 3) ? atoi(argv[3]) : 500;
Remeshing_triangulation t3;
Constraints_set constraints;
CGAL::Tetrahedral_remeshing::generate_input_cube(nbv, t3, constraints);
make_constraints_from_cube_edges(t3, constraints);
CGAL::tetrahedral_isotropic_remeshing(t3, target_edge_length,
CGAL::parameters::edge_is_constrained_map(Constraints_pmap(constraints))
.number_of_iterations(nb_iter));
return EXIT_SUCCESS;
}

It is also possible to define the polyline features as the ones stored as complex edges in a Mesh_complex_3_in_triangulation_3 (e.g. generated by the 3D Mesh Generation package mesh generation algorithms).


File Tetrahedral_remeshing/mesh_and_remesh_c3t3.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Mesh_triangulation_3.h>
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
#include <CGAL/Mesh_criteria_3.h>
#include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include <CGAL/IO/File_medit.h>
// Domain
#ifdef CGAL_CONCURRENT_MESH_3
using Concurrency_tag = CGAL::Parallel_tag;
#else
using Concurrency_tag = CGAL::Sequential_tag;
#endif
// Triangulation for Meshing
Tr, Mesh_domain::Corner_index, Mesh_domain::Curve_index>;
// Criteria
using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
// Triangulation for Remeshing
using Triangulation_3 = CGAL::Triangulation_3<Tr::Geom_traits,
Tr::Triangulation_data_structure>;
using Vertex_handle = Triangulation_3::Vertex_handle;
using Vertex_pair = std::pair<Vertex_handle, Vertex_handle>;
using Constraints_set = std::unordered_set<Vertex_pair, boost::hash<Vertex_pair>>;
using Constraints_pmap = CGAL::Boolean_property_map<Constraints_set>;
using Corners_set = std::unordered_set<Vertex_handle, boost::hash<Vertex_handle>>;
using Corners_pmap = CGAL::Boolean_property_map<Corners_set>;
// To avoid verbose function and named parameters call
using namespace CGAL::parameters;
int main(int argc, char* argv[])
{
const std::string fname = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/fandisk.off");
std::ifstream input(fname);
Polyhedron polyhedron;
input >> polyhedron;
if (input.fail() || !CGAL::is_triangle_mesh(polyhedron)) {
std::cerr << "Error: Input invalid " << fname << std::endl;
return EXIT_FAILURE;
}
// Create domain
Mesh_domain domain(polyhedron);
// Get sharp features
domain.detect_features();
// Mesh criteria
Mesh_criteria criteria(edge_size = 0.025,
facet_angle = 25, facet_size = 0.05, facet_distance = 0.005,
cell_radius_edge_ratio = 3, cell_size = 0.05);
// Mesh generation
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria);
// Property map of constraints
Constraints_set constraints;
Constraints_pmap constraints_pmap(constraints);
Corners_set corners;
Corners_pmap corners_pmap(corners);
CGAL::parameters::edge_is_constrained_map(constraints_pmap).
vertex_is_constrained_map(corners_pmap));
//note we use the move semantic, with std::move(c3t3),
// to avoid a copy of the triangulation by the function
// `CGAL::convert_to_triangulation_3()`
// After the call to this function, c3t3 is an empty and valid C3t3.
//It is possible to use : CGAL::convert_to_triangulation_3(c3t3),
// Then the triangulation is copied and duplicated, and c3t3 remains as is.
const double target_edge_length = 0.1;//coarsen the mesh
CGAL::tetrahedral_isotropic_remeshing(tr, target_edge_length,
CGAL::parameters::number_of_iterations(3)
.edge_is_constrained_map(constraints_pmap));
std::ofstream out("out_remeshed.mesh");
out.close();
return EXIT_SUCCESS;
}
Triangulation_data_structure::Vertex_handle Vertex_handle
bool is_triangle_mesh(const FaceGraph &g)
CGAL::Triangulation_3< typename Tr::Geom_traits, typename Tr::Triangulation_data_structure > convert_to_triangulation_3(CGAL::Mesh_complex_3_in_triangulation_3< Tr, CornerIndex, CurveIndex > c3t3, const NamedParameters &np=parameters::default_values())
converts the triangulation contained in the input to a Triangulation_3.
Definition: tetrahedral_remeshing.h:371

Tetrahedral Remeshing After Mesh Generation

The tetrahedral remeshing algorithm is designed as a post-processing for mesh generation algorithms. The API allows to generate a tetrahedral mesh with the CGAL 3D Mesh Generation package, and further improve it with the tetrahedral remeshing algorithm. This example shows how to use tetrahedral mesh generation and remeshing in sequence, from a polyhedral domain with features.


File Tetrahedral_remeshing/mesh_and_remesh_polyhedral_domain_with_features.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Mesh_triangulation_3.h>
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
#include <CGAL/Mesh_criteria_3.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/tetrahedral_remeshing.h>
// Domain
using Polyhedron = CGAL::Surface_mesh<K::Point_3>;
#ifdef CGAL_CONCURRENT_MESH_3
using Concurrency_tag = CGAL::Parallel_tag;
#else
using Concurrency_tag = CGAL::Sequential_tag;
#endif
// Triangulation
Tr, Mesh_domain::Corner_index, Mesh_domain::Curve_index>;
// Criteria
using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
// Triangulation for Remeshing
using Triangulation_3 = CGAL::Triangulation_3<Tr::Geom_traits,
Tr::Triangulation_data_structure>;
using Vertex_handle = Triangulation_3::Vertex_handle;
using Vertex_pair = std::pair<Vertex_handle, Vertex_handle>;
using Constraints_set = std::unordered_set<Vertex_pair, boost::hash<Vertex_pair>>;
using Constraints_pmap = CGAL::Boolean_property_map<Constraints_set>;
// To avoid verbose function and named parameters call
using namespace CGAL::parameters;
int main(int argc, char* argv[])
{
const std::string fname = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/anchor.off");
const int nb_iter = (argc > 2) ? atoi(argv[2]) : 5;
std::ifstream input(fname);
Polyhedron polyhedron;
std::string filename(fname);
input >> polyhedron;
if (input.fail()) {
std::cerr << "Error: Cannot read file " << fname << std::endl;
return EXIT_FAILURE;
}
if (!CGAL::is_triangle_mesh(polyhedron)) {
std::cerr << "Input geometry is not triangulated." << std::endl;
return EXIT_FAILURE;
}
// Create domain
Mesh_domain domain(polyhedron);
// Get sharp features
domain.detect_features();
// Mesh criteria
const double size = 0.072;
Mesh_criteria criteria(edge_size = size,
facet_angle = 25,
facet_size = size,
cell_radius_edge_ratio = 2,
cell_size = size);
// Mesh generation
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria, no_perturb(), no_exude());
Constraints_set constraints;
Constraints_pmap constraints_pmap(constraints);
CGAL::parameters::edge_is_constrained_map(constraints_pmap));
// Remeshing
CGAL::parameters::number_of_iterations(nb_iter));
std::ofstream out("out_remeshed.mesh");
out.close();
return EXIT_SUCCESS;
}
unspecified_type no_perturb()
unspecified_type no_exude()

Tetrahedral Remeshing from Any Tetrahedral Mesh

The following example shows how to read a mesh from a triangulation stored in a Medit file, perform tetrahedral remeshing, and save the output triangulation. The input triangulation should follow the validity requirements of a CGAL::Triangulation_3 (valid connectivity, positive orientation of the cells, and coverage of the convex hull of the vertices).


File Tetrahedral_remeshing/tetrahedral_remeshing_from_mesh.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Tetrahedral_remeshing/Remeshing_triangulation_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include <CGAL/IO/File_medit.h>
#include <fstream>
int main(int argc, char* argv[])
{
const char* filename = (argc > 1) ? argv[1] : "data/sphere.mesh";
const double target_edge_length = (argc > 2) ? atof(argv[2]) : 0.1;
Remeshing_triangulation tr;
std::ifstream is(filename, std::ios_base::in);
CGAL::tetrahedral_isotropic_remeshing(tr, target_edge_length);
std::ofstream os("after_remeshing.mesh");
return EXIT_SUCCESS;
}
bool read_MEDIT(std::istream &in, T3 &t3, const NamedParameters &np=parameters::default_values())

Implementation History

This package implements the uniform version of the "Multi-Material Adaptive Volume Remesher" algorithm for quality tetrahedral remeshing, described by Noura Faraj et al. in [1].

A first version of the code was written by Noura Faraj, Jean-Marc Thiery, and Tamy Boubekeur. Jane Tournois worked on the finalization of the code, the API, and documentation.

It was initially published in CGAL-5.1.