Ponca  aa50bfdf187919869239c5b44b748842569114c1
Point Cloud Analysis library
Loading...
Searching...
No Matches
Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T > Class Template Reference

Algebraic Sphere fitting procedure on point sets with non-oriented normals. More...

#include <unorientedSphereFit.h>

+ Inheritance diagram for Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >:
+ Collaboration diagram for Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >:

Public Types

using Scalar = typename DataPoint::Scalar
 Alias to scalar type.
 
using VectorType = typename Base::VectorType
 Alias to vector type.
 
using WFunctor = typename Base::WFunctor
 Alias to weight function.
 
using Solver = Eigen::EigenSolver< MatrixBB >
 

Public Member Functions

UnorientedSphereFitImpl< DataPoint, _WFunctor, T > & unorientedSphereFit ()
 Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.
 
const UnorientedSphereFitImpl< DataPoint, _WFunctor, T > & unorientedSphereFit () const
 Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.
 
void init (const VectorType &_evalPos)
 Set the evaluation position and reset the internal states.
 
bool addLocalNeighbor (Scalar w, const VectorType &localQ, const DataPoint &attributes)
 Add a neighbor to perform the fit.
 
FIT_RESULT finalize ()
 Finalize the procedure.
 
constexpr bool isSigned ()
 Is scalar field signed.
 

Public Attributes

MatrixBB m_matA {MatrixBB::Zero()}
 The accumulated covariance matrix.
 
MatrixBB m_matQ {MatrixBB::Zero()}
 The constraint matrix.
 
Scalar m_sumDotPP {0}
 Sum of the squared relative positions.
 
Solver m_solver
 

Protected Types

enum  { Check = Base::PROVIDES_ALGEBRAIC_SPHERE && Base::PROVIDES_MEAN_POSITION }
 
using Base = T
 Base class of the procedure.
 
typedef Eigen::Matrix< Scalar, DataPoint::Dim+1, 1 > VectorB
 
typedef Eigen::Matrix< Scalar, DataPoint::Dim+1, DataPoint::Dim+1 > MatrixBB
 

Detailed Description

template<class DataPoint, class _WFunctor, typename T>
class Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >

Algebraic Sphere fitting procedure on point sets with non-oriented normals.

This method published in [4] maximizes the sum of squared dot product between the input normal vectors and the gradient of the algebraic sphere. The maximization is done under the constraint that the norm of the gradient is unitary on average. In practice, it amounts to solve the generalized eigenvalue problem

\[ A \mathbf{u} = \lambda Q \mathbf{u} \]

where

\[ \mathbf{u} = \begin{bmatrix} u_l \\ u_q \end{bmatrix} \]

\[ A = \sum_i w_i \begin{bmatrix} \mathbf{n}_i \\ \mathbf{n}_i^T\mathbf{p}_i \end{bmatrix}\begin{bmatrix} \mathbf{n}_i^T & \mathbf{n}_i^T\mathbf{p}_i \end{bmatrix} \]

\[ Q = \frac{1}{\sum_i w_i} \begin{bmatrix} I & \sum_i w_i \mathbf{p}_i \\ \sum_i w_i \mathbf{p}_i^T & \sum_i w_i \mathbf{p}_i^T\mathbf{p}_i \end{bmatrix} \]

The constant coefficient \(u_c\) is computed as in OrientedSphereFitImpl by minimizing the sum of squared potential evaluated at the points \(\mathbf{p}_i\). This fitting method corresponds to the case where the scalar field is \( f^* \) in the Section 3.3 of [4].

See also
class AlgebraicSphere, class OrientedSphereFit

Definition at line 45 of file unorientedSphereFit.h.

Member Typedef Documentation

◆ Base

template<class DataPoint , class _WFunctor , typename T >
using Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::Base = T
protected

Base class of the procedure.

Definition at line 47 of file unorientedSphereFit.h.

◆ MatrixBB

template<class DataPoint , class _WFunctor , typename T >
typedef Eigen::Matrix<Scalar, DataPoint::Dim+1, DataPoint::Dim+1> Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::MatrixBB
protected

