Ponca  e26a0e88a45818354616c1a7180bcd203aecad3c
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
11#include "./algebraicSphere.h"
12#include "./mean.h" // used to define UnorientedSphereFit
13
14#include <Eigen/Dense>
15
16namespace Ponca
17{
44template < class DataPoint, class _WFunctor, typename T >
46{
47PONCA_FITTING_DECLARE_DEFAULT_TYPES
48
49protected:
50 enum { Check = Base::PROVIDES_ALGEBRAIC_SPHERE && Base::PROVIDES_MEAN_POSITION };
51
52 typedef Eigen::Matrix<Scalar, DataPoint::Dim+1, 1> VectorB;
53 typedef Eigen::Matrix<Scalar, DataPoint::Dim+1, DataPoint::Dim+1> MatrixBB;
54
55public:
56 using Solver = Eigen::EigenSolver<MatrixBB>;
57
58 MatrixBB m_matA {MatrixBB::Zero()};
59 MatrixBB m_matQ {MatrixBB::Zero()};
62 Solver m_solver;
63
64public:
65 PONCA_EXPLICIT_CAST_OPERATORS(UnorientedSphereFitImpl,unorientedSphereFit)
66 PONCA_FITTING_DECLARE_INIT_ADD_FINALIZE
67 PONCA_FITTING_IS_SIGNED(false)
68
69}; // class UnorientedSphereFitImpl
70
72template < class DataPoint, class _WFunctor, typename T>
74UnorientedSphereFitImpl<DataPoint, _WFunctor,
75 MeanPosition<DataPoint, _WFunctor,
77
78
79
80template < class DataPoint, class _WFunctor, int DiffType, typename T>
82{
83protected:
84 PONCA_FITTING_DECLARE_DEFAULT_TYPES
85 PONCA_FITTING_DECLARE_DEFAULT_DER_TYPES
86
87 using VectorB = typename Base::VectorB;
88 using MatrixBB = typename Base::MatrixBB;
89
90protected:
91 enum
92 {
93 Check = Base::PROVIDES_ALGEBRAIC_SPHERE &
94 Base::PROVIDES_MEAN_POSITION_DERIVATIVE &
95 Base::PROVIDES_PRIMITIVE_DERIVATIVE,
96 PROVIDES_ALGEBRAIC_SPHERE_DERIVATIVE,
97 PROVIDES_NORMAL_DERIVATIVE
98 };
99
100protected:
101 // computation data
102 MatrixBB m_dmatA[Base::NbDerivatives];
103 ScalarArray m_dSumDotPP;
104
105public:
106 // results
107 ScalarArray m_dUc;
108 VectorArray m_dUl;
109 ScalarArray m_dUq;
110
111public:
112 PONCA_EXPLICIT_CAST_OPERATORS_DER(UnorientedSphereDerImpl,unorientedSphereDer)
113 PONCA_FITTING_DECLARE_INIT_ADDDER_FINALIZE
114
115 PONCA_MULTIARCH inline ScalarArray dPotential() const;
116 PONCA_MULTIARCH inline VectorArray dNormal() const;
117
118}; //class UnorientedSphereDerImpl
119
120
121template < class DataPoint, class _WFunctor, int DiffType, typename T>
123 UnorientedSphereDerImpl<DataPoint, _WFunctor, DiffType,
125
126
127} //namespace Ponca
128
129#include "unorientedSphereFit.hpp"
Algebraic Sphere primitive.
Compute the barycenter of the input points.
Definition: mean.h:24
UnorientedSphereDerImpl< DataPoint, _WFunctor, DiffType, T > & unorientedSphereDer()
Explicit conversion to UnorientedSphereDerImpl , to access methods potentially hidden by heritage.
typename Base::ScalarArray ScalarArray
Alias to scalar derivatives array.
typename Base::VectorArray VectorArray
Alias to vector derivatives array.
Algebraic Sphere fitting procedure on point sets with non-oriented normals.
MatrixBB m_matA
The accumulated covariance matrix.
Scalar m_sumDotPP
Sum of the squared relative positions.
typename DataPoint::Scalar Scalar
Alias to scalar type.
UnorientedSphereFitImpl< DataPoint, _WFunctor, T > & unorientedSphereFit()
Explicit conversion to UnorientedSphereFitImpl , to access methods potentially hidden by heritage.
MatrixBB m_matQ
The constraint matrix.
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