00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #define EIGEN_STACK_ALLOCATION_LIMIT 0
00028 #define EIGEN_RUNTIME_NO_MALLOC
00029 #include "main.h"
00030 #include <Eigen/SVD>
00031
00032 template<typename MatrixType, int QRPreconditioner>
00033 void jacobisvd_check_full(const MatrixType& m, const JacobiSVD<MatrixType, QRPreconditioner>& svd)
00034 {
00035 typedef typename MatrixType::Index Index;
00036 Index rows = m.rows();
00037 Index cols = m.cols();
00038
00039 enum {
00040 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00041 ColsAtCompileTime = MatrixType::ColsAtCompileTime
00042 };
00043
00044 typedef typename MatrixType::Scalar Scalar;
00045 typedef typename NumTraits<Scalar>::Real RealScalar;
00046 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
00047 typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
00048 typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
00049 typedef Matrix<Scalar, ColsAtCompileTime, 1> InputVectorType;
00050
00051 MatrixType sigma = MatrixType::Zero(rows,cols);
00052 sigma.diagonal() = svd.singularValues().template cast<Scalar>();
00053 MatrixUType u = svd.matrixU();
00054 MatrixVType v = svd.matrixV();
00055
00056 VERIFY_IS_APPROX(m, u * sigma * v.adjoint());
00057 VERIFY_IS_UNITARY(u);
00058 VERIFY_IS_UNITARY(v);
00059 }
00060
00061 template<typename MatrixType, int QRPreconditioner>
00062 void jacobisvd_compare_to_full(const MatrixType& m,
00063 unsigned int computationOptions,
00064 const JacobiSVD<MatrixType, QRPreconditioner>& referenceSvd)
00065 {
00066 typedef typename MatrixType::Index Index;
00067 Index rows = m.rows();
00068 Index cols = m.cols();
00069 Index diagSize = std::min(rows, cols);
00070
00071 JacobiSVD<MatrixType, QRPreconditioner> svd(m, computationOptions);
00072
00073 VERIFY_IS_APPROX(svd.singularValues(), referenceSvd.singularValues());
00074 if(computationOptions & ComputeFullU)
00075 VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU());
00076 if(computationOptions & ComputeThinU)
00077 VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize));
00078 if(computationOptions & ComputeFullV)
00079 VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV());
00080 if(computationOptions & ComputeThinV)
00081 VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize));
00082 }
00083
00084 template<typename MatrixType, int QRPreconditioner>
00085 void jacobisvd_solve(const MatrixType& m, unsigned int computationOptions)
00086 {
00087 typedef typename MatrixType::Scalar Scalar;
00088 typedef typename MatrixType::Index Index;
00089 Index rows = m.rows();
00090 Index cols = m.cols();
00091
00092 enum {
00093 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00094 ColsAtCompileTime = MatrixType::ColsAtCompileTime
00095 };
00096
00097 typedef Matrix<Scalar, RowsAtCompileTime, Dynamic> RhsType;
00098 typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
00099
00100 RhsType rhs = RhsType::Random(rows, internal::random<Index>(1, cols));
00101 JacobiSVD<MatrixType, QRPreconditioner> svd(m, computationOptions);
00102 SolutionType x = svd.solve(rhs);
00103
00104 VERIFY_IS_APPROX(m.adjoint()*m*x,m.adjoint()*rhs);
00105 }
00106
00107 template<typename MatrixType, int QRPreconditioner>
00108 void jacobisvd_test_all_computation_options(const MatrixType& m)
00109 {
00110 if (QRPreconditioner == NoQRPreconditioner && m.rows() != m.cols())
00111 return;
00112 JacobiSVD<MatrixType, QRPreconditioner> fullSvd(m, ComputeFullU|ComputeFullV);
00113
00114 jacobisvd_check_full(m, fullSvd);
00115 jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeFullU | ComputeFullV);
00116
00117 if(QRPreconditioner == FullPivHouseholderQRPreconditioner)
00118 return;
00119
00120 jacobisvd_compare_to_full(m, ComputeFullU, fullSvd);
00121 jacobisvd_compare_to_full(m, ComputeFullV, fullSvd);
00122 jacobisvd_compare_to_full(m, 0, fullSvd);
00123
00124 if (MatrixType::ColsAtCompileTime == Dynamic) {
00125
00126 jacobisvd_compare_to_full(m, ComputeFullU|ComputeThinV, fullSvd);
00127 jacobisvd_compare_to_full(m, ComputeThinV, fullSvd);
00128 jacobisvd_compare_to_full(m, ComputeThinU|ComputeFullV, fullSvd);
00129 jacobisvd_compare_to_full(m, ComputeThinU , fullSvd);
00130 jacobisvd_compare_to_full(m, ComputeThinU|ComputeThinV, fullSvd);
00131 jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeFullU | ComputeThinV);
00132 jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeThinU | ComputeFullV);
00133 jacobisvd_solve<MatrixType, QRPreconditioner>(m, ComputeThinU | ComputeThinV);
00134 }
00135 }
00136
00137 template<typename MatrixType>
00138 void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
00139 {
00140 MatrixType m = pickrandom ? MatrixType::Random(a.rows(), a.cols()) : a;
00141
00142 jacobisvd_test_all_computation_options<MatrixType, FullPivHouseholderQRPreconditioner>(m);
00143 jacobisvd_test_all_computation_options<MatrixType, ColPivHouseholderQRPreconditioner>(m);
00144 jacobisvd_test_all_computation_options<MatrixType, HouseholderQRPreconditioner>(m);
00145 jacobisvd_test_all_computation_options<MatrixType, NoQRPreconditioner>(m);
00146 }
00147
00148 template<typename MatrixType> void jacobisvd_verify_assert(const MatrixType& m)
00149 {
00150 typedef typename MatrixType::Scalar Scalar;
00151 typedef typename MatrixType::Index Index;
00152 Index rows = m.rows();
00153 Index cols = m.cols();
00154
00155 enum {
00156 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00157 ColsAtCompileTime = MatrixType::ColsAtCompileTime
00158 };
00159
00160 typedef Matrix<Scalar, RowsAtCompileTime, 1> RhsType;
00161
00162 RhsType rhs(rows);
00163
00164 JacobiSVD<MatrixType> svd;
00165 VERIFY_RAISES_ASSERT(svd.matrixU())
00166 VERIFY_RAISES_ASSERT(svd.singularValues())
00167 VERIFY_RAISES_ASSERT(svd.matrixV())
00168 VERIFY_RAISES_ASSERT(svd.solve(rhs))
00169
00170 MatrixType a = MatrixType::Zero(rows, cols);
00171 a.setZero();
00172 svd.compute(a, 0);
00173 VERIFY_RAISES_ASSERT(svd.matrixU())
00174 VERIFY_RAISES_ASSERT(svd.matrixV())
00175 svd.singularValues();
00176 VERIFY_RAISES_ASSERT(svd.solve(rhs))
00177
00178 if (ColsAtCompileTime == Dynamic)
00179 {
00180 svd.compute(a, ComputeThinU);
00181 svd.matrixU();
00182 VERIFY_RAISES_ASSERT(svd.matrixV())
00183 VERIFY_RAISES_ASSERT(svd.solve(rhs))
00184
00185 svd.compute(a, ComputeThinV);
00186 svd.matrixV();
00187 VERIFY_RAISES_ASSERT(svd.matrixU())
00188 VERIFY_RAISES_ASSERT(svd.solve(rhs))
00189
00190 JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> svd_fullqr;
00191 VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeFullU|ComputeThinV))
00192 VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeThinV))
00193 VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeFullV))
00194 }
00195 else
00196 {
00197 VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinU))
00198 VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinV))
00199 }
00200 }
00201
00202 template<typename MatrixType>
00203 void jacobisvd_method()
00204 {
00205 enum { Size = MatrixType::RowsAtCompileTime };
00206 typedef typename MatrixType::RealScalar RealScalar;
00207 typedef Matrix<RealScalar, Size, 1> RealVecType;
00208 MatrixType m = MatrixType::Identity();
00209 VERIFY_IS_APPROX(m.jacobiSvd().singularValues(), RealVecType::Ones());
00210 VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixU());
00211 VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixV());
00212 VERIFY_IS_APPROX(m.jacobiSvd(ComputeFullU|ComputeFullV).solve(m), m);
00213 }
00214
00215
00216
00217 template<typename Scalar>
00218 EIGEN_DONT_INLINE Scalar zero() { return Scalar(0); }
00219
00220
00221 template<typename T> EIGEN_DONT_INLINE T sub(T a, T b) { return a - b; }
00222
00223 template<typename MatrixType>
00224 void jacobisvd_inf_nan()
00225 {
00226
00227
00228 JacobiSVD<MatrixType> svd;
00229 typedef typename MatrixType::Scalar Scalar;
00230 Scalar some_inf = Scalar(1) / zero<Scalar>();
00231 VERIFY(sub(some_inf, some_inf) != sub(some_inf, some_inf));
00232 svd.compute(MatrixType::Constant(10,10,some_inf), ComputeFullU | ComputeFullV);
00233
00234 Scalar some_nan = zero<Scalar>() / zero<Scalar>();
00235 VERIFY(some_nan != some_nan);
00236 svd.compute(MatrixType::Constant(10,10,some_nan), ComputeFullU | ComputeFullV);
00237
00238 MatrixType m = MatrixType::Zero(10,10);
00239 m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_inf;
00240 svd.compute(m, ComputeFullU | ComputeFullV);
00241
00242 m = MatrixType::Zero(10,10);
00243 m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_nan;
00244 svd.compute(m, ComputeFullU | ComputeFullV);
00245 }
00246
00247 void jacobisvd_preallocate()
00248 {
00249 Vector3f v(3.f, 2.f, 1.f);
00250 MatrixXf m = v.asDiagonal();
00251
00252 internal::set_is_malloc_allowed(false);
00253 VERIFY_RAISES_ASSERT(VectorXf v(10);)
00254 JacobiSVD<MatrixXf> svd;
00255 internal::set_is_malloc_allowed(true);
00256 svd.compute(m);
00257 VERIFY_IS_APPROX(svd.singularValues(), v);
00258
00259 JacobiSVD<MatrixXf> svd2(3,3);
00260 internal::set_is_malloc_allowed(false);
00261 svd2.compute(m);
00262 internal::set_is_malloc_allowed(true);
00263 VERIFY_IS_APPROX(svd2.singularValues(), v);
00264 VERIFY_RAISES_ASSERT(svd2.matrixU());
00265 VERIFY_RAISES_ASSERT(svd2.matrixV());
00266 svd2.compute(m, ComputeFullU | ComputeFullV);
00267 VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
00268 VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
00269 internal::set_is_malloc_allowed(false);
00270 svd2.compute(m);
00271 internal::set_is_malloc_allowed(true);
00272
00273 JacobiSVD<MatrixXf> svd3(3,3,ComputeFullU|ComputeFullV);
00274 internal::set_is_malloc_allowed(false);
00275 svd2.compute(m);
00276 internal::set_is_malloc_allowed(true);
00277 VERIFY_IS_APPROX(svd2.singularValues(), v);
00278 VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
00279 VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
00280 internal::set_is_malloc_allowed(false);
00281 svd2.compute(m, ComputeFullU|ComputeFullV);
00282 internal::set_is_malloc_allowed(true);
00283
00284
00285 }
00286
00287 void test_jacobisvd()
00288 {
00289 CALL_SUBTEST_3(( jacobisvd_verify_assert(Matrix3f()) ));
00290 CALL_SUBTEST_4(( jacobisvd_verify_assert(Matrix4d()) ));
00291 CALL_SUBTEST_7(( jacobisvd_verify_assert(MatrixXf(10,12)) ));
00292 CALL_SUBTEST_8(( jacobisvd_verify_assert(MatrixXcd(7,5)) ));
00293
00294 for(int i = 0; i < g_repeat; i++) {
00295 Matrix2cd m;
00296 m << 0, 1,
00297 0, 1;
00298 CALL_SUBTEST_1(( jacobisvd(m, false) ));
00299 m << 1, 0,
00300 1, 0;
00301 CALL_SUBTEST_1(( jacobisvd(m, false) ));
00302
00303 Matrix2d n;
00304 n << 0, 0,
00305 0, 0;
00306 CALL_SUBTEST_2(( jacobisvd(n, false) ));
00307 n << 0, 0,
00308 0, 1;
00309 CALL_SUBTEST_2(( jacobisvd(n, false) ));
00310
00311 CALL_SUBTEST_3(( jacobisvd<Matrix3f>() ));
00312 CALL_SUBTEST_4(( jacobisvd<Matrix4d>() ));
00313 CALL_SUBTEST_5(( jacobisvd<Matrix<float,3,5> >() ));
00314 CALL_SUBTEST_6(( jacobisvd<Matrix<double,Dynamic,2> >(Matrix<double,Dynamic,2>(10,2)) ));
00315
00316 int r = internal::random<int>(1, 30),
00317 c = internal::random<int>(1, 30);
00318 CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(r,c)) ));
00319 CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(r,c)) ));
00320 (void) r;
00321 (void) c;
00322
00323
00324 CALL_SUBTEST_7( jacobisvd_inf_nan<MatrixXf>() );
00325 }
00326
00327 CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(internal::random<int>(100, 150), internal::random<int>(100, 150))) ));
00328 CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(internal::random<int>(80, 100), internal::random<int>(80, 100))) ));
00329
00330
00331 CALL_SUBTEST_1(( jacobisvd_method<Matrix2cd>() ));
00332 CALL_SUBTEST_3(( jacobisvd_method<Matrix3f>() ));
00333
00334
00335 CALL_SUBTEST_7( JacobiSVD<MatrixXf>(10,10) );
00336
00337
00338 CALL_SUBTEST_9( jacobisvd_preallocate() );
00339 }