Ponca  40f245e28b920cbb763a1c6282156c87c626f24c
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 "./defines.h"
13#include "./plane.h"
14#include "./primitive.h"
15#include "./mean.h" // used to define CovarianceLineFit
16#include "./covarianceFit.h" // use to define CovariancePlaneFit
17
18#include <Eigen/Eigenvalues>
19
20namespace Ponca
21{
22
33template < class DataPoint, class _WFunctor, typename T >
34class CovariancePlaneFitImpl : public T
35{
36PONCA_FITTING_DECLARE_DEFAULT_TYPES
37PONCA_FITTING_DECLARE_MATRIX_TYPE
38
39protected:
40 enum
41 {
42 Check = Base::PROVIDES_PLANE &&
43 Base::PROVIDES_POSITION_COVARIANCE,
51 };
52
53public:
54 PONCA_EXPLICIT_CAST_OPERATORS(CovariancePlaneFitImpl,covariancePlaneFit)
55 PONCA_FITTING_DECLARE_FINALIZE
56 PONCA_FITTING_IS_SIGNED(false)
57
58 /**************************************************************************/
59 /* Results */
60 /**************************************************************************/
61
62
71 template <bool ignoreTranslation = false>
72 PONCA_MULTIARCH inline VectorType worldToTangentPlane(const VectorType &_q) const;
73
81 template <bool ignoreTranslation = false>
82 PONCA_MULTIARCH inline VectorType tangentPlaneToWorld(const VectorType &_q) const;
83}; //class CovariancePlaneFitImpl
84
87 template < class DataPoint, class _WFunctor, typename T>
89 CovariancePlaneFitImpl<DataPoint, _WFunctor,
90 CovarianceFitBase<DataPoint, _WFunctor,
91 MeanPosition<DataPoint, _WFunctor,
92 Plane<DataPoint, _WFunctor,T>>>>;
94
101template < class DataPoint, class _WFunctor, int DiffType, typename T>
103{
104 PONCA_FITTING_DECLARE_DEFAULT_TYPES
105 PONCA_FITTING_DECLARE_MATRIX_TYPE
106 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
107 static_assert ( DataPoint::Dim == 3, "CovariancePlaneDer is only valid in 3D");
108
109protected:
110 enum
111 {
112 Check = Base::PROVIDES_PLANE &
113 Base::PROVIDES_POSITION_COVARIANCE_DERIVATIVE,
115 PROVIDES_NORMAL_DERIVATIVE
116 };
117
118private:
119 VectorArray m_dNormal {VectorArray::Zero()};
120 ScalarArray m_dDist {ScalarArray::Zero()};
122public:
123 PONCA_EXPLICIT_CAST_OPERATORS_DER(CovariancePlaneDerImpl,covariancePlaneDer)
124
125
126 PONCA_MULTIARCH FIT_RESULT finalize();
127
128 /**************************************************************************/
129 /* Use results */
130 /**************************************************************************/
131
132 // \brief Returns the derivatives of the scalar field at the evaluation point
134 PONCA_MULTIARCH inline ScalarArray dPotential() const { return m_dDist; }
135
137 PONCA_MULTIARCH inline VectorArray dNormal() const { return m_dNormal; }
138
139}; //class CovariancePlaneDer
140
141
143template < class DataPoint, class _WFunctor, int DiffType, typename T>
145CovariancePlaneDerImpl<DataPoint, _WFunctor, DiffType,
146 CovarianceFitDer<DataPoint, _WFunctor, DiffType,
148
149#include "covariancePlaneFit.hpp"
150
151} //namespace Ponca
Procedure that compute and decompose the covariance matrix of the neighbors positions in .
Internal generic class computing the derivatives of covariance matrix computed by CovarianceFitBase.
[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