Ponca  7d8ac87a7de01d881c9fde3c42e397b44bffb901
Point Cloud Analysis library
Loading...
Searching...
No Matches
orientedSphereFit.h
1/*
2This Source Code Form is subject to the terms of the Mozilla Public
3License, v. 2.0. If a copy of the MPL was not distributed with this
4file, You can obtain one at http://mozilla.org/MPL/2.0/.
5*/
6
7#pragma once
8
9#include "./algebraicSphere.h"
10#include "./mean.h" // used to define OrientedSphereFit
11
12namespace Ponca
13{
14
24 template <class DataPoint, class _NFilter, typename T>
25 class OrientedSphereFitImpl : public T
26 {
27 PONCA_FITTING_DECLARE_DEFAULT_TYPES
28
29 protected:
30 enum
31 {
32 Check = Base::PROVIDES_ALGEBRAIC_SPHERE && Base::PROVIDES_MEAN_NORMAL && Base::PROVIDES_MEAN_POSITION
33 };
34
35 // computation data
38 m_nume{0},
39 m_deno{0};
41 public:
42 PONCA_EXPLICIT_CAST_OPERATORS(OrientedSphereFitImpl, orientedSphereFit)
43 PONCA_FITTING_DECLARE_INIT_ADD_FINALIZE
44 PONCA_FITTING_IS_SIGNED(true)
45 }; // class OrientedSphereFitImpl
46
49 template <class DataPoint, class _NFilter, typename T>
51 DataPoint, _NFilter,
54
58 template <class DataPoint, class _NFilter, int DiffType, typename T>
59 class OrientedSphereDerImpl : public T
60 {
61 PONCA_FITTING_DECLARE_DEFAULT_TYPES
62 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
63
64 protected:
65 enum
66 {
67 Check = Base::PROVIDES_ALGEBRAIC_SPHERE & Base::PROVIDES_MEAN_POSITION_DERIVATIVE &
68 Base::PROVIDES_PRIMITIVE_DERIVATIVE,
69 PROVIDES_ALGEBRAIC_SPHERE_DERIVATIVE,
70 PROVIDES_NORMAL_DERIVATIVE
71 };
72
73 // computation data
74 VectorArray m_dSumN{VectorArray::Zero()};
75 ScalarArray m_dSumDotPN{ScalarArray::Zero()},
78 ScalarArray::Zero()},
79 m_dNume{ScalarArray::Zero()},
80 m_dDeno{ScalarArray::Zero()};
82 public:
83 // results
84 ScalarArray m_dUc{ScalarArray::Zero()},
85 m_dUq{ScalarArray::Zero()};
86 VectorArray m_dUl{VectorArray::Zero()};
88 public:
89 PONCA_EXPLICIT_CAST_OPERATORS_DER(OrientedSphereDerImpl, orientedSphereDer)
90
91 PONCA_FITTING_DECLARE_INIT_ADDDER_FINALIZE
92
94 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dPotential() const;
95
97 PONCA_MULTIARCH [[nodiscard]] inline VectorArray dNormal() const;
98
100 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dprattNorm2() const
101 {
102 return Scalar(2.) * Base::m_ul.transpose() * m_dUl - Scalar(4.) * Base::m_uq * m_dUc -
103 Scalar(4.) * Base::m_uc * m_dUq;
104 }
105
107 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm2(unsigned int _d) const
108 {
109 return Scalar(2.) * m_dUl.col(_d).dot(Base::m_ul) - Scalar(4.) * m_dUc.col(_d)[0] * Base::m_uq -
110 Scalar(4.) * m_dUq.col(_d)[0] * Base::m_uc;
111 }
112
114 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm(unsigned int _d) const
115 {
116 PONCA_MULTIARCH_STD_MATH(sqrt);
117 return sqrt(dprattNorm2(_d));
118 }
119
121 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm() const
122 {
123 PONCA_MULTIARCH_STD_MATH(sqrt);
124 return dprattNorm2().array().sqrt();
125 }
127
131 PONCA_MULTIARCH [[nodiscard]] inline bool applyPrattNorm();
132
133 }; // class OrientedSphereDerImpl
134
136 template <class DataPoint, class _NFilter, int DiffType, typename T>
139
140#include "orientedSphereFit.hpp"
141
142} // namespace Ponca
Aggregator class used to declare specialized structures using CRTP.
Definition basket.h:318
[OrientedSphereFit Definition]
ScalarArray m_dSumDotPP
Sum of the squared relative positions with weight derivatives.
ScalarArray m_dUc
Derivatives of the hyper-sphere constant term
ScalarArray m_dSumDotPN
Sum of the dot product between relative positions and normals with weight derivatives.
VectorArray m_dSumN
Sum of the normal vectors with weight derivatives.
ScalarArray m_dNume
Derivatives of the numerator of the quadratic parameter
VectorArray m_dUl
Derivatives of the hyper-sphere linear term
Scalar dprattNorm2(unsigned int _d) const
compute the square of the Pratt norm derivative for dimension _d
ScalarArray m_dUq
Derivatives of the hyper-sphere quadratic term.
VectorArray dNormal() const
Returns the derivatives of the primitive normal.
ScalarArray dprattNorm2() const
compute the square of the Pratt norm derivative
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
ScalarArray dPotential() const
Returns the derivatives of the scalar field at the evaluation point.
typename DataPoint::Scalar Scalar
Alias to scalar type.
bool applyPrattNorm()
Normalize the scalar field by the Pratt norm.
Scalar dprattNorm() const
compute the Pratt norm derivative
OrientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & orientedSphereDer()
Explicit conversion to OrientedSphereDerImpl , to access methods potentially hidden by heritage.
ScalarArray m_dDeno
Derivatives of the denominator of the quadratic parameter.
Scalar dprattNorm(unsigned int _d) const
compute the Pratt norm derivative for the dimension _d
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
Algebraic Sphere fitting procedure on oriented point sets.
Scalar m_sumDotPP
Sum of the squared relative positions.
OrientedSphereFitImpl< DataPoint, _NFilter, T > & orientedSphereFit()
Explicit conversion to OrientedSphereFitImpl , to access methods potentially hidden by heritage.
typename DataPoint::Scalar Scalar
Alias to scalar type.
Scalar m_sumDotPN
Sum of the dot product between relative positions and normals.
Scalar m_nume
Numerator of the quadratic parameter (excluding the 0.5 coefficient)
Scalar m_deno
Denominator of the quadratic parameter (excluding the 0.5 coefficient)
This Source Code Form is subject to the terms of the Mozilla Public License, v.