#include <CGAL/Simple_cartesian.h>
#include <CGAL/Isosurfacing_3/Cartesian_grid_3.h>
#include <CGAL/Isosurfacing_3/dual_contouring_3.h>
#include <CGAL/Isosurfacing_3/Dual_contouring_domain_3.h>
#include <CGAL/Isosurfacing_3/Finite_difference_gradient_3.h>
#include <CGAL/Isosurfacing_3/Interpolated_discrete_values_3.h>
#include <CGAL/Isosurfacing_3/marching_cubes_3.h>
#include <CGAL/Isosurfacing_3/Marching_cubes_domain_3.h>
#include <CGAL/Image_3.h>
#include <CGAL/Isosurfacing_3/IO/Image_3.h>
#include <CGAL/IO/polygon_soup_io.h>
#include <iostream>
#include <vector>
using FT = typename Kernel::FT;
using Point_range = std::vector<Point>;
using Polygon_range = std::vector<std::vector<std::size_t> >;
void run_marching_cubes(const Grid& grid,
const FT isovalue,
const Values& values)
{
using Domain = IS::Marching_cubes_domain_3<Grid, Values, IS::Linear_interpolation_edge_intersection>;
std::cout << "\n ---- " << std::endl;
std::cout << "Running Marching Cubes with isovalue = " << isovalue << std::endl;
Domain domain { grid, values };
Point_range points;
Polygon_range triangles;
IS::marching_cubes<CGAL::Parallel_if_available_tag>(domain, isovalue, points, triangles);
std::cout << "Output #vertices: " << points.size() << std::endl;
std::cout << "Output #triangles: " << triangles.size() << std::endl;
}
void run_dual_contouring(const Grid& grid,
const FT isovalue,
const Values& values)
{
using Domain = IS::Dual_contouring_domain_3<Grid, Values, Gradients, IS::Linear_interpolation_edge_intersection>;
std::cout << "\n ---- " << std::endl;
std::cout << "Running Dual Contouring with isovalue = " << isovalue << std::endl;
const FT step = CGAL::approximate_sqrt(grid.spacing().squared_length()) * 0.01;
Gradients gradients { values, step };
Domain domain { grid, values, gradients };
Point_range points;
Polygon_range triangles;
IS::dual_contouring<CGAL::Parallel_if_available_tag>(domain, isovalue, points, triangles);
std::cout << "Output #vertices: " << points.size() << std::endl;
std::cout << "Output #triangles: " << triangles.size() << std::endl;
}
int main(int argc, char* argv[])
{
const FT isovalue = (argc > 2) ? std::stod(argv[2]) : - 2.9;
{
std::cerr << "Error: Cannot read image file " << fname << std::endl;
return EXIT_FAILURE;
}
Grid grid;
Values values { grid };
if(!IS::IO::convert_image_to_grid(image, grid, values))
{
std::cerr << "Error: Cannot convert image to Cartesian grid" << std::endl;
return EXIT_FAILURE;
}
std::cout << "Span: " << grid.span() << std::endl;
std::cout << "Cell dimensions: " << grid.spacing()[0] << " " << grid.spacing()[1] << " " << grid.spacing()[2] << std::endl;
std::cout << "Cell #: " << grid.xdim() << ", " << grid.ydim() << ", " << grid.zdim() << std::endl;
run_marching_cubes(grid, isovalue, values);
run_dual_contouring(grid, isovalue, values);
return EXIT_SUCCESS;
}
bool read(const char *file)
The class Cartesian_grid_3 represents a 3D Cartesian grid, that is the partition of an iso-cuboid int...
Definition: Cartesian_grid_3.h:148
Class template for a gradient that is calculated using finite differences.
Definition: Finite_difference_gradient_3.h:39
Class template for a field of values that are calculated using discrete values and interpolation.
Definition: Interpolated_discrete_values_3.h:41
bool write_polygon_soup(const std::string &fname, const PointRange &points, const PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
Definition: partition_traits.h:2
std::string data_file_path(const std::string &filename)