CGAL 6.0 - 3D Isosurfacing
Loading...
Searching...
No Matches
Isosurfacing_3/marching_cubes.cpp
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Isosurfacing_3/Cartesian_grid_3.h>
#include <CGAL/Isosurfacing_3/Value_function_3.h>
#include <CGAL/Isosurfacing_3/marching_cubes_3.h>
#include <CGAL/Isosurfacing_3/Marching_cubes_domain_3.h>
#include <CGAL/IO/polygon_soup_io.h>
#include <vector>
using FT = typename Kernel::FT;
using Point = typename Kernel::Point_3;
using Vector = typename Kernel::Vector_3;
using Point_range = std::vector<Point>;
using Polygon_range = std::vector<std::vector<std::size_t> >;
int main(int argc, char** argv)
{
const FT isovalue = (argc > 1) ? std::stod(argv[1]) : 0.8;
// create bounding box and grid
const CGAL::Bbox_3 bbox { -1., -1., -1., 1., 1., 1. };
Grid grid { bbox, CGAL::make_array<std::size_t>(30, 30, 30) };
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;
// fill up values
auto sphere_value_fn = [](const Point& p) -> FT
{
return sqrt(p.x()*p.x() + p.y()*p.y() + p.z()*p.z());
};
Values values { sphere_value_fn, grid };
// Below is equivalent to:
// Domain domain { grid, values };
Point_range points;
Polygon_range triangles;
// run marching cubes isosurfacing
std::cout << "Running Marching Cubes with isovalue = " << isovalue << std::endl;
CGAL::Isosurfacing::marching_cubes<CGAL::Parallel_if_available_tag>(domain, isovalue, points, triangles,
CGAL::parameters::use_topologically_correct_marching_cubes(true));
std::cout << "Output #vertices: " << points.size() << std::endl;
std::cout << "Output #triangles: " << triangles.size() << std::endl;
CGAL::IO::write_polygon_soup("marching_cubes.off", points, triangles);
std::cout << "Done" << std::endl;
return EXIT_SUCCESS;
}
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
A domain that can be used with the Marching Cubes algorithm.
Definition: Marching_cubes_domain_3.h:50
The class Value_function_3 represents a field of scalars computed using a user-provided unary functio...
Definition: Value_function_3.h:40
bool write_polygon_soup(const std::string &fname, const PointRange &points, const PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
Marching_cubes_domain_3< Partition, ValueField, EdgeIntersectionOracle > create_marching_cubes_domain_3(const Partition &partition, const ValueField &values, const EdgeIntersectionOracle &intersection_oracle=EdgeIntersectionOracle())
creates a new instance of a domain that can be used with the Marching Cubes algorithm.
Definition: Marching_cubes_domain_3.h:92