Julian Stahl and Daniel Zint
This package implements isosurfacing algorithms to generate a mesh from a scalar field. The algorithms provide different guarantees for the output and should be chosen depending on the use-case. The representations of input data is flexible and includes explicit and implicit formats.
This package provides algorithms to extract isosurfaces from different inputs. The input is represented as a domain and can be an implicit function, a cartesion grid, or an octree. The output is an indexed face set that stores an isosurface in the form of a surface mesh. Available algorithms include Marching Cubes, topologically correct Marching Cubes, and Dual Contouring.
Concepts
Classes
Cartesian_grid_3
Zero_gradient
Finite_difference_gradient
Explicit_cartesian_grid_gradient
Free Functions
|
template<typename Concurrency_tag = Sequential_tag, class Domain_ , class PointRange , class PolygonRange , class Positioning = internal::Positioning::QEM_SVD<true>> |
void | CGAL::Isosurfacing::dual_contouring (const Domain_ &domain, const typename Domain_::FT iso_value, PointRange &points, PolygonRange &polygons, const Positioning &positioning=Positioning()) |
| Creates an indexed face set that represents an isosurface using the Dual Contouring algorithm. More...
|
|
template<class GeomTraits , typename Gradient_ = Zero_gradient<GeomTraits>> |
Explicit_cartesian_grid_domain< GeomTraits, Gradient_ > | CGAL::Isosurfacing::create_explicit_cartesian_grid_domain (const std::shared_ptr< Cartesian_grid_3< GeomTraits >> grid, const Gradient_ &gradient=Gradient_()) |
| Creates a domain from a Cartesian_grid_3 that can be used as input for isosurfacing algorithms. More...
|
|
template<class GeomTraits , typename PointFunction , typename Gradient_ = Zero_gradient<GeomTraits>> |
Implicit_cartesian_grid_domain< GeomTraits, PointFunction, Gradient_ > | CGAL::Isosurfacing::create_implicit_cartesian_grid_domain (const Bbox_3 &bbox, const typename GeomTraits::Vector_3 &spacing, const PointFunction &point_function, const Gradient_ &gradient=Gradient_()) |
| Creates a domain from a Cartesian_grid_3 that can be used as input for isosurfacing algorithms. More...
|
|
template<typename Concurrency_tag = Sequential_tag, class Domain_ , class PointRange , class TriangleRange > |
void | CGAL::Isosurfacing::marching_cubes (const Domain_ &domain, const typename Domain_::FT iso_value, PointRange &points, TriangleRange &polygons, bool topologically_correct=true) |
| Creates a polygon soup that represents an isosurface using the marching cubes algorithm. More...
|
|
◆ create_explicit_cartesian_grid_domain()
template<class GeomTraits , typename Gradient_ = Zero_gradient<GeomTraits>>
#include <CGAL/Explicit_cartesian_grid_domain.h>
Creates a domain from a Cartesian_grid_3 that can be used as input for isosurfacing algorithms.
- Template Parameters
-
GeomTraits | the traits type |
- Parameters
-
grid | the cartesian grid containing input data |
gradient | a function that describes the gradient of the data |
◆ create_implicit_cartesian_grid_domain()
template<class GeomTraits , typename PointFunction , typename Gradient_ = Zero_gradient<GeomTraits>>
Implicit_cartesian_grid_domain<GeomTraits, PointFunction, Gradient_> CGAL::Isosurfacing::create_implicit_cartesian_grid_domain |
( |
const Bbox_3 & |
bbox, |
|
|
const typename GeomTraits::Vector_3 & |
spacing, |
|
|
const PointFunction & |
point_function, |
|
|
const Gradient_ & |
gradient = Gradient_() |
|
) |
| |
#include <CGAL/Implicit_cartesian_grid_domain.h>
Creates a domain from a Cartesian_grid_3 that can be used as input for isosurfacing algorithms.
- Template Parameters
-
GeomTraits | the traits type |
PointFunction | the type of the implicit function |
- Parameters
-
bbox | a bounding box that specifies the size of the functions domain |
spacing | the distance between discretized points on the function |
point_function | the function with a point as argument |
gradient | a function that describes the gradient of the data |
◆ dual_contouring()
template<typename Concurrency_tag = Sequential_tag, class Domain_ , class PointRange , class PolygonRange , class Positioning = internal::Positioning::QEM_SVD<true>>
void CGAL::Isosurfacing::dual_contouring |
( |
const Domain_ & |
domain, |
|
|
const typename Domain_::FT |
iso_value, |
|
|
PointRange & |
points, |
|
|
PolygonRange & |
polygons, |
|
|
const Positioning & |
positioning = Positioning() |
|
) |
| |
◆ marching_cubes()
template<typename Concurrency_tag = Sequential_tag, class Domain_ , class PointRange , class TriangleRange >
void CGAL::Isosurfacing::marching_cubes |
( |
const Domain_ & |
domain, |
|
|
const typename Domain_::FT |
iso_value, |
|
|
PointRange & |
points, |
|
|
TriangleRange & |
polygons, |
|
|
bool |
topologically_correct = true |
|
) |
| |