aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/jacobisvd.cpp
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-08-31 22:26:15 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-08-31 22:26:15 -0400
commit6e4e94ff3266f8d85adbfe6556e612a2ff9a7e68 (patch)
tree33a1923d88f095df3f98edd80cb1b02fc7799b91 /test/jacobisvd.cpp
parent29c6b2452dbe82cd49aa701921f2fa5a20017cc0 (diff)
* JacobiSVD:
- support complex numbers - big rewrite of the 2x2 kernel, much more robust * Jacobi: - fix weirdness in initial design, e.g. applyJacobiOnTheRight actually did the inverse transformation - fully support complex numbers - fix logic to decide whether to vectorize - remove several clumsy methods fix for complex numbers
Diffstat (limited to 'test/jacobisvd.cpp')
-rw-r--r--test/jacobisvd.cpp105
1 files changed, 105 insertions, 0 deletions
diff --git a/test/jacobisvd.cpp b/test/jacobisvd.cpp
new file mode 100644
index 000000000..9a4d79e45
--- /dev/null
+++ b/test/jacobisvd.cpp
@@ -0,0 +1,105 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
+// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
+//
+// 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 <http://www.gnu.org/licenses/>.
+
+#include "main.h"
+#include <Eigen/SVD>
+#include <Eigen/LU>
+
+template<typename MatrixType> void svd(const MatrixType& m, bool pickrandom = true)
+{
+ int rows = m.rows();
+ int cols = m.cols();
+
+ enum {
+ RowsAtCompileTime = MatrixType::RowsAtCompileTime,
+ ColsAtCompileTime = MatrixType::ColsAtCompileTime
+ };
+
+ typedef typename MatrixType::Scalar Scalar;
+ typedef typename NumTraits<Scalar>::Real RealScalar;
+ typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
+ typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
+ typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
+ typedef Matrix<Scalar, ColsAtCompileTime, 1> InputVectorType;
+
+ MatrixType a;
+ if(pickrandom) a = MatrixType::Random(rows,cols);
+ else a = m;
+
+ JacobiSVD<MatrixType> svd(a);
+ MatrixType sigma = MatrixType::Zero(rows,cols);
+ sigma.diagonal() = svd.singularValues().template cast<Scalar>();
+ MatrixUType u = svd.matrixU();
+ MatrixVType v = svd.matrixV();
+
+ VERIFY_IS_APPROX(a, u * sigma * v.adjoint());
+ VERIFY_IS_UNITARY(u);
+ VERIFY_IS_UNITARY(v);
+}
+
+template<typename MatrixType> void svd_verify_assert()
+{
+ MatrixType tmp;
+
+ SVD<MatrixType> 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( svd(m, false) );
+ m << 1, 0,
+ 1, 0;
+ CALL_SUBTEST( svd(m, false) );
+ Matrix2d n;
+ n << 1, 1,
+ 1, -1;
+ CALL_SUBTEST( svd(n, false) );
+ CALL_SUBTEST( svd(Matrix3f()) );
+ CALL_SUBTEST( svd(Matrix4d()) );
+ CALL_SUBTEST( svd(MatrixXf(50,50)) );
+// CALL_SUBTEST( svd(MatrixXd(14,7)) );
+ CALL_SUBTEST( svd(MatrixXcf(3,3)) );
+ CALL_SUBTEST( svd(MatrixXd(30,30)) );
+ }
+ CALL_SUBTEST( svd(MatrixXf(200,200)) );
+ CALL_SUBTEST( svd(MatrixXcd(100,100)) );
+
+ CALL_SUBTEST( svd_verify_assert<Matrix3f>() );
+ CALL_SUBTEST( svd_verify_assert<Matrix3d>() );
+ CALL_SUBTEST( svd_verify_assert<MatrixXf>() );
+ CALL_SUBTEST( svd_verify_assert<MatrixXd>() );
+}