Point Cloud Library (PCL)  1.11.0
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
vector_average.hpp
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2010-2012, Willow Garage, Inc.
6  * Copyright (c) 2012-, Open Perception, Inc.
7  *
8  * All rights reserved.
9  *
10  * Redistribution and use in source and binary forms, with or without
11  * modification, are permitted provided that the following conditions
12  * are met:
13  *
14  * * Redistributions of source code must retain the above copyright
15  * notice, this list of conditions and the following disclaimer.
16  * * Redistributions in binary form must reproduce the above
17  * copyright notice, this list of conditions and the following
18  * disclaimer in the documentation and/or other materials provided
19  * with the distribution.
20  * * Neither the name of the copyright holder(s) nor the names of its
21  * contributors may be used to endorse or promote products derived
22  * from this software without specific prior written permission.
23  *
24  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
28  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
29  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
30  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
31  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
32  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
33  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
34  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
35  * POSSIBILITY OF SUCH DAMAGE.
36  */
37 
38 #pragma once
39 
40 #include <pcl/common/vector_average.h>
41 
42 
43 namespace pcl
44 {
45  template <typename real, int dimension>
47  {
48  reset();
49  }
50 
51  template <typename real, int dimension>
53  {
54  noOfSamples_ = 0;
55  accumulatedWeight_ = 0.0;
56  mean_.fill(0);
57  covariance_.fill(0);
58  }
59 
60  template <typename real, int dimension>
61  inline void VectorAverage<real, dimension>::add(const Eigen::Matrix<real, dimension, 1>& sample, real weight) {
62  if (weight == 0.0f)
63  return;
64 
65  ++noOfSamples_;
66  accumulatedWeight_ += weight;
67  real alpha = weight/accumulatedWeight_;
68 
69  Eigen::Matrix<real, dimension, 1> diff = sample - mean_;
70  covariance_ = (covariance_ + (diff * diff.transpose())*alpha)*(1.0f-alpha);
71 
72  mean_ += (diff)*alpha;
73 
74  //if (std::isnan(covariance_(0,0)))
75  //{
76  //std::cout << PVARN(weight);
77  //exit(0);
78  //}
79  }
80 
81  template <typename real, int dimension>
82  inline void VectorAverage<real, dimension>::doPCA(Eigen::Matrix<real, dimension, 1>& eigen_values, Eigen::Matrix<real, dimension, 1>& eigen_vector1,
83  Eigen::Matrix<real, dimension, 1>& eigen_vector2, Eigen::Matrix<real, dimension, 1>& eigen_vector3) const
84  {
85  // The following step is necessary for cases where the values in the covariance matrix are small
86  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
87  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
88  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance);
89  //eigen_values = ei_symm.eigenvalues().template cast<real>();
90  //Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors().template cast<real>();
91 
92  //std::cout << "My covariance is \n"<<covariance_<<"\n";
93  //std::cout << "My mean is \n"<<mean_<<"\n";
94  //std::cout << "My Eigenvectors \n"<<eigen_vectors<<"\n";
95 
96  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_);
97  eigen_values = ei_symm.eigenvalues();
98  Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors();
99 
100  eigen_vector1 = eigen_vectors.col(0);
101  eigen_vector2 = eigen_vectors.col(1);
102  eigen_vector3 = eigen_vectors.col(2);
103  }
104 
105  template <typename real, int dimension>
106  inline void VectorAverage<real, dimension>::doPCA(Eigen::Matrix<real, dimension, 1>& eigen_values) const
107  {
108  // The following step is necessary for cases where the values in the covariance matrix are small
109  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
110  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
111  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance, false);
112  //eigen_values = ei_symm.eigenvalues().template cast<real>();
113 
114  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_, false);
115  eigen_values = ei_symm.eigenvalues();
116  }
117 
118  template <typename real, int dimension>
119  inline void VectorAverage<real, dimension>::getEigenVector1(Eigen::Matrix<real, dimension, 1>& eigen_vector1) const
120  {
121  // The following step is necessary for cases where the values in the covariance matrix are small
122  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
123  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
124  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance);
125  //eigen_values = ei_symm.eigenvalues().template cast<real>();
126  //Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors().template cast<real>();
127 
128  //std::cout << "My covariance is \n"<<covariance_<<"\n";
129  //std::cout << "My mean is \n"<<mean_<<"\n";
130  //std::cout << "My Eigenvectors \n"<<eigen_vectors<<"\n";
131 
132  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_);
133  Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors();
134  eigen_vector1 = eigen_vectors.col(0);
135  }
136 
137 
138  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
139  // Special cases for real=float & dimension=3 -> Partial specialization does not work with class templates. :( //
140  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
141  ///////////
142  // float //
143  ///////////
144  template <>
145  inline void VectorAverage<float, 3>::doPCA(Eigen::Matrix<float, 3, 1>& eigen_values, Eigen::Matrix<float, 3, 1>& eigen_vector1,
146  Eigen::Matrix<float, 3, 1>& eigen_vector2, Eigen::Matrix<float, 3, 1>& eigen_vector3) const
147  {
148  //std::cout << "Using specialized 3x3 version of doPCA!\n";
149  Eigen::Matrix<float, 3, 3> eigen_vectors;
150  eigen33(covariance_, eigen_vectors, eigen_values);
151  eigen_vector1 = eigen_vectors.col(0);
152  eigen_vector2 = eigen_vectors.col(1);
153  eigen_vector3 = eigen_vectors.col(2);
154  }
155  template <>
156  inline void VectorAverage<float, 3>::doPCA(Eigen::Matrix<float, 3, 1>& eigen_values) const
157  {
158  //std::cout << "Using specialized 3x3 version of doPCA!\n";
159  computeRoots (covariance_, eigen_values);
160  }
161  template <>
162  inline void VectorAverage<float, 3>::getEigenVector1(Eigen::Matrix<float, 3, 1>& eigen_vector1) const
163  {
164  //std::cout << "Using specialized 3x3 version of doPCA!\n";
165  Eigen::Vector3f::Scalar eigen_value;
166  Eigen::Vector3f eigen_vector;
167  eigen33(covariance_, eigen_value, eigen_vector);
168  eigen_vector1 = eigen_vector;
169  }
170 
171  ////////////
172  // double //
173  ////////////
174  template <>
175  inline void VectorAverage<double, 3>::doPCA(Eigen::Matrix<double, 3, 1>& eigen_values, Eigen::Matrix<double, 3, 1>& eigen_vector1,
176  Eigen::Matrix<double, 3, 1>& eigen_vector2, Eigen::Matrix<double, 3, 1>& eigen_vector3) const
177  {
178  //std::cout << "Using specialized 3x3 version of doPCA!\n";
179  Eigen::Matrix<double, 3, 3> eigen_vectors;
180  eigen33(covariance_, eigen_vectors, eigen_values);
181  eigen_vector1 = eigen_vectors.col(0);
182  eigen_vector2 = eigen_vectors.col(1);
183  eigen_vector3 = eigen_vectors.col(2);
184  }
185  template <>
186  inline void VectorAverage<double, 3>::doPCA(Eigen::Matrix<double, 3, 1>& eigen_values) const
187  {
188  //std::cout << "Using specialized 3x3 version of doPCA!\n";
189  computeRoots (covariance_, eigen_values);
190  }
191  template <>
192  inline void VectorAverage<double, 3>::getEigenVector1(Eigen::Matrix<double, 3, 1>& eigen_vector1) const
193  {
194  //std::cout << "Using specialized 3x3 version of doPCA!\n";
195  Eigen::Vector3d::Scalar eigen_value;
196  Eigen::Vector3d eigen_vector;
197  eigen33(covariance_, eigen_value, eigen_vector);
198  eigen_vector1 = eigen_vector;
199  }
200 } // namespace pcl
void reset()
Reset the object to work with a new data set.
void computeRoots(const Matrix &m, Roots &roots)
computes the roots of the characteristic polynomial of the input matrix m, which are the eigenvalues ...
Definition: eigen.hpp:68
VectorAverage()
Constructor - dimension gives the size of the vectors to work with.
void doPCA(VectorType &eigen_values, VectorType &eigen_vector1, VectorType &eigen_vector2, VectorType &eigen_vector3) const
Do Principal component analysis.
void getEigenVector1(VectorType &eigen_vector1) const
Get the eigenvector corresponding to the smallest eigenvalue.
void eigen33(const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector)
determines the eigenvector and eigenvalue of the smallest eigenvalue of the symmetric positive semi d...
Definition: eigen.hpp:296
void add(const VectorType &sample, real weight=1.0)
Add a new sample.