Tetrahedral remeshing - add a test (#8953)

## Summary of Changes

Add a test that does
```
  c3t3 = make_mesh_3(domain);
  write_MEDIT(file, c3t3);
  read_MEDIT(file, tr);
  tetrahedral_isotropic_remeshing(tr);
```

to test the pipeline.

Adding this test is triggered by issue #8948, but it does not fix it.


## Release Management

* Affected package(s): Tetrahedral_remeshing
* License and copyright ownership:
This commit is contained in:
Sebastien Loriot 2025-10-17 10:32:30 +02:00 committed by GitHub
commit 6f4bb669d6
No known key found for this signature in database
GPG Key ID: B5690EEEBB952194
3 changed files with 84 additions and 3 deletions

View File

@ -566,6 +566,7 @@ bool build_triangulation_from_file(std::istream& is,
{
using Point_3 = typename Tr::Point;
using Subdomain_index = typename Tr::Cell::Subdomain_index;
using Surface_patch_index = typename Tr::Cell::Surface_patch_index;
using Facet = std::array<int, 3>; // 3 = id
using Tet_with_ref = std::array<int, 4>; // 4 = id
@ -576,7 +577,7 @@ bool build_triangulation_from_file(std::istream& is,
std::vector<Tet_with_ref> finite_cells;
std::vector<Subdomain_index> subdomains;
std::vector<Point_3> points;
boost::unordered_map<Facet, typename Tr::Cell::Surface_patch_index> border_facets;
boost::unordered_map<Facet, Surface_patch_index> border_facets;
int dim;
int nv, nf, ntet, ref;
@ -632,7 +633,7 @@ bool build_triangulation_from_file(std::istream& is,
if(line.find("Triangles") != std::string::npos)
{
bool has_negative_surface_patch_ids = false;
typename Tr::Cell::Surface_patch_index max_surface_patch_id = 0;
Surface_patch_index max_surface_patch_id{0};
is >> nf;
if(verbose)
@ -641,7 +642,7 @@ bool build_triangulation_from_file(std::istream& is,
for(int i=0; i<nf; ++i)
{
int n[3];
typename Tr::Cell::Surface_patch_index surface_patch_id;
Surface_patch_index surface_patch_id;
if(!(is >> n[0] >> n[1] >> n[2] >> surface_patch_id))
{
if(verbose)

View File

@ -40,6 +40,9 @@ if(TARGET CGAL::Eigen3_support)
create_single_source_cgal_program("features_and_adaptive_sizing.cpp")
target_link_libraries(features_and_adaptive_sizing PRIVATE CGAL::Eigen3_support)
create_single_source_cgal_program("test_mesh_and_remesh_with_io.cpp")
target_link_libraries(test_mesh_and_remesh_with_io PRIVATE CGAL::Eigen3_support)
# with MLS projection
add_executable(test_mesh_and_remesh_polyhedron_with_features_mls
"test_mesh_and_remesh_polyhedron_with_features.cpp")

View File

@ -0,0 +1,77 @@
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
#include <CGAL/Mesh_criteria_3.h>
#include <CGAL/Mesh_triangulation_3.h>
#include <CGAL/Image_3.h>
#include <CGAL/Labeled_mesh_domain_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/Tetrahedral_remeshing/Remeshing_triangulation_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include <CGAL/IO/File_medit.h>
#include <string>
// Domain
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Image_3 Image;
typedef CGAL::Labeled_mesh_domain_3<K> Mesh_domain;
// Triangulation
typedef CGAL::Mesh_triangulation_3<Mesh_domain>::type Tr;
typedef CGAL::Mesh_complex_3_in_triangulation_3<Tr> C3t3;
// Mesh Criteria
typedef CGAL::Mesh_criteria_3<Tr> Mesh_criteria;
typedef Mesh_criteria::Facet_criteria Facet_criteria;
typedef Mesh_criteria::Cell_criteria Cell_criteria;
typedef Tr::Geom_traits Gt;
typedef CGAL::Triangulation_3<Gt, typename Tr::Triangulation_data_structure> T3;
typedef CGAL::Tetrahedral_remeshing::Remeshing_triangulation_3<Gt> Remeshing_triangulation;
using namespace CGAL::parameters;
int main()
{
const std::string filename = CGAL::data_file_path("images/liver.inr.gz");
CGAL::Image_3 image;
if(!image.read(filename)) {
std::cerr << "Error: Cannot read file " << filename << std::endl;
return EXIT_FAILURE;
}
Mesh_domain domain = Mesh_domain::create_labeled_image_mesh_domain(image, relative_error_bound = 1e-9);
// Mesh criteria
Facet_criteria facet_criteria(25, 20, 2); // angle, size, approximation
Cell_criteria cell_criteria(3, 20); // radius-edge ratio, size
Mesh_criteria criteria(facet_criteria, cell_criteria);
// Mesh
C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain, criteria, no_exude(), no_perturb());
std::cout << "Meshing done." << std::endl;
// Write
std::ofstream os("after_meshing_io.mesh");
CGAL::IO::write_MEDIT(os, c3t3);
os.close();
// Read
Remeshing_triangulation tr;
std::ifstream is("after_meshing_io.mesh");
CGAL::IO::read_MEDIT(is, tr);
// Remesh
double target_edge_length = 10.;
CGAL::tetrahedral_isotropic_remeshing(tr, target_edge_length);
std::cout << "Remeshing done." << std::endl;
return 0;
}