Helios++
Helios software for LiDAR simulations
SurfaceInspector::maths::PlaneFitter Class Reference

Handle plane fitting operations. More...

#include <PlaneFitter.hpp>

Inheritance diagram for SurfaceInspector::maths::PlaneFitter:
Collaboration diagram for SurfaceInspector::maths::PlaneFitter:

Static Public Member Functions

template<typename T >
static vector< T > centerCoordinatesMatrix (Mat< T > &M)
 Modify coordinates at M so they are centered at origin. More...
 
template<typename T >
static vector< T > translateToOrigin (Mat< T > &M, vector< T > center)
 Translate M to origin transposing by center point. More...
 
template<typename T >
static Plane< T > bestFittingPlaneSVD (Mat< T > &M)
 Compute the best fitting plane for given Matrix of coordinates through singular value decomposition. More...
 
template<typename T >
static Plane< T > bestFittingPlanePCA (Mat< T > &M)
 Compute the best fitting plane for given Matrix of coordinates through principal component analysis. More...
 
template<typename T >
static Plane< T > bestFittingPlaneFromCovariances (Mat< T > &M)
 Compute the best fitting plane for given Matrix of covariances through eigen value decomposition. More...
 
template<typename T >
static DetailedPlane< T > bestFittingDetailedPlaneSVD (Mat< T > &M, vector< T > center)
 Like bestFittingPlaneSVD but returning a DetailedPlane which provides more information than a simple Plane object. More...
 
template<typename T >
static DetailedPlane< T > bestFittingDetailedPlanePCA (Mat< T > &M, vector< T > center)
 Like bestFittingPlanePCA but returning a DetailedPlane which provides more information than a simple Plane object. More...
 

Static Protected Member Functions

template<typename T >
static void extractCommonPlaneComponents (arma::mat const &A, arma::vec const &v, T &sum, Plane< T > &p)
 Extract common plane components. More...
 
template<typename T >
static void extractDetailedPlaneComponents (Mat< T > const &M, arma::mat const &A, arma::vec const &v, T const &sum, DetailedPlane< T > &p)
 Extract detailed plane components. It must be invoked after common plane components are extracted into given plane. More...
 

Detailed Description

Handle plane fitting operations.

Author
Alberto M. Esmoris Pena
Version
1.0

Member Function Documentation

◆ bestFittingDetailedPlanePCA()

template<typename T >
static DetailedPlane<T> SurfaceInspector::maths::PlaneFitter::bestFittingDetailedPlanePCA ( Mat< T > &  M,
vector< T >  center 
)
static

Like bestFittingPlanePCA but returning a DetailedPlane which provides more information than a simple Plane object.

Parameters
centerWhen center is a non empty vector, it will be used as the origin of the vector subspace. In consequence, the matrix M will be treated as an affine space such that when transposed by center it will be a vector subspace
See also
PlaneFitter::bestFittingPlanePCA
SurfaceInspector::maths::DetailedPlane

◆ bestFittingDetailedPlaneSVD()

template<typename T >
static DetailedPlane<T> SurfaceInspector::maths::PlaneFitter::bestFittingDetailedPlaneSVD ( Mat< T > &  M,
vector< T >  center 
)
static

Like bestFittingPlaneSVD but returning a DetailedPlane which provides more information than a simple Plane object.

Parameters
centerWhen center is a non empty vector, it will be used as the origin of the vector subspace. In consequence, the matrix M will be treated as an affine space such that when transposed by center it will be a vector subspace
See also
PlaneFitter::bestFittingPlaneSVD
SurfaceInspector::maths::DetailedPlane

◆ bestFittingPlaneFromCovariances()

template<typename T >
static Plane<T> SurfaceInspector::maths::PlaneFitter::bestFittingPlaneFromCovariances ( Mat< T > &  M)
static

Compute the best fitting plane for given Matrix of covariances through eigen value decomposition.

Template Parameters
TType of number
Parameters
MMatrix of covariances
Returns
Best fitting plane

◆ bestFittingPlanePCA()

template<typename T >
static Plane<T> SurfaceInspector::maths::PlaneFitter::bestFittingPlanePCA ( Mat< T > &  M)
static