Definition at line 53 of file unorientedSphereFit.h.

◆ Scalar

template<class DataPoint , class _WFunctor , typename T >
using Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::Scalar = typename DataPoint::Scalar

Alias to scalar type.

Definition at line 47 of file unorientedSphereFit.h.

◆ Solver

template<class DataPoint , class _WFunctor , typename T >
using Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::Solver = Eigen::EigenSolver<MatrixBB>

Definition at line 56 of file unorientedSphereFit.h.

◆ VectorB

template<class DataPoint , class _WFunctor , typename T >
typedef Eigen::Matrix<Scalar, DataPoint::Dim+1, 1> Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::VectorB
protected

Definition at line 52 of file unorientedSphereFit.h.

◆ VectorType

template<class DataPoint , class _WFunctor , typename T >
using Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::VectorType = typename Base::VectorType

Alias to vector type.

Definition at line 47 of file unorientedSphereFit.h.

◆ WFunctor

template<class DataPoint , class _WFunctor , typename T >
using Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::WFunctor = typename Base::WFunctor

Alias to weight function.

Definition at line 47 of file unorientedSphereFit.h.

Member Enumeration Documentation

◆ anonymous enum

template<class DataPoint , class _WFunctor , typename T >
anonymous enum
protected

Definition at line 50 of file unorientedSphereFit.h.

Member Function Documentation

◆ addLocalNeighbor()

template<class DataPoint , class _WFunctor , typename T >
bool Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::addLocalNeighbor ( Scalar  w,
const VectorType localQ,
const DataPoint &  attributes 
)
inline

Add a neighbor to perform the fit.

Returns
false if param nei is not a valid neighbour (weight = 0)

Definition at line 22 of file unorientedSphereFit.hpp.

◆ finalize()

template<class DataPoint , class _WFunctor , typename T >
FIT_RESULT Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::finalize
inline

Finalize the procedure.

Returns
Fitting Status
Warning
Must be called be for any use of the fitting output

Definition at line 41 of file unorientedSphereFit.hpp.

◆ init()

template<class DataPoint , class _WFunctor , typename T >
void Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::init ( const VectorType _evalPos)
inline

Set the evaluation position and reset the internal states.

Warning
Must be called be for any computation (but after #setWeightFunc)

Definition at line 12 of file unorientedSphereFit.hpp.

◆ isSigned()

template<class DataPoint , class _WFunctor , typename T >
constexpr bool Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::isSigned ( )
inlineconstexpr

Is scalar field signed.

If not, the method the sign of potential() must be ignored

Definition at line 67 of file unorientedSphereFit.h.

◆ unorientedSphereFit() [1/2]

template<class DataPoint , class _WFunctor , typename T >
UnorientedSphereFitImpl< DataPoint, _WFunctor, T > & Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::unorientedSphereFit ( )
inline

Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.

Definition at line 65 of file unorientedSphereFit.h.

◆ unorientedSphereFit() [2/2]

template<class DataPoint , class _WFunctor , typename T >
const UnorientedSphereFitImpl< DataPoint, _WFunctor, T > & Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::unorientedSphereFit ( ) const
inline

Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.

Definition at line 65 of file unorientedSphereFit.h.

Member Data Documentation

◆ m_matA

template<class DataPoint , class _WFunctor , typename T >
MatrixBB Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::m_matA {MatrixBB::Zero()}

The accumulated covariance matrix.

Definition at line 58 of file unorientedSphereFit.h.

◆ m_matQ

template<class DataPoint , class _WFunctor , typename T >
MatrixBB Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::m_matQ {MatrixBB::Zero()}

The constraint matrix.

Definition at line 59 of file unorientedSphereFit.h.

◆ m_solver

template<class DataPoint , class _WFunctor , typename T >
Solver Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::m_solver

Definition at line 62 of file unorientedSphereFit.h.

◆ m_sumDotPP

template<class DataPoint , class _WFunctor , typename T >
Scalar Ponca::UnorientedSphereFitImpl< DataPoint, _WFunctor, T >::m_sumDotPP {0}

Sum of the squared relative positions.

Definition at line 60 of file unorientedSphereFit.h.