svd_common.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
00006 //
00007 // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
00008 // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
00009 // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
00010 // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
00011 //
00012 // This Source Code Form is subject to the terms of the Mozilla
00013 // Public License v. 2.0. If a copy of the MPL was not distributed
00014 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00015 
00016 // discard stack allocation as that too bypasses malloc
00017 #define EIGEN_STACK_ALLOCATION_LIMIT 0
00018 #define EIGEN_RUNTIME_NO_MALLOC
00019 
00020 #include "main.h"
00021 #include <unsupported/Eigen/SVD>
00022 #include <Eigen/LU>
00023 
00024 
00025 // check if "svd" is the good image of "m"  
00026 template<typename MatrixType, typename SVD>
00027 void svd_check_full(const MatrixType& m, const SVD& svd)
00028 {
00029   typedef typename MatrixType::Index Index;
00030   Index rows = m.rows();
00031   Index cols = m.cols();
00032   enum {
00033     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00034     ColsAtCompileTime = MatrixType::ColsAtCompileTime
00035   };
00036 
00037   typedef typename MatrixType::Scalar Scalar;
00038   typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
00039   typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
00040 
00041   
00042   MatrixType sigma = MatrixType::Zero(rows, cols);
00043   sigma.diagonal() = svd.singularValues().template cast<Scalar>();
00044   MatrixUType u = svd.matrixU();
00045   MatrixVType v = svd.matrixV();
00046   VERIFY_IS_APPROX(m, u * sigma * v.adjoint());
00047   VERIFY_IS_UNITARY(u);
00048   VERIFY_IS_UNITARY(v);
00049 } // end svd_check_full
00050 
00051 
00052 
00053 // Compare to a reference value
00054 template<typename MatrixType, typename SVD>
00055 void svd_compare_to_full(const MatrixType& m,
00056                          unsigned int computationOptions,
00057                          const SVD& referenceSvd)
00058 {
00059   typedef typename MatrixType::Index Index;
00060   Index rows = m.rows();
00061   Index cols = m.cols();
00062   Index diagSize = (std::min)(rows, cols);
00063 
00064   SVD svd(m, computationOptions);
00065 
00066   VERIFY_IS_APPROX(svd.singularValues(), referenceSvd.singularValues());
00067   if(computationOptions & ComputeFullU)
00068     VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU());
00069   if(computationOptions & ComputeThinU)
00070     VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize));
00071   if(computationOptions & ComputeFullV)
00072     VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV());
00073   if(computationOptions & ComputeThinV)
00074     VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize));
00075 } // end svd_compare_to_full
00076 
00077 
00078 
00079 template<typename MatrixType, typename SVD>
00080 void svd_solve(const MatrixType& m, unsigned int computationOptions)
00081 {
00082   typedef typename MatrixType::Scalar Scalar;
00083   typedef typename MatrixType::Index Index;
00084   Index rows = m.rows();
00085   Index cols = m.cols();
00086 
00087   enum {
00088     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00089     ColsAtCompileTime = MatrixType::ColsAtCompileTime
00090   };
00091 
00092   typedef Matrix<Scalar, RowsAtCompileTime, Dynamic> RhsType;
00093   typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
00094 
00095   RhsType rhs = RhsType::Random(rows, internal::random<Index>(1, cols));
00096   SVD svd(m, computationOptions);
00097   SolutionType x = svd.solve(rhs);
00098   // evaluate normal equation which works also for least-squares solutions
00099   VERIFY_IS_APPROX(m.adjoint()*m*x,m.adjoint()*rhs);
00100 } // end svd_solve
00101 
00102 
00103 // test computations options
00104 // 2 functions because Jacobisvd can return before the second function
00105 template<typename MatrixType, typename SVD>
00106 void svd_test_computation_options_1(const MatrixType& m, const SVD& fullSvd)
00107 {
00108   svd_check_full< MatrixType, SVD >(m, fullSvd);
00109   svd_solve< MatrixType, SVD >(m, ComputeFullU | ComputeFullV);
00110 }
00111 
00112 
00113 template<typename MatrixType, typename SVD>
00114 void svd_test_computation_options_2(const MatrixType& m, const SVD& fullSvd)
00115 {
00116   svd_compare_to_full< MatrixType, SVD >(m, ComputeFullU, fullSvd);
00117   svd_compare_to_full< MatrixType, SVD >(m, ComputeFullV, fullSvd);
00118   svd_compare_to_full< MatrixType, SVD >(m, 0, fullSvd);
00119 
00120   if (MatrixType::ColsAtCompileTime == Dynamic) {
00121     // thin U/V are only available with dynamic number of columns
00122  
00123     svd_compare_to_full< MatrixType, SVD >(m, ComputeFullU|ComputeThinV, fullSvd);
00124     svd_compare_to_full< MatrixType, SVD >(m,              ComputeThinV, fullSvd);
00125     svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU|ComputeFullV, fullSvd);
00126     svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU             , fullSvd);
00127     svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU|ComputeThinV, fullSvd);
00128     svd_solve<MatrixType, SVD>(m, ComputeFullU | ComputeThinV);
00129     svd_solve<MatrixType, SVD>(m, ComputeThinU | ComputeFullV);
00130     svd_solve<MatrixType, SVD>(m, ComputeThinU | ComputeThinV);
00131     
00132     typedef typename MatrixType::Index Index;
00133     Index diagSize = (std::min)(m.rows(), m.cols());
00134     SVD svd(m, ComputeThinU | ComputeThinV);
00135     VERIFY_IS_APPROX(m, svd.matrixU().leftCols(diagSize) * svd.singularValues().asDiagonal() * svd.matrixV().leftCols(diagSize).adjoint());
00136   }
00137 }
00138 
00139 template<typename MatrixType, typename SVD> 
00140 void svd_verify_assert(const MatrixType& m)
00141 {
00142   typedef typename MatrixType::Scalar Scalar;
00143   typedef typename MatrixType::Index Index;
00144   Index rows = m.rows();
00145   Index cols = m.cols();
00146 
00147   enum {
00148     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00149     ColsAtCompileTime = MatrixType::ColsAtCompileTime
00150   };
00151 
00152   typedef Matrix<Scalar, RowsAtCompileTime, 1> RhsType;
00153   RhsType rhs(rows);
00154   SVD svd;
00155   VERIFY_RAISES_ASSERT(svd.matrixU())
00156   VERIFY_RAISES_ASSERT(svd.singularValues())
00157   VERIFY_RAISES_ASSERT(svd.matrixV())
00158   VERIFY_RAISES_ASSERT(svd.solve(rhs))
00159   MatrixType a = MatrixType::Zero(rows, cols);
00160   a.setZero();
00161   svd.compute(a, 0);
00162   VERIFY_RAISES_ASSERT(svd.matrixU())
00163   VERIFY_RAISES_ASSERT(svd.matrixV())
00164   svd.singularValues();
00165   VERIFY_RAISES_ASSERT(svd.solve(rhs))
00166     
00167   if (ColsAtCompileTime == Dynamic)
00168   {
00169     svd.compute(a, ComputeThinU);
00170     svd.matrixU();
00171     VERIFY_RAISES_ASSERT(svd.matrixV())
00172     VERIFY_RAISES_ASSERT(svd.solve(rhs))
00173     svd.compute(a, ComputeThinV);
00174     svd.matrixV();
00175     VERIFY_RAISES_ASSERT(svd.matrixU())
00176     VERIFY_RAISES_ASSERT(svd.solve(rhs))
00177   }
00178   else
00179   {
00180     VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinU))
00181     VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinV))
00182   }
00183 }
00184 
00185 // work around stupid msvc error when constructing at compile time an expression that involves
00186 // a division by zero, even if the numeric type has floating point
00187 template<typename Scalar>
00188 EIGEN_DONT_INLINE Scalar zero() { return Scalar(0); }
00189 
00190 // workaround aggressive optimization in ICC
00191 template<typename T> EIGEN_DONT_INLINE  T sub(T a, T b) { return a - b; }
00192 
00193 
00194 template<typename MatrixType, typename SVD>
00195 void svd_inf_nan()
00196 {
00197   // all this function does is verify we don't iterate infinitely on nan/inf values
00198 
00199   SVD svd;
00200   typedef typename MatrixType::Scalar Scalar;
00201   Scalar some_inf = Scalar(1) / zero<Scalar>();
00202   VERIFY(sub(some_inf, some_inf) != sub(some_inf, some_inf));
00203   svd.compute(MatrixType::Constant(10,10,some_inf), ComputeFullU | ComputeFullV);
00204 
00205   Scalar some_nan = zero<Scalar> () / zero<Scalar> ();
00206   VERIFY(some_nan != some_nan);
00207   svd.compute(MatrixType::Constant(10,10,some_nan), ComputeFullU | ComputeFullV);
00208 
00209   MatrixType m = MatrixType::Zero(10,10);
00210   m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_inf;
00211   svd.compute(m, ComputeFullU | ComputeFullV);
00212 
00213   m = MatrixType::Zero(10,10);
00214   m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_nan;
00215   svd.compute(m, ComputeFullU | ComputeFullV);
00216 }
00217 
00218 
00219 template<typename SVD>
00220 void svd_preallocate()
00221 {
00222   Vector3f v(3.f, 2.f, 1.f);
00223   MatrixXf m = v.asDiagonal();
00224 
00225   internal::set_is_malloc_allowed(false);
00226   VERIFY_RAISES_ASSERT(VectorXf v(10);)
00227     SVD svd;
00228   internal::set_is_malloc_allowed(true);
00229   svd.compute(m);
00230   VERIFY_IS_APPROX(svd.singularValues(), v);
00231 
00232   SVD svd2(3,3);
00233   internal::set_is_malloc_allowed(false);
00234   svd2.compute(m);
00235   internal::set_is_malloc_allowed(true);
00236   VERIFY_IS_APPROX(svd2.singularValues(), v);
00237   VERIFY_RAISES_ASSERT(svd2.matrixU());
00238   VERIFY_RAISES_ASSERT(svd2.matrixV());
00239   svd2.compute(m, ComputeFullU | ComputeFullV);
00240   VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
00241   VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
00242   internal::set_is_malloc_allowed(false);
00243   svd2.compute(m);
00244   internal::set_is_malloc_allowed(true);
00245 
00246   SVD svd3(3,3,ComputeFullU|ComputeFullV);
00247   internal::set_is_malloc_allowed(false);
00248   svd2.compute(m);
00249   internal::set_is_malloc_allowed(true);
00250   VERIFY_IS_APPROX(svd2.singularValues(), v);
00251   VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
00252   VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
00253   internal::set_is_malloc_allowed(false);
00254   svd2.compute(m, ComputeFullU|ComputeFullV);
00255   internal::set_is_malloc_allowed(true);
00256 }
00257 
00258 
00259 
00260 
00261 


shape_reconstruction
Author(s): Roberto Martín-Martín
autogenerated on Sat Jun 8 2019 18:36:57