\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.13 - Estimation of Local Differential Properties of Point-Sampled Surfaces
CGAL::Monge_via_jet_fitting< DataKernel, LocalKernel, SvdTraits > Class Template Reference

#include <CGAL/Monge_via_jet_fitting.h>

Definition

The class Monge_via_jet_fitting is designed to perform the estimation of the local differential quantities at a given point.

The point range is given by a pair of input iterators, and it is assumed that the point where the calculation is carried out is the point that the begin iterator refers to. The results are stored in an instance of the nested class Monge_form, the particular information returned depending on the degrees specified for the polynomial fitting and for the Monge form.

If CGAL_EIGEN3_ENABLED is defined, LocalKernel and SvdTraits template parameters have defaults, Simple_cartesian<double> and Eigen_svd respectively.

Template Parameters
DataKernelprovides the geometric classes and tools corresponding to the input points, and also members of the Monge_form class.
LocalKernelprovides the geometric classes and tools required by local computations.
SvdTraitsfeatures the linear algebra algorithm required by the fitting method. The scalar type, SvdTraits::FT, must be the same as that of the LocalKernel concept : LocalKernel::FT.
See also
Eigen_svd
Monge_form
Note
This class requires the Eigen library.
Examples:
Jet_fitting_3/Mesh_estimation.cpp, and Jet_fitting_3/Single_estimation.cpp.

Classes

class  Monge_form
 The class Monge_form stores the Monge representation, i.e., the Monge coordinate system and the coefficients of the Monge form in this system. More...
 

Types

typedef DataKernel Data_kernel
 
typedef LocalKernel Local_kernel
 
typedef Local_kernel::FT FT
 
typedef Local_kernel::Vector_3 Vector_3
 
typedef unspecified_type Monge_form
 see the page Monge_via_jet_fitting::Monge_form.
 

Creation

 Monge_via_jet_fitting ()
 default constructor
 

Operations

template<class InputIterator >
Monge_form operator() (InputIterator begin, InputIterator end, size_t d, size_t d')
 This operator performs all the computations. More...
 
FT condition_number ()
 condition number of the linear fitting system.
 
std::pair< FT, Vector_3pca_basis (size_t i)
 pca eigenvalues and eigenvectors, the pca_basis has always 3 such pairs. More...
 

Member Function Documentation

◆ operator()()

template<typename DataKernel , typename LocalKernel , typename SvdTraits >
template<class InputIterator >
Monge_form CGAL::Monge_via_jet_fitting< DataKernel, LocalKernel, SvdTraits >::operator() ( InputIterator  begin,
InputIterator  end,
size_t  d,
size_t d'   
)

This operator performs all the computations.

The \( N\) input points are given by the InputIterator parameters which value-type are Data_kernel::Point_3, d is the degree of the fitted polynomial, d' is the degree of the expected Monge coefficients.

Precondition
\( N \geq N_{d}:=(d+1)(d+2)/2\), \( 1 \leq d' \leq\min(d,4) \).

◆ pca_basis()

template<typename DataKernel , typename LocalKernel , typename SvdTraits >
std::pair<FT, Vector_3> CGAL::Monge_via_jet_fitting< DataKernel, LocalKernel, SvdTraits >::pca_basis ( size_t  i)

pca eigenvalues and eigenvectors, the pca_basis has always 3 such pairs.

Precondition
\( i\) ranges from 0 to 2.