Ponca  7b3f8ad3fde25a027e6452783ccee143798a71b8
Point Cloud Analysis library
Loading...
Searching...
No Matches
unorientedSphereFit.h
1/*
2 Copyright (C) 2013 Gael Guennebaud <gael.guennebaud@inria.fr>
3
4 This Source Code Form is subject to the terms of the Mozilla Public
5 License, v. 2.0. If a copy of the MPL was not distributed with this
6 file, You can obtain one at http://mozilla.org/MPL/2.0/.
7*/
8
9#pragma once
10#include "../../defines.h"
11#include "../concepts.h"
12
13#include "../primitive/algebraicSphere.h" // used to define UnorientedSphereFit
14#include "../tools/mean.h" // used to define UnorientedSphereFit
15
16#include <Eigen/Dense>
17
18#define UNORIENTED_SPHERE_FIT_REQUIREMENTS ProvidesAlgebraicSphere<T>&& ProvidesMeanPosition<T>
19#define UNORIENTED_SPHERE_DER_REQUIREMENTS \
20 ProvidesBasketDiffUnitBase<T>&& ProvidesAlgebraicSphere<T>&& ProvidesMeanPositionDerivative<T>
21
22namespace Ponca
23{
53 template <class DataPoint, class _NFilter, typename T>
54 requires UNORIENTED_SPHERE_FIT_REQUIREMENTS
55 class UnorientedSphereFitImpl : public T
56 {
57 PONCA_FITTING_DECLARE_DEFAULT_TYPES
58
59 protected:
60 using VectorB = Eigen::Matrix<Scalar, DataPoint::Dim + 1, 1>;
61 using MatrixBB = Eigen::Matrix<Scalar, DataPoint::Dim + 1, DataPoint::Dim + 1>;
62
63 public:
64 using Solver = Eigen::EigenSolver<MatrixBB>;
65
66 MatrixBB m_matA{MatrixBB::Zero()};
67 MatrixBB m_matQ{MatrixBB::Zero()};
70 Solver m_solver;
71
72 public:
73 PONCA_EXPLICIT_CAST_OPERATORS(UnorientedSphereFitImpl, unorientedSphereFit)
74 PONCA_FITTING_DECLARE_INIT_ADD_FINALIZE
75 PONCA_FITTING_IS_SIGNED(false)
76
77 }; // class UnorientedSphereFitImpl
78
80 template <class DataPoint, class _NFilter, typename T>
84
85 template <class DataPoint, class _NFilter, int DiffType, typename T>
87 class UnorientedSphereDerImpl : public T
88 {
89 protected:
90 PONCA_FITTING_DECLARE_DEFAULT_TYPES
91 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
92
93 using VectorB = typename Base::VectorB;
94 using MatrixBB = typename Base::MatrixBB;
95
96 protected:
97 // computation data
98 MatrixBB m_dmatA[Base::NbDerivatives];
99 ScalarArray m_dSumDotPP{ScalarArray::Zero()};
100
101 public:
102 // results
103 ScalarArray m_dUc{ScalarArray::Zero()};
104 VectorArray m_dUl{VectorArray::Zero()};
105 ScalarArray m_dUq{ScalarArray::Zero()};
106
107 public:
108 PONCA_EXPLICIT_CAST_OPERATORS_DER(UnorientedSphereDerImpl, unorientedSphereDer)
109 PONCA_EXPLICIT_CAST_OPERATORS_DER(UnorientedSphereDerImpl, algebraicSphereDer)
110 PONCA_EXPLICIT_CAST_OPERATORS_DER(UnorientedSphereDerImpl, implicitPrimitiveDer)
111
112 PONCA_FITTING_DECLARE_INIT_ADDDER_FINALIZE
113
114 PONCA_MULTIARCH inline ScalarArray dPotential() const;
115 PONCA_MULTIARCH inline VectorArray dNormal() const;
116
118 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dprattNorm2() const
119 {
120 return Scalar(2.) * Base::m_ul.transpose() * m_dUl - Scalar(4.) * Base::m_uq * m_dUc -
121 Scalar(4.) * Base::m_uc * m_dUq;
122 }
123
125 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm2(unsigned int _d) const
126 {
127 return Scalar(2.) * m_dUl.col(_d).dot(Base::m_ul) - Scalar(4.) * m_dUc.col(_d)[0] * Base::m_uq -
128 Scalar(4.) * m_dUq.col(_d)[0] * Base::m_uc;
129 }
130
132 PONCA_MULTIARCH [[nodiscard]] inline Scalar dprattNorm(unsigned int _d) const
133 {
134 PONCA_MULTIARCH_STD_MATH(sqrt);
135 return sqrt(dprattNorm2(_d));
136 }
137
139 PONCA_MULTIARCH [[nodiscard]] inline ScalarArray dprattNorm() const
140 {
141 PONCA_MULTIARCH_STD_MATH(sqrt);
142 return dprattNorm2().array().sqrt();
143 }
145
149 PONCA_MULTIARCH [[nodiscard]] inline bool applyPrattNorm();
150
151 }; // class UnorientedSphereDerImpl
152
153 template <class DataPoint, class _NFilter, int DiffType, typename T>
154 using UnorientedSphereDer =
155 UnorientedSphereDerImpl<DataPoint, _NFilter, DiffType, MeanPositionDer<DataPoint, _NFilter, DiffType, T>>;
156
157} // namespace Ponca
158
159#include "unorientedSphereFit.hpp"
Aggregator class used to declare specialized structures using CRTP.
Definition basket.h:294
UnorientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & implicitPrimitiveDer()
Explicit conversion to UnorientedSphereDerImpl , to access methods potentially hidden by heritage.
ScalarArray dprattNorm() const
compute the Pratt norm derivative
Scalar dprattNorm2(unsigned int _d) const
compute the square of the Pratt norm derivative for dimension _d
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
UnorientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & algebraicSphereDer()
Explicit conversion to UnorientedSphereDerImpl , to access methods potentially hidden by heritage.
Scalar dprattNorm(unsigned int _d) const
compute the Pratt norm derivative for the dimension _d
ScalarArray dprattNorm2() const
compute the square of the Pratt norm derivative
typename DataPoint::Scalar Scalar
Alias to scalar type.
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
bool applyPrattNorm()
Normalize the scalar field by the Pratt norm.
UnorientedSphereDerImpl< DataPoint, _NFilter, DiffType, T > & unorientedSphereDer()
Explicit conversion to UnorientedSphereDerImpl , to access methods potentially hidden by heritage.
Algebraic Sphere fitting procedure on point sets with non-oriented normals.
Scalar m_sumDotPP
Sum of the squared relative positions.
UnorientedSphereFitImpl< DataPoint, _NFilter, T > & unorientedSphereFit()
Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.
MatrixBB m_matA
The accumulated covariance matrix.
MatrixBB m_matQ
The constraint matrix.
typename DataPoint::Scalar Scalar
Alias to scalar type.
This concept ensures that the default types and accessors in a BasketDiff are well-formed.
Definition concepts.h:40
This Source Code Form is subject to the terms of the Mozilla Public License, v.
Definition concepts.h:11