Point Cloud Library (PCL)  1.7.1
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 #ifndef PCL_COMMON_VECTOR_AVERAGE_IMPL_HPP_
39 #define PCL_COMMON_VECTOR_AVERAGE_IMPL_HPP_
40 
41 namespace pcl
42 {
43  template <typename real, int dimension>
45  noOfSamples_ (0), accumulatedWeight_ (0),
46  mean_ (Eigen::Matrix<real, dimension, 1>::Identity ()),
47  covariance_ (Eigen::Matrix<real, dimension, dimension>::Identity ())
48  {
49  reset();
50  }
51 
52  template <typename real, int dimension>
54  {
55  noOfSamples_ = 0;
56  accumulatedWeight_ = 0.0;
57  mean_.fill(0);
58  covariance_.fill(0);
59  }
60 
61  template <typename real, int dimension>
62  inline void VectorAverage<real, dimension>::add(const Eigen::Matrix<real, dimension, 1>& sample, real weight) {
63  if (weight == 0.0f)
64  return;
65 
66  ++noOfSamples_;
67  accumulatedWeight_ += weight;
68  real alpha = weight/accumulatedWeight_;
69 
70  Eigen::Matrix<real, dimension, 1> diff = sample - mean_;
71  covariance_ = (covariance_ + (diff * diff.transpose())*alpha)*(1.0f-alpha);
72 
73  mean_ += (diff)*alpha;
74 
75  //if (pcl_isnan(covariance_(0,0)))
76  //{
77  //cout << PVARN(weight);
78  //exit(0);
79  //}
80  }
81 
82  template <typename real, int dimension>
83  inline void VectorAverage<real, dimension>::doPCA(Eigen::Matrix<real, dimension, 1>& eigen_values, Eigen::Matrix<real, dimension, 1>& eigen_vector1,
84  Eigen::Matrix<real, dimension, 1>& eigen_vector2, Eigen::Matrix<real, dimension, 1>& eigen_vector3) const
85  {
86  // The following step is necessary for cases where the values in the covariance matrix are small
87  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
88  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
89  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance);
90  //eigen_values = ei_symm.eigenvalues().template cast<real>();
91  //Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors().template cast<real>();
92 
93  //cout << "My covariance is \n"<<covariance_<<"\n";
94  //cout << "My mean is \n"<<mean_<<"\n";
95  //cout << "My Eigenvectors \n"<<eigen_vectors<<"\n";
96 
97  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_);
98  eigen_values = ei_symm.eigenvalues();
99  Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors();
100 
101  eigen_vector1 = eigen_vectors.col(0);
102  eigen_vector2 = eigen_vectors.col(1);
103  eigen_vector3 = eigen_vectors.col(2);
104  }
105 
106  template <typename real, int dimension>
107  inline void VectorAverage<real, dimension>::doPCA(Eigen::Matrix<real, dimension, 1>& eigen_values) const
108  {
109  // The following step is necessary for cases where the values in the covariance matrix are small
110  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
111  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
112  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance, false);
113  //eigen_values = ei_symm.eigenvalues().template cast<real>();
114 
115  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_, false);
116  eigen_values = ei_symm.eigenvalues();
117  }
118 
119  template <typename real, int dimension>
120  inline void VectorAverage<real, dimension>::getEigenVector1(Eigen::Matrix<real, dimension, 1>& eigen_vector1) const
121  {
122  // The following step is necessary for cases where the values in the covariance matrix are small
123  // In this case float accuracy is nor enough to calculate the eigenvalues and eigenvectors.
124  //Eigen::Matrix<double, dimension, dimension> tmp_covariance = covariance_.template cast<double>();
125  //Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, dimension, dimension> > ei_symm(tmp_covariance);
126  //eigen_values = ei_symm.eigenvalues().template cast<real>();
127  //Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors().template cast<real>();
128 
129  //cout << "My covariance is \n"<<covariance_<<"\n";
130  //cout << "My mean is \n"<<mean_<<"\n";
131  //cout << "My Eigenvectors \n"<<eigen_vectors<<"\n";
132 
133  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<real, dimension, dimension> > ei_symm(covariance_);
134  Eigen::Matrix<real, dimension, dimension> eigen_vectors = ei_symm.eigenvectors();
135  eigen_vector1 = eigen_vectors.col(0);
136  }
137 
138 
139  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
140  // Special cases for real=float & dimension=3 -> Partial specialization does not work with class templates. :( //
141  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
142  ///////////
143  // float //
144  ///////////
145  template <>
146  inline void VectorAverage<float, 3>::doPCA(Eigen::Matrix<float, 3, 1>& eigen_values, Eigen::Matrix<float, 3, 1>& eigen_vector1,
147  Eigen::Matrix<float, 3, 1>& eigen_vector2, Eigen::Matrix<float, 3, 1>& eigen_vector3) const
148  {
149  //cout << "Using specialized 3x3 version of doPCA!\n";
150  Eigen::Matrix<float, 3, 3> eigen_vectors;
151  eigen33(covariance_, eigen_vectors, eigen_values);
152  eigen_vector1 = eigen_vectors.col(0);
153  eigen_vector2 = eigen_vectors.col(1);
154  eigen_vector3 = eigen_vectors.col(2);
155  }
156  template <>
157  inline void VectorAverage<float, 3>::doPCA(Eigen::Matrix<float, 3, 1>& eigen_values) const
158  {
159  //cout << "Using specialized 3x3 version of doPCA!\n";
160  computeRoots (covariance_, eigen_values);
161  }
162  template <>
163  inline void VectorAverage<float, 3>::getEigenVector1(Eigen::Matrix<float, 3, 1>& eigen_vector1) const
164  {
165  //cout << "Using specialized 3x3 version of doPCA!\n";
166  Eigen::Vector3f::Scalar eigen_value;
167  Eigen::Vector3f eigen_vector;
168  eigen33(covariance_, eigen_value, eigen_vector);
169  eigen_vector1 = eigen_vector;
170  }
171 
172  ////////////
173  // double //
174  ////////////
175  template <>
176  inline void VectorAverage<double, 3>::doPCA(Eigen::Matrix<double, 3, 1>& eigen_values, Eigen::Matrix<double, 3, 1>& eigen_vector1,
177  Eigen::Matrix<double, 3, 1>& eigen_vector2, Eigen::Matrix<double, 3, 1>& eigen_vector3) const
178  {
179  //cout << "Using specialized 3x3 version of doPCA!\n";
180  Eigen::Matrix<double, 3, 3> eigen_vectors;
181  eigen33(covariance_, eigen_vectors, eigen_values);
182  eigen_vector1 = eigen_vectors.col(0);
183  eigen_vector2 = eigen_vectors.col(1);
184  eigen_vector3 = eigen_vectors.col(2);
185  }
186  template <>
187  inline void VectorAverage<double, 3>::doPCA(Eigen::Matrix<double, 3, 1>& eigen_values) const
188  {
189  //cout << "Using specialized 3x3 version of doPCA!\n";
190  computeRoots (covariance_, eigen_values);
191  }
192  template <>
193  inline void VectorAverage<double, 3>::getEigenVector1(Eigen::Matrix<double, 3, 1>& eigen_vector1) const
194  {
195  //cout << "Using specialized 3x3 version of doPCA!\n";
196  Eigen::Vector3d::Scalar eigen_value;
197  Eigen::Vector3d eigen_vector;
198  eigen33(covariance_, eigen_value, eigen_vector);
199  eigen_vector1 = eigen_vector;
200  }
201 } // END namespace
202 
203 #endif
204