Ponca  e26a0e88a45818354616c1a7180bcd203aecad3c
Point Cloud Analysis library
Loading...
Searching...
No Matches
covariancePlaneFit.h
1/*
2 Copyright (C) 2014 Nicolas Mellado <nmellado0@gmail.com>
3 Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
4
5 This Source Code Form is subject to the terms of the Mozilla Public
6 License, v. 2.0. If a copy of the MPL was not distributed with this
7 file, You can obtain one at http://mozilla.org/MPL/2.0/.
8*/
9
10#pragma once
11
12#include "./plane.h"
13#include "./primitive.h"
14#include "./mean.h" // used to define CovarianceLineFit
15#include "./covarianceFit.h" // use to define CovariancePlaneFit
16
17#include <Eigen/Eigenvalues>
18
19namespace Ponca
20{
21
32template < class DataPoint, class _WFunctor, typename T >
33class CovariancePlaneFitImpl : public T
34{
35PONCA_FITTING_DECLARE_DEFAULT_TYPES
36PONCA_FITTING_DECLARE_MATRIX_TYPE
37
38protected:
39 enum
40 {
41 Check = Base::PROVIDES_PLANE &&
42 Base::PROVIDES_POSITION_COVARIANCE,
50 };
51
52public:
53 PONCA_EXPLICIT_CAST_OPERATORS(CovariancePlaneFitImpl,covariancePlaneFit)
54 PONCA_FITTING_DECLARE_FINALIZE
55 PONCA_FITTING_IS_SIGNED(false)
56
57 /**************************************************************************/
58 /* Results */
59 /**************************************************************************/
60
61
70 template <bool ignoreTranslation = false>
71 PONCA_MULTIARCH inline VectorType worldToTangentPlane(const VectorType &_q) const;
72
80 template <bool ignoreTranslation = false>
81 PONCA_MULTIARCH inline VectorType tangentPlaneToWorld(const VectorType &_q) const;
82}; //class CovariancePlaneFitImpl
83
86 template < class DataPoint, class _WFunctor, typename T>
88 CovariancePlaneFitImpl<DataPoint, _WFunctor,
89 CovarianceFitBase<DataPoint, _WFunctor,
90 MeanPosition<DataPoint, _WFunctor,
91 Plane<DataPoint, _WFunctor,T>>>>;
93
100template < class DataPoint, class _WFunctor, int DiffType, typename T>
102{
103 PONCA_FITTING_DECLARE_DEFAULT_TYPES
104 PONCA_FITTING_DECLARE_MATRIX_TYPE
105 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
106 static_assert ( DataPoint::Dim == 3, "CovariancePlaneDer is only valid in 3D");
107
108protected:
109 enum
110 {
111 Check = Base::PROVIDES_PLANE &
112 Base::PROVIDES_POSITION_COVARIANCE_DERIVATIVE,
114 PROVIDES_NORMAL_DERIVATIVE
115 };
116
117private:
118 VectorArray m_dNormal {VectorArray::Zero()};
119 ScalarArray m_dDist {ScalarArray::Zero()};
121public:
122 PONCA_EXPLICIT_CAST_OPERATORS_DER(CovariancePlaneDerImpl,covariancePlaneDer)
123
124
125 PONCA_MULTIARCH FIT_RESULT finalize();
126
127 /**************************************************************************/
128 /* Use results */
129 /**************************************************************************/
130
131 // \brief Returns the derivatives of the scalar field at the evaluation point
133 PONCA_MULTIARCH inline ScalarArray dPotential() const { return m_dDist; }
134
136 PONCA_MULTIARCH inline VectorArray dNormal() const { return m_dNormal; }
137
138}; //class CovariancePlaneDer
139
140
142template < class DataPoint, class _WFunctor, int DiffType, typename T>
144CovariancePlaneDerImpl<DataPoint, _WFunctor, DiffType,
145 CovarianceFitDer<DataPoint, _WFunctor, DiffType,
147
148#include "covariancePlaneFit.hpp"
149
150} //namespace Ponca
Line fitting procedure that minimize the orthogonal distance between the samples and the fitted primi...
Definition: covarianceFit.h:33
Internal generic class computing the derivatives of covariance matrix computed by CovarianceFitBase.
Definition: covarianceFit.h:75
[CovariancePlaneFit Definition]
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
@ PROVIDES_COVARIANCE_PLANE_DERIVATIVE
Provides derivatives for hyper-planes.
VectorArray dNormal() const
Returns the derivatives of the primitive normal.
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
Plane fitting procedure using only points position.
FIT_RESULT finalize()
Finalize the procedure.
CovariancePlaneFitImpl< DataPoint, _WFunctor, T > & covariancePlaneFit()
Explicit conversion to CovariancePlaneFitImpl , to access methods potentially hidden by heritage.
@ PROVIDES_TANGENT_PLANE_BASIS
Expose a method worldToTangentPlane(VectorType), which turns a point in ambient 3D space to the tange...
VectorType worldToTangentPlane(const VectorType &_q) const
Express a point in ambient space relatively to the tangent plane.
VectorType tangentPlaneToWorld(const VectorType &_q) const
Transform a point from the tangent plane [h, u, v]^T to ambient space.
typename Base::VectorType VectorType
Alias to vector type.
Compute the barycenter of the input points.
Definition: mean.h:24
Implicit hyperplane defined by an homogeneous vector .
Definition: plane.h:37
This Source Code Form is subject to the terms of the Mozilla Public License, v.
DiffType
Flags defining which derivatives need to be computed.
Definition: enums.h:34
FIT_RESULT
Enum corresponding to the state of a fitting method (and what the finalize function returns)
Definition: enums.h:15