mirror of https://github.com/CGAL/cgal
extend APSS user manual (including a full example)
This commit is contained in:
parent
d11c7dd142
commit
a1c1269720
|
|
@ -4010,6 +4010,7 @@ Surface_reconstruction_points_3/doc_tex/Surface_reconstruction_points_3/output.t
|
|||
Surface_reconstruction_points_3/doc_tex/Surface_reconstruction_points_3/pipeline.eps -text svneol=unset#application/postscript
|
||||
Surface_reconstruction_points_3/doc_tex/Surface_reconstruction_points_3/pipeline.jpg -text svneol=unset#image/jpeg
|
||||
Surface_reconstruction_points_3/examples/Surface_reconstruction_points_3/APSS_reconstruction.cmd eol=lf
|
||||
Surface_reconstruction_points_3/examples/Surface_reconstruction_points_3/APSS_reconstruction_example.cpp -text
|
||||
Surface_reconstruction_points_3/examples/Surface_reconstruction_points_3/cgal_test_with_cmake -text
|
||||
Surface_reconstruction_points_3/examples/Surface_reconstruction_points_3/cgal_test_with_cmake.bat eol=crlf
|
||||
Surface_reconstruction_points_3/examples/Surface_reconstruction_points_3/data/camel-5kv.xyz -text svneol=unset#application/octet-stream
|
||||
|
|
|
|||
|
|
@ -1,29 +1,106 @@
|
|||
\section{Algebraic Point Set Surfaces}
|
||||
|
||||
The Algebraic Point Set Surfaces \cite{Guennebaud07} falls into the category of so-called point set surfaces, where the approximate signed distance function to the inferred surface can be evaluated at any query point, on the fly during iso-contouring: \\
|
||||
\ccc{CGAL::APSS_reconstruction_function<GeomTraits>} \\
|
||||
|
||||
% Insert image APSS.jpg/eps
|
||||
\begin{center}
|
||||
\label{Surface_reconstruction_points_3-fig-APSS}
|
||||
% Image
|
||||
\begin{ccTexOnly}
|
||||
\includegraphics[width=0.5\textwidth]{Surface_reconstruction_points_3/APSS} % omit .eps suffix
|
||||
\end{ccTexOnly}
|
||||
\begin{ccHtmlOnly}
|
||||
<img width="50%" border=0 src="./APSS.jpg"><P>
|
||||
\end{ccHtmlOnly}
|
||||
% Title
|
||||
\begin{figure}[h]
|
||||
\caption{APSS surface reconstruction from 10K
|
||||
points sampled with a Minolta laser scanner.}
|
||||
% later: add computed function
|
||||
\end{figure}
|
||||
\end{center}
|
||||
|
||||
|
||||
\emph{more to come here about interface and example}
|
||||
|
||||
%\subsection{Interface}
|
||||
%\subsection{Example}
|
||||
|
||||
|
||||
\section{Algebraic Point Set Surfaces}
|
||||
|
||||
The Algebraic Point Set Surfaces \cite{Guennebaud07} method defines an approximate signed distance field from a set of oriented points sampled onto a surface, such that the zero-isosurface best approximates the input points with respect to some given criteria. It uses local Moving Least Squares (MLS) fitting of algebaic spheres. The implicit function can be evaluated for any point in space with almost no preprocessing.
|
||||
|
||||
More precisely, the value of the implicit function $f(x)$ is obtained by fitting an algebraic sphere to the neighbors of the point $x$ in a weighted least square sense where the weights depends on both the Euclidean distance between $x$ and the neighbor points. This weight function behave like a low pass filter allowing to control the degree of smoothing. Intuitively, the fitting algorithm first tries to minimize the distance between the gradient of the algebraic sphere and the input normals, and then it minimizes the algebraic distance between the 0-isosphere and the neighbors to somehow {\em set} the radius of the sphere. Therefore the method is quite sensible to the quality of the input normals, and sampling.
|
||||
|
||||
Even though the \ccc{CGAL::APSS_reconstruction_function} class can be directly used, it is aimed to be used with the \cgal\ surface mesh generator to extract the $0$-isosurface.
|
||||
|
||||
|
||||
% Insert image APSS.jpg/eps
|
||||
\begin{center}
|
||||
\label{Surface_reconstruction_points_3-fig-APSS}
|
||||
% Image
|
||||
\begin{ccTexOnly}
|
||||
\includegraphics[width=0.5\textwidth]{Surface_reconstruction_points_3/APSS} % omit .eps suffix
|
||||
\end{ccTexOnly}
|
||||
\begin{ccHtmlOnly}
|
||||
<img width="50%" border=0 src="./APSS.jpg"><P>
|
||||
\end{ccHtmlOnly}
|
||||
% Title
|
||||
\begin{figure}[h]
|
||||
\caption{APSS surface reconstruction from 10K
|
||||
points sampled with a Minolta laser scanner.}
|
||||
% later: add computed function
|
||||
\end{figure}
|
||||
\end{center}
|
||||
|
||||
\subsection{Interface}
|
||||
|
||||
The class template declaration is:
|
||||
|
||||
template$<$ \\
|
||||
class Gt$>$ \\
|
||||
class \ccc{APSS_reconstruction_function};
|
||||
|
||||
with \\
|
||||
\ccc{Gt}: Geometric traits class.
|
||||
|
||||
|
||||
Creation:
|
||||
|
||||
% Reduce left margin
|
||||
\ccTwo{1234567890123456789012}{}
|
||||
|
||||
\ccFunction{template<typename InputIterator, typename PointPMap, typename NormalPMap, typename RadiusPMap> APSS_reconstruction_function(InputIterator first, InputIterator beyond, PointPMap point_pmap, NormalPMap normal_pmap, RadiusPMap radius_pmap, FT smoothness);}
|
||||
{
|
||||
Creates a APSS implicit function from the [first, beyond) range of points.
|
||||
\ccCommentHeading{Template Parameters} \\
|
||||
\ccc{InputIterator}: iterator over input points.
|
||||
|
||||
\ccc{PointPMap}: is a model of \ccc{boost::ReadablePropertyMap} with a \ccc{value_type} = \ccc{Point_3}. It can be omitted if \ccc{InputIterator} \ccc{value_type} is convertible to \ccc{Point_3}.
|
||||
|
||||
\ccc{NormalPMap}: is a model of \ccc{boost::ReadablePropertyMap} with a \ccc{value_type} = \ccc{Vector_3}.
|
||||
|
||||
\ccc{RadiusPMap}: is a model of \ccc{boost::ReadablePropertyMap} with a \ccc{value_type} = \ccc{FT}. It is optional (see below).
|
||||
|
||||
\ccCommentHeading{Parameters} \\
|
||||
\ccc{first}: iterator over the first input point.
|
||||
|
||||
\ccc{beyond}: past-the-end iterator over the input points.
|
||||
|
||||
\ccc{point_pmap}: property map to access the position of an input point.
|
||||
|
||||
\ccc{normal_pmap}: property map to access the {\bf oriented} normal of an input point.
|
||||
|
||||
\ccc{radius_pmap}: property map to access the influence radius of an input point. It can be omitted, in which case the influence radius of each point is automatically computed from a basic estimate of the local density using the 16 nearest neighbors.
|
||||
}
|
||||
|
||||
|
||||
The main operations are:
|
||||
|
||||
\ccFunction{Sphere bounding_sphere() const;}
|
||||
{
|
||||
Returns a sphere bounding the inferred surface.
|
||||
}
|
||||
\ccGlue
|
||||
\ccFunction{FT value(const Point& p, bool* ok = 0) const;}
|
||||
{
|
||||
Evaluates the implicit function at a given 3D query point. The optional parameter \ccc{ok} returns true if the evaluation succeeded and false otherwise.
|
||||
}
|
||||
% \ccGlue
|
||||
% \ccFunction{void project(Point& p, Vector& n, unsigned int maxNofIterations = 20) const;}
|
||||
% {
|
||||
% Iteratively projects a point onto the underlying implicit surface until convergence or a maximal number of iteration (\ccc{maxNofIterations}). The optional parameter \ccc{n} returns the surface normal at the projection point.
|
||||
% }
|
||||
\ccGlue
|
||||
\ccFunction{Point get_inner_point() const;}
|
||||
{
|
||||
Returns a point located inside the inferred surface.
|
||||
}
|
||||
|
||||
|
||||
See details in \\
|
||||
\ccc{CGAL::APSS_reconstruction_function<GeomTraits>}
|
||||
|
||||
|
||||
% falls into the category of so-called point set surfaces, where the approximate signed distance function to the inferred surface can be evaluated at any query point, on the fly during iso-contouring: \\
|
||||
% \ccc{CGAL::APSS_reconstruction_function<GeomTraits>} \\
|
||||
|
||||
|
||||
\subsection{Example}
|
||||
|
||||
\ccc{APSS_reconstruction_example.cpp} reads a point set, creates a APSS implicit function and extract the $0$-isosurface as a triangle mesh.
|
||||
|
||||
\ccIncludeExampleCode{Surface_reconstruction_points_3/APSS_reconstruction_example.cpp}
|
||||
|
|
|
|||
|
|
@ -0,0 +1,83 @@
|
|||
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
|
||||
#include <CGAL/Surface_mesh_default_triangulation_3.h>
|
||||
#include <CGAL/make_surface_mesh.h>
|
||||
#include <CGAL/Implicit_surface_3.h>
|
||||
#include <CGAL/IO/Complex_2_in_triangulation_3_file_writer.h>
|
||||
#include <CGAL/APSS_reconstruction_function.h>
|
||||
#include <CGAL/Point_with_normal_3.h>
|
||||
#include <CGAL/property_map.h>
|
||||
#include <CGAL/IO/read_xyz_points.h>
|
||||
|
||||
// Types
|
||||
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
|
||||
typedef Kernel::FT FT;
|
||||
typedef Kernel::Point_3 Point;
|
||||
typedef Kernel::Vector_3 Vector;
|
||||
typedef CGAL::Point_with_normal_3<Kernel> Point_with_normal;
|
||||
typedef Kernel::Sphere_3 Sphere;
|
||||
typedef std::deque<Point_with_normal> PointList;
|
||||
typedef CGAL::Polyhedron_3<Kernel> Polyhedron;
|
||||
typedef CGAL::APSS_reconstruction_function<Kernel> APSS_reconstruction_function;
|
||||
typedef CGAL::Surface_mesh_default_triangulation_3 STr;
|
||||
typedef CGAL::Surface_mesh_complex_2_in_triangulation_3<STr> C2t3;
|
||||
typedef CGAL::Implicit_surface_3<Kernel, APSS_reconstruction_function> Surface_3;
|
||||
|
||||
int main(int argc, char * argv[])
|
||||
{
|
||||
// Reads the point set file in points[].
|
||||
// Note: read_xyz_points_and_normals() requires an iterator over points
|
||||
// + property maps to access each point's position and normal.
|
||||
// The position property map can be omitted here as we use iterators over Point_3 elements.
|
||||
PointList points;
|
||||
std::ifstream stream("data/oni.xyz");
|
||||
if (!stream ||
|
||||
!CGAL::read_xyz_points_and_normals(
|
||||
stream,
|
||||
std::back_inserter(points),
|
||||
CGAL::make_normal_of_point_with_normal_pmap(std::back_inserter(points))))
|
||||
{
|
||||
return EXIT_FAILURE;
|
||||
}
|
||||
|
||||
// Creates implicit function from the read points.
|
||||
// Note: APSS_reconstruction_function() requires an iterator over points
|
||||
// + property maps to access each point's position and normal.
|
||||
// The position property map can be omitted here as we use iterators over Point_3 elements.
|
||||
APSS_reconstruction_function implicit_function(
|
||||
points.begin(), points.end(),
|
||||
CGAL::make_normal_of_point_with_normal_pmap(points.begin()),
|
||||
2.0); // Smoothness factor. In the rangcgal_manuale 2 (clean datasets) and 8 (noisy datasets).
|
||||
|
||||
// Gets point inside the implicit surface
|
||||
Point inner_point = implicit_function.get_inner_point();
|
||||
|
||||
// Gets implicit function's radius
|
||||
Sphere bsphere = implicit_function.bounding_sphere();
|
||||
FT size = sqrt(bsphere.squared_radius());
|
||||
|
||||
// defining the implicit surface = implicit function + bounding sphere centered at inner_point
|
||||
Point sm_sphere_center = inner_point;
|
||||
FT sm_sphere_radius = size + std::sqrt(CGAL::squared_distance(bsphere.center(),inner_point));
|
||||
sm_sphere_radius *= 1.01; // make sure that the bounding sphere contains the surface
|
||||
Surface_3 surface(implicit_function,
|
||||
Sphere(sm_sphere_center,sm_sphere_radius*sm_sphere_radius));
|
||||
|
||||
// defining meshing criteria
|
||||
CGAL::Surface_mesh_default_criteria_3<STr> criteria(20.0, // Min triangle angle (degrees). 20=fast, 30 guaranties convergence.
|
||||
0.1*size, // Max triangle radius w.r.t. point set radius. 0.1 is fine.
|
||||
0.003*size); // Approximation error w.r.t. p.s.r.. For APSS: 0.015=fast, 0.003=smooth.
|
||||
|
||||
// meshing surface
|
||||
STr tr; // 3D-Delaunay triangulation for Surface Mesher
|
||||
C2t3 surface_mesher_c2t3 (tr); // 2D-complex in 3D-Delaunay triangulation
|
||||
CGAL::make_surface_mesh(surface_mesher_c2t3, // reconstructed mesh
|
||||
surface, // implicit surface
|
||||
criteria, // meshing criteria
|
||||
CGAL::Manifold_with_boundary_tag()); // require manifold mesh
|
||||
|
||||
// save the mesh
|
||||
std::ofstream out("oni_apss.off");
|
||||
CGAL::output_surface_facets_to_off(out, surface_mesher_c2t3);
|
||||
|
||||
return EXIT_SUCCESS;
|
||||
}
|
||||
|
|
@ -58,6 +58,7 @@ if ( CGAL_FOUND )
|
|||
|
||||
# Executables that do *not* require BLAS, LAPACK nor TAUCS
|
||||
create_single_source_cgal_program( "APSS_reconstruction.cpp" )
|
||||
create_single_source_cgal_program( "APSS_reconstruction_example.cpp" )
|
||||
|
||||
# Link with BLAS, LAPACK and TAUCS (optional)
|
||||
find_package(TAUCS)
|
||||
|
|
|
|||
|
|
@ -198,8 +198,6 @@ private:
|
|||
Point p = get(point_pmap,it);
|
||||
m->bounding_box = m->bounding_box + BBox(p.x(), p.y(), p.z(), p.x(), p.y(), p.z());
|
||||
}
|
||||
std::cerr << m->bounding_box.xmin() << " " << m->bounding_box.ymin() << " " << m->bounding_box.zmin() << " "
|
||||
<< m->bounding_box.xmax() << " " << m->bounding_box.ymax() << " " << m->bounding_box.zmax() << "\n";
|
||||
|
||||
// Find a point inside the surface.
|
||||
find_inner_point();
|
||||
|
|
|
|||
Loading…
Reference in New Issue