Ponca  e26a0e88a45818354616c1a7180bcd203aecad3c
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
24template < class DataPoint, class _WFunctor, typename T >
25class OrientedSphereFitImpl : public T
26{
27 PONCA_FITTING_DECLARE_DEFAULT_TYPES
28
29protected:
30 enum
31 {
32 Check = Base::PROVIDES_ALGEBRAIC_SPHERE &&
33 Base::PROVIDES_MEAN_NORMAL &&
34 Base::PROVIDES_MEAN_POSITION
35 };
36
37 // computation data
43public:
44 PONCA_EXPLICIT_CAST_OPERATORS(OrientedSphereFitImpl,orientedSphereFit)
45 PONCA_FITTING_DECLARE_INIT_ADD_FINALIZE
46 PONCA_FITTING_IS_SIGNED(true)
47}; //class OrientedSphereFitImpl
48
51template < class DataPoint, class _WFunctor, typename T>
53OrientedSphereFitImpl<DataPoint, _WFunctor,
54 MeanPosition<DataPoint, _WFunctor,
55 MeanNormal<DataPoint, _WFunctor,
58
62template < class DataPoint, class _WFunctor, int DiffType, typename T>
63class OrientedSphereDerImpl : public T
64{
65 PONCA_FITTING_DECLARE_DEFAULT_TYPES
66 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
67
68protected:
69 enum
70 {
71 Check = Base::PROVIDES_ALGEBRAIC_SPHERE &
72 Base::PROVIDES_MEAN_POSITION_DERIVATIVE &
73 Base::PROVIDES_PRIMITIVE_DERIVATIVE,
74 PROVIDES_ALGEBRAIC_SPHERE_DERIVATIVE,
75 PROVIDES_NORMAL_DERIVATIVE
76 };
77
78 // computation data
85public:
86 // results
91public:
92 PONCA_EXPLICIT_CAST_OPERATORS_DER(OrientedSphereDerImpl,orientedSphereDer)
93 PONCA_FITTING_DECLARE_INIT_ADDDER_FINALIZE
94
96 PONCA_MULTIARCH inline ScalarArray dPotential() const;
97
99 PONCA_MULTIARCH inline VectorArray dNormal() const;
100
102 PONCA_MULTIARCH inline ScalarArray dprattNorm2() const
103 {
104 return Scalar(2.) * Base::m_ul.transpose() * m_dUl
105 - Scalar(4.) * Base::m_uq * m_dUc
106 - Scalar(4.) * Base::m_uc * m_dUq;
107 }
108
110 PONCA_MULTIARCH inline Scalar dprattNorm2(unsigned int _d) const
111 {
112 return Scalar(2.) * m_dUl.col(_d).dot(Base::m_ul)
113 - Scalar(4.) * m_dUc.col(_d)[0]*Base::m_uq
114 - Scalar(4.) * m_dUq.col(_d)[0]*Base::m_uc;
115 }
116
118 PONCA_MULTIARCH inline Scalar dprattNorm(unsigned int _d) const
119 {
120 PONCA_MULTIARCH_STD_MATH(sqrt);
121 return sqrt(dprattNorm2(_d));
122 }
123
125 PONCA_MULTIARCH inline Scalar dprattNorm() const
126 {
127 PONCA_MULTIARCH_STD_MATH(sqrt);
128 return dprattNorm2().array().sqrt();
129 }
131
135 PONCA_MULTIARCH inline bool applyPrattNorm();
136
137}; //class OrientedSphereDerImpl
138
140 template < class DataPoint, class _WFunctor, int DiffType, typename T>
142 OrientedSphereDerImpl<DataPoint, _WFunctor, DiffType,
144
145#include "orientedSphereFit.hpp"
146
147} //namespace Ponca
Algebraic Sphere primitive.
Compute the barycenter of the input points + their normals.
Definition: mean.h:65
Compute the barycenter of the input points.
Definition: mean.h:24
[OrientedSphereFit Definition]
OrientedSphereDerImpl< DataPoint, _WFunctor, DiffType, T > & orientedSphereDer()
Explicit conversion to OrientedSphereDerImpl , to access methods potentially hidden by heritage.
VectorArray m_dSumN
Sum of the normal vectors with differenciated weights.
Scalar dprattNorm2(unsigned int _d) const
compute the square of the Pratt norm derivative for dimension _d
ScalarArray m_dUc
Derivative of the hyper-sphere constant term
Scalar dprattNorm(unsigned int _d) const
compute the Pratt norm derivative for the dimension _d
typename DataPoint::Scalar Scalar
Alias to scalar type.
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
ScalarArray m_dDeno
Differenciation of the denominator of the quadratic parameter.
ScalarArray m_dSumDotPN
Sum of the dot product betwen relative positions and normals with differenciated weights.
ScalarArray dprattNorm2() const
compute the square of the Pratt norm derivative
VectorArray dNormal() const
Returns the derivatives of the primitive normal.
VectorArray m_dUl
Derivative of the hyper-sphere linear term
Scalar dprattNorm() const
compute the Pratt norm derivative
ScalarArray dPotential() const
Returns the derivatives of the scalar field at the evaluation point.
ScalarArray m_dNume
Differenciation of the numerator of the quadratic parameter
bool applyPrattNorm()
Normalize the scalar field by the Pratt norm.
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
ScalarArray m_dSumDotPP
Sum of the squared relative positions with differenciated weights.
ScalarArray m_dUq
Derivative of the hyper-sphere quadratic term.
Algebraic Sphere fitting procedure on oriented point sets.
typename DataPoint::Scalar Scalar
Alias to scalar type.
OrientedSphereFitImpl< DataPoint, _WFunctor, T > & orientedSphereFit()
Explicit conversion to OrientedSphereFitImpl , to access methods potentially hidden by heritage.
Scalar m_deno
Denominator of the quadratic parameter (excluding the 0.5 coefficient)
Scalar m_sumDotPN
Sum of the dot product betwen relative positions and normals.
Scalar m_sumDotPP
Sum of the squared relative positions.
Scalar m_nume
Numerator of the quadratic parameter (excluding the 0.5 coefficient)
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