// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Gael Guennebaud // Copyright (C) 2009 Benoit Jacob // // Eigen is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 3 of the License, or (at your option) any later version. // // Alternatively, you can redistribute it and/or // modify it under the terms of the GNU General Public License as // published by the Free Software Foundation; either version 2 of // the License, or (at your option) any later version. // // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the // GNU General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License and a copy of the GNU General Public License along with // Eigen. If not, see . #include "main.h" #include #include template void jacobisvd_check_full(const MatrixType& m, const JacobiSVD& svd) { typedef typename MatrixType::Index Index; Index rows = m.rows(); Index cols = m.cols(); enum { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime }; typedef typename MatrixType::Scalar Scalar; typedef typename NumTraits::Real RealScalar; typedef Matrix MatrixUType; typedef Matrix MatrixVType; typedef Matrix ColVectorType; typedef Matrix InputVectorType; MatrixType sigma = MatrixType::Zero(rows,cols); sigma.diagonal() = svd.singularValues().template cast(); MatrixUType u = svd.matrixU(); MatrixVType v = svd.matrixV(); VERIFY_IS_APPROX(m, u * sigma * v.adjoint()); VERIFY_IS_UNITARY(u); VERIFY_IS_UNITARY(v); } template void jacobisvd_compare_to_full(const MatrixType& m, unsigned int computationOptions, const JacobiSVD& referenceSvd) { typedef typename MatrixType::Index Index; Index rows = m.rows(); Index cols = m.cols(); Index diagSize = std::min(rows, cols); JacobiSVD svd(m, computationOptions); VERIFY_IS_EQUAL(svd.singularValues(), referenceSvd.singularValues()); if(computationOptions & ComputeFullU) VERIFY_IS_EQUAL(svd.matrixU(), referenceSvd.matrixU()); if(computationOptions & ComputeThinU) VERIFY_IS_EQUAL(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize)); if(computationOptions & ComputeFullV) VERIFY_IS_EQUAL(svd.matrixV(), referenceSvd.matrixV()); if(computationOptions & ComputeThinV) VERIFY_IS_EQUAL(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize)); } template void jacobisvd_test_all_computation_options(const MatrixType& m) { if (QRPreconditioner == NoQRPreconditioner && m.rows() != m.cols()) return; JacobiSVD fullSvd(m, ComputeFullU|ComputeFullV); jacobisvd_check_full(m, fullSvd); if(QRPreconditioner == FullPivHouseholderQRPreconditioner) return; jacobisvd_compare_to_full(m, ComputeFullU, fullSvd); jacobisvd_compare_to_full(m, ComputeFullV, fullSvd); jacobisvd_compare_to_full(m, 0, fullSvd); if (MatrixType::ColsAtCompileTime == Dynamic) { // thin U/V are only available with dynamic number of columns jacobisvd_compare_to_full(m, ComputeFullU|ComputeThinV, fullSvd); jacobisvd_compare_to_full(m, ComputeThinV, fullSvd); jacobisvd_compare_to_full(m, ComputeThinU|ComputeFullV, fullSvd); jacobisvd_compare_to_full(m, ComputeThinU , fullSvd); jacobisvd_compare_to_full(m, ComputeThinU|ComputeThinV, fullSvd); } } template void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true) { MatrixType m = pickrandom ? MatrixType::Random(a.rows(), a.cols()) : a; jacobisvd_test_all_computation_options(m); jacobisvd_test_all_computation_options(m); jacobisvd_test_all_computation_options(m); jacobisvd_test_all_computation_options(m); } template void jacobisvd_verify_assert() { MatrixType tmp; JacobiSVD svd; //VERIFY_RAISES_ASSERT(svd.solve(tmp, &tmp)) VERIFY_RAISES_ASSERT(svd.matrixU()) VERIFY_RAISES_ASSERT(svd.singularValues()) VERIFY_RAISES_ASSERT(svd.matrixV()) /*VERIFY_RAISES_ASSERT(svd.computeUnitaryPositive(&tmp,&tmp)) VERIFY_RAISES_ASSERT(svd.computePositiveUnitary(&tmp,&tmp)) VERIFY_RAISES_ASSERT(svd.computeRotationScaling(&tmp,&tmp)) VERIFY_RAISES_ASSERT(svd.computeScalingRotation(&tmp,&tmp))*/ } void test_jacobisvd() { for(int i = 0; i < g_repeat; i++) { Matrix2cd m; m << 0, 1, 0, 1; CALL_SUBTEST_1(( jacobisvd(m, false) )); m << 1, 0, 1, 0; CALL_SUBTEST_1(( jacobisvd(m, false) )); Matrix2d n; n << 1, 1, 1, -1; CALL_SUBTEST_2(( jacobisvd(n, false) )); CALL_SUBTEST_3(( jacobisvd() )); CALL_SUBTEST_4(( jacobisvd() )); CALL_SUBTEST_5(( jacobisvd >() )); CALL_SUBTEST_6(( jacobisvd >(Matrix(10,2)) )); CALL_SUBTEST_7(( jacobisvd(MatrixXf(50,50)) )); CALL_SUBTEST_8(( jacobisvd(MatrixXcd(14,7)) )); } CALL_SUBTEST_9(( jacobisvd(MatrixXf(300,200)) )); CALL_SUBTEST_10(( jacobisvd(MatrixXcd(100,150)) )); CALL_SUBTEST_3(( jacobisvd_verify_assert() )); CALL_SUBTEST_3(( jacobisvd_verify_assert() )); CALL_SUBTEST_9(( jacobisvd_verify_assert() )); CALL_SUBTEST_11(( jacobisvd_verify_assert() )); // Test problem size constructors CALL_SUBTEST_12( JacobiSVD(10, 20) ); }