Helios++
Helios software for LiDAR simulations
PlaneFitter.hpp
1 #ifndef _SURFACEINSPECTOR_MATHS_PLANEFITTER_HPP_
2 #define _SURFACEINSPECTOR_MATHS_PLANEFITTER_HPP_
3 
4 #include <vector>
5 
6 #include <armadillo>
7 
8 #include <surfaceinspector/util/Object.hpp>
9 #include <surfaceinspector/maths/Plane.hpp>
10 #include <surfaceinspector/maths/DetailedPlane.hpp>
11 
12 using std::vector;
13 
14 using arma::Mat;
15 
19 
20 namespace SurfaceInspector { namespace maths{
21 
28 class PlaneFitter : public Object{
29 public:
30  // *** STATIC FUNCTIONS *** //
31  // ************************** //
47  template <typename T>
48  static vector<T> centerCoordinatesMatrix(Mat<T> & M);
49 
59  template <typename T>
60  static vector<T> translateToOrigin(Mat<T> &M, vector<T> center);
61 
87  template <typename T>
88  static Plane<T> bestFittingPlaneSVD(Mat<T> & M);
106  template <typename T>
107  static Plane<T> bestFittingPlanePCA(Mat<T> &M);
108 
116  template <typename T>
118 
129  template <typename T>
131  Mat<T> &M,
132  vector<T> center
133  );
144  template <typename T>
146  Mat<T> &M,
147  vector<T> center
148  );
149 
150 protected:
151  // *** INNER FUNCTIONS *** //
152  // ************************* //
161  template <typename T>
163  arma::mat const &A,
164  arma::vec const &v,
165  T &sum,
166  Plane<T> &p
167  );
168 
180  template <typename T>
182  Mat<T> const &M,
183  arma::mat const &A,
184  arma::vec const &v,
185  T const &sum,
187  );
188 
189 };
190 
191 }}
192 
193 #include <surfaceinspector/maths/PlaneFitter.tpp>
194 
195 #endif
Definition: DetailedPlane.hpp:8
Handle plane fitting operations.
Definition: PlaneFitter.hpp:28
static DetailedPlane< T > bestFittingDetailedPlanePCA(Mat< T > &M, vector< T > center)
Like bestFittingPlanePCA but returning a DetailedPlane which provides more information than a simple ...
static vector< T > translateToOrigin(Mat< T > &M, vector< T > center)
Translate M to origin transposing by center point.
static void extractDetailedPlaneComponents(Mat< T > const &M, arma::mat const &A, arma::vec const &v, T const &sum, DetailedPlane< T > &p)
Extract detailed plane components. It must be invoked after common plane components are extracted int...
static DetailedPlane< T > bestFittingDetailedPlaneSVD(Mat< T > &M, vector< T > center)
Like bestFittingPlaneSVD but returning a DetailedPlane which provides more information than a simple ...
static vector< T > centerCoordinatesMatrix(Mat< T > &M)
Modify coordinates at M so they are centered at origin.
static Plane< T > bestFittingPlaneSVD(Mat< T > &M)
Compute the best fitting plane for given Matrix of coordinates through singular value decomposition.
static Plane< T > bestFittingPlaneFromCovariances(Mat< T > &M)
Compute the best fitting plane for given Matrix of covariances through eigen value decomposition.
static Plane< T > bestFittingPlanePCA(Mat< T > &M)
Compute the best fitting plane for given Matrix of coordinates through principal component analysis.
static void extractCommonPlaneComponents(arma::mat const &A, arma::vec const &v, T &sum, Plane< T > &p)
Extract common plane components.
Class representing a plane.
Definition: Plane.hpp:22
Class representing an object. All surface inspector classes must extend Object.
Definition: Object.hpp:12