aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/svd.cpp
blob: e6a32bd3f85f801b9976616cf757443ad4c03afa (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// 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)
{
  /* this test covers the following files:
     SVD.h
  */
  int rows = m.rows();
  int cols = m.cols();

  typedef typename MatrixType::Scalar Scalar;
  typedef typename NumTraits<Scalar>::Real RealScalar;
  MatrixType a = MatrixType::Random(rows,cols);
  Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> b =
    Matrix<Scalar, MatrixType::RowsAtCompileTime, 1>::Random(rows,1);
  Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> x(cols,1), x2(cols,1);

  {
    SVD<MatrixType> svd(a);
    MatrixType sigma = MatrixType::Zero(rows,cols);
    MatrixType matU  = MatrixType::Zero(rows,rows);
    sigma.diagonal() = svd.singularValues();
    matU = svd.matrixU();
    VERIFY_IS_APPROX(a, matU * sigma * svd.matrixV().transpose());
  }


  if (rows==cols)
  {
    if (ei_is_same_type<RealScalar,float>::ret)
    {
      MatrixType a1 = MatrixType::Random(rows,cols);
      a += a * a.adjoint() + a1 * a1.adjoint();
    }
    SVD<MatrixType> svd(a);
    svd.solve(b, &x);
    VERIFY_IS_APPROX(a * x,b);
  }


  if(rows==cols)
  {
    SVD<MatrixType> svd(a);
    MatrixType unitary, positive;
    svd.computeUnitaryPositive(&unitary, &positive);
    VERIFY_IS_APPROX(unitary * unitary.adjoint(), MatrixType::Identity(unitary.rows(),unitary.rows()));
    VERIFY_IS_APPROX(positive, positive.adjoint());
    for(int i = 0; i < rows; i++) VERIFY(positive.diagonal()[i] >= 0); // cheap necessary (not sufficient) condition for positivity
    VERIFY_IS_APPROX(unitary*positive, a);

    svd.computePositiveUnitary(&positive, &unitary);
    VERIFY_IS_APPROX(unitary * unitary.adjoint(), MatrixType::Identity(unitary.rows(),unitary.rows()));
    VERIFY_IS_APPROX(positive, positive.adjoint());
    for(int i = 0; i < rows; i++) VERIFY(positive.diagonal()[i] >= 0); // cheap necessary (not sufficient) condition for positivity
    VERIFY_IS_APPROX(positive*unitary, a);
  }
}

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_svd()
{
  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST( svd(Matrix3f()) );
    CALL_SUBTEST( svd(Matrix4d()) );
    CALL_SUBTEST( svd(MatrixXf(7,7)) );
    CALL_SUBTEST( svd(MatrixXd(14,7)) );
    // complex are not implemented yet
//     CALL_SUBTEST( svd(MatrixXcd(6,6)) );
//     CALL_SUBTEST( svd(MatrixXcf(3,3)) );
  }

  CALL_SUBTEST( svd_verify_assert<Matrix3f>() );
  CALL_SUBTEST( svd_verify_assert<Matrix3d>() );
  CALL_SUBTEST( svd_verify_assert<MatrixXf>() );
  CALL_SUBTEST( svd_verify_assert<MatrixXd>() );
}