Compute the best fitting plane for given Matrix of coordinates through principal component analysis.

Template Parameters
TType of number
Parameters
MMatrix of coordinates. Row-coordinates format is mandatory:
M[0] := First X coordinate
M[1] := Second X coordinate
M[n-1] := Last X coordinate
M[n] := First Y coordinate
M[n+1] := Second Y coordinate
M[2n-1] := Last Y coordinate
M[2n] := First Z coordinate
M[2n+1] := Second Z coordinate
M[3n-1] := Last Z coordinate
Returns
Best fitting plane

◆ bestFittingPlaneSVD()

template<typename T >
static Plane<T> SurfaceInspector::maths::PlaneFitter::bestFittingPlaneSVD ( Mat< T > &  M)
static

Compute the best fitting plane for given Matrix of coordinates through singular value decomposition.

\[ M_{n_{x}3} = U_{n_{x}r} \Sigma_{r_{x}r} V^{T}_{r_{x}3} \]

Then \((V^{T}_{1,3} , V^{T}_{2,3} , V^{T}_{3,3})\) it is the orthonormal of best fitting plane.

Template Parameters
TType of number
Parameters
MMatrix of coordinates. Row-coordinates format is mandatory:
M[0] := First X coordinate
M[1] := Second X coordinate
M[n-1] := Last X coordinate
M[n] := First Y coordinate
M[n+1] := Second Y coordinate
M[2n-1] := Last Y coordinate
M[2n] := First Z coordinate
M[2n+1] := Second Z coordinate
M[3n-1] := Last Z coordinate
Returns
Best fitting plane

◆ centerCoordinatesMatrix()

template<typename T >
static vector<T> SurfaceInspector::maths::PlaneFitter::centerCoordinatesMatrix ( Mat< T > &  M)
static

Modify coordinates at M so they are centered at origin.

Template Parameters
TType of number
Parameters
MMatrix of coordinates. Following format is mandatory:
M[0] := First X coordinate
M[1] := Second X coordinate
M[n-1] := Last X coordinate
M[n] := First Y coordinate
M[n+1] := Second Y coordinate
M[2n-1] := Last Y coordinate
M[2n] := First Z coordinate
M[2n+1] := Second Z coordinate
M[3n-1] := Last Z coordinate
Returns
Original centroid

◆ extractCommonPlaneComponents()

template<typename T >
static void SurfaceInspector::maths::PlaneFitter::extractCommonPlaneComponents ( arma::mat const &  A,
arma::vec const &  v,
T &  sum,
Plane< T > &  p 
)
staticprotected

Extract common plane components.

Parameters
[in]AMatrix containing eigen vectors or singular vectors
[in]vVector of singular values or eigen values
[out]sumWhere summation of singular or eigen values will be stored. It must be passed with 0 as initial value.
[out]pPlane where extracted components will be stored

◆ extractDetailedPlaneComponents()

template<typename T >
static void SurfaceInspector::maths::PlaneFitter::extractDetailedPlaneComponents ( Mat< T > const &  M,
arma::mat const &  A,
arma::vec const &  v,
T const &  sum,
DetailedPlane< T > &  p 
)
staticprotected

Extract detailed plane components. It must be invoked after common plane components are extracted into given plane.

Parameters
[in]MMatrix of coordinates
[in]AMatrix containing eigen vectors or singular vectors
[in]vVector of singular values or eigen values
[in]sumIst must contain summation of singular or eigen values
pPlane where common components are stored and where detailed components will be placed. It is both an input/output argument
See also
PlaneFitter::extractCommonPlaneComponents

◆ translateToOrigin()

template<typename T >
static vector<T> SurfaceInspector::maths::PlaneFitter::translateToOrigin ( Mat< T > &  M,
vector< T >  center 
)
static

Translate M to origin transposing by center point.

Template Parameters
TType of number
Parameters
MMatrix of coordinates to be translated. It must have the same format than the one specified for centerCoordinatesMatrix function
centerReference point to compute the origin of vector subspace
Returns
Center as centroid by convenience
See also
PlaneFitter::centerCoordinatesMatrix(Mat<T> &M)

The documentation for this class was generated from the following file: