Ponca  7b3f8ad3fde25a027e6452783ccee143798a71b8
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#include "../../defines.h"
9#include "../concepts.h"
10
11#include "../primitive/algebraicSphere.h" // used to define OrientedSphereFit
12#include "../tools/mean.h" // used to define OrientedSphereFit
13
14#define ORIENTED_SPHERE_FIT_REQUIREMENTS ProvidesAlgebraicSphere<T>&& ProvidesMeanPosition<T>&& ProvidesMeanNormal<T>
15#define ORIENTED_SPHERE_DER_REQUIREMENTS \
16 ProvidesBasketDiffUnitBase<T>&& ProvidesAlgebraicSphere<T>&& ProvidesMeanPositionDerivative<T>
17
18namespace Ponca
19{
20
30 template <class DataPoint, class _NFilter, typename T>
31 requires ORIENTED_SPHERE_FIT_REQUIREMENTS
32 class OrientedSphereFitImpl : public T
33 {
34 PONCA_FITTING_DECLARE_DEFAULT_TYPES
35
36 protected:
37 // computation data
40 m_nume{0},
41 m_deno{0};
43 public:
44 PONCA_EXPLICIT_CAST_OPERATORS(OrientedSphereFitImpl, orientedSphereFit)
45 PONCA_FITTING_DECLARE_INIT_ADD_FINALIZE
46 PONCA_FITTING_IS_SIGNED(true)
47 }; // class OrientedSphereFitImpl
48
51 template <class DataPoint, class _NFilter, typename T>
53 DataPoint, _NFilter,
56
60 template <class DataPoint, class _NFilter, int DiffType, typename T>
61 requires ORIENTED_SPHERE_DER_REQUIREMENTS
62 class OrientedSphereDerImpl : public T
63 {
64 PONCA_FITTING_DECLARE_DEFAULT_TYPES
65 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
66
67 protected:
68 // computation data
69 VectorArray m_dSumN{VectorArray::Zero()};
70 ScalarArray m_dSumDotPN{ScalarArray::Zero()},
73 ScalarArray::Zero()},
74 m_dNume{ScalarArray::Zero()},
75 m_dDeno{ScalarArray::Zero()};
77 public:
78 // results
79 ScalarArray m_dUc{ScalarArray::Zero()},
80 m_dUq{ScalarArray::Zero()};
81 VectorArray m_dUl{VectorArray::Zero()};
83 public:
84 PONCA_EXPLICIT_CAST_OPERATORS_DER(OrientedSphereDerImpl, orientedSphereDer)
85 PONCA_EXPLICIT_CAST_OPERATORS_DER(OrientedSphereDerImpl, algebraicSphereDer)
86 PONCA_EXPLICIT_CAST_OPERATORS_DER(OrientedSphereDerImpl, implicitPrimitiveDer)
87
88 PONCA_FITTING_DECLARE_INIT_ADDDER_FINALIZE
89
91 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dPotential() const;
92
94 PONCA_MULTIARCH [[nodiscard]] inline VectorArray dNormal() const;
95
97 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dprattNorm2() const
98 {
99 return Scalar(2.) * Base::m_ul.transpose() * m_dUl - Scalar(4.) * Base::m_uq * m_dUc -
100 Scalar(4.) * Base::m_uc * m_dUq;
101 }
102
104 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm2(unsigned int _d) const
105 {
106 return Scalar(2.) * m_dUl.col(_d).dot(Base::m_ul) - Scalar(4.) * m_dUc.col(_d)[0] * Base::m_uq -
107 Scalar(4.) * m_dUq.col(_d)[0] * Base::m_uc;
108 }
109
111 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm(unsigned int _d) const
112 {
113 PONCA_MULTIARCH_STD_MATH(sqrt);
114 return sqrt(dprattNorm2(_d));
115 }
116
118 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dprattNorm() const
119 {
120 PONCA_MULTIARCH_STD_MATH(sqrt);
121 return dprattNorm2().array().sqrt();
122 }
124
128 PONCA_MULTIARCH [[nodiscard]] inline bool applyPrattNorm();
129
130 }; // class OrientedSphereDerImpl
131
133 template <class DataPoint, class _NFilter, int DiffType, typename T>
136
137#include "orientedSphereFit.hpp"
138
139} // namespace Ponca
Aggregator class used to declare specialized structures using CRTP.
Definition basket.h:294
[OrientedSphereFit Definition]
ScalarArray m_dSumDotPP
Sum of the squared relative positions with weight derivatives.
OrientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & implicitPrimitiveDer()
Explicit conversion to OrientedSphereDerImpl , to access methods potentially hidden by heritage.
ScalarArray dPotential() const
Returns the derivatives of the scalar field at the evaluation point.
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.
ScalarArray dprattNorm2() const
compute the square of the Pratt norm derivative
ScalarArray dprattNorm() const
compute the Pratt norm derivative
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
VectorArray dNormal() const
Returns the derivatives of the primitive normal.
typename DataPoint::Scalar Scalar
Alias to scalar type.
bool applyPrattNorm()
Normalize the scalar field by the Pratt norm.
OrientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & algebraicSphereDer()
Explicit conversion to OrientedSphereDerImpl , to access methods potentially hidden by heritage.
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.
Definition concepts.h:11