Ponca  94e4a36411c3364f6192f97adfa8ceee67834d6e
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 _NFilter, 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
72 PONCA_MULTIARCH inline VectorType worldToTangentPlane(const VectorType &_q,
73 bool _isPositionVector = true) const;
74
83 PONCA_MULTIARCH inline VectorType tangentPlaneToWorld(const VectorType &_q,
84 bool _isPositionVector = true) const;
85}; //class CovariancePlaneFitImpl
86
89 template < class DataPoint, class _NFilter, typename T>
91 CovariancePlaneFitImpl<DataPoint, _NFilter,
92 CovarianceFitBase<DataPoint, _NFilter,
93 MeanPosition<DataPoint, _NFilter,
94 Plane<DataPoint, _NFilter, T>>>>;
96
103template < class DataPoint, class _NFilter, int DiffType, typename T>
105{
106 PONCA_FITTING_DECLARE_DEFAULT_TYPES
107 PONCA_FITTING_DECLARE_MATRIX_TYPE
108 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
109 static_assert ( DataPoint::Dim == 3, "CovariancePlaneDer is only valid in 3D");
110
111protected:
112 enum
113 {
114 Check = Base::PROVIDES_PLANE &
115 Base::PROVIDES_POSITION_COVARIANCE_DERIVATIVE,
117 PROVIDES_NORMAL_DERIVATIVE
118 };
119
120private:
121 VectorArray m_dNormal {VectorArray::Zero()};
122 ScalarArray m_dDist {ScalarArray::Zero()};
124public:
125 PONCA_EXPLICIT_CAST_OPERATORS_DER(CovariancePlaneDerImpl,covariancePlaneDer)
126
127
128 PONCA_MULTIARCH FIT_RESULT finalize();
129
130 /**************************************************************************/
131 /* Use results */
132 /**************************************************************************/
133
134 // \brief Returns the derivatives of the scalar field at the evaluation point
136 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dPotential() const { return m_dDist; }
137
139 PONCA_MULTIARCH [[nodiscard]] inline VectorArray dNormal() const { return m_dNormal; }
140
141}; //class CovariancePlaneDer
142
143
145template < class DataPoint, class _NFilter, int DiffType, typename T>
147CovariancePlaneDerImpl<DataPoint, _NFilter, DiffType,
148 CovarianceFitDer<DataPoint, _NFilter, DiffType,
150
151#include "covariancePlaneFit.hpp"
152
153} //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.
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
@ PROVIDES_COVARIANCE_PLANE_DERIVATIVE
Provides derivatives for hyper-planes.
VectorArray dNormal() const
Returns the derivatives of the primitive normal.
Plane fitting procedure using only points position.
VectorType tangentPlaneToWorld(const VectorType &_q, bool _isPositionVector=true) const
Transform a point from the tangent plane [h, u, v]^T to ambient space.
FIT_RESULT finalize()
Finalize the procedure.
CovariancePlaneFitImpl< DataPoint, _NFilter, T > & covariancePlaneFit()
Explicit conversion to CovariancePlaneFitImpl , to access methods potentially hidden by heritage.
VectorType worldToTangentPlane(const VectorType &_q, bool _isPositionVector=true) const
Express a point in ambient space relatively to the tangent plane.
@ PROVIDES_TANGENT_PLANE_BASIS
Expose a method worldToTangentPlane(VectorType), which turns a point in ambient 3D space to the tange...
typename Base::VectorType VectorType
Alias to vector type.
Compute the derivatives of the input points barycenter.
Definition mean.h:114
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