aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SVD/JacobiSquareSVD.h
blob: 18c3db980a91d024ff38a3d1020c2e44a82f56bb (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
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// 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/>.

#ifndef EIGEN_JACOBISQUARESVD_H
#define EIGEN_JACOBISQUARESVD_H

/** \ingroup SVD_Module
  * \nonstableyet
  *
  * \class JacobiSquareSVD
  *
  * \brief Jacobi SVD decomposition of a square matrix
  *
  * \param MatrixType the type of the matrix of which we are computing the SVD decomposition
  * \param ComputeU whether the U matrix should be computed
  * \param ComputeV whether the V matrix should be computed
  *
  * \sa MatrixBase::jacobiSvd()
  */
template<typename MatrixType, bool ComputeU, bool ComputeV> class JacobiSquareSVD
{
  private:
    typedef typename MatrixType::Scalar Scalar;
    typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
    enum {
      RowsAtCompileTime = MatrixType::RowsAtCompileTime,
      MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
      Options = MatrixType::Options
    };
    
    typedef Matrix<Scalar, Dynamic, Dynamic, Options> DummyMatrixType;
    typedef typename ei_meta_if<ComputeU,
                                Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
                                       Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime>,
                                DummyMatrixType>::ret MatrixUType;
    typedef typename Diagonal<MatrixType,0>::PlainMatrixType SingularValuesType;
    typedef Matrix<Scalar, 1, RowsAtCompileTime, Options, 1, MaxRowsAtCompileTime> RowType;
    typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> ColType;

  public:

    JacobiSquareSVD() : m_isInitialized(false) {}

    JacobiSquareSVD(const MatrixType& matrix) : m_isInitialized(false) 
    {
      compute(matrix);
    }
    
    void compute(const MatrixType& matrix);
    
    const MatrixUType& matrixU() const
    {
      ei_assert(m_isInitialized && "SVD is not initialized.");
      return m_matrixU;
    }

    const SingularValuesType& singularValues() const
    {
      ei_assert(m_isInitialized && "SVD is not initialized.");
      return m_singularValues;
    }

    const MatrixUType& matrixV() const
    {
      ei_assert(m_isInitialized && "SVD is not initialized.");
      return m_matrixV;
    }

  protected:
    MatrixUType m_matrixU;
    MatrixUType m_matrixV;
    SingularValuesType m_singularValues;
    bool m_isInitialized;
};

template<typename MatrixType, bool ComputeU, bool ComputeV>
void JacobiSquareSVD<MatrixType, ComputeU, ComputeV>::compute(const MatrixType& matrix)
{
  MatrixType work_matrix(matrix);
  int size = matrix.rows();
  if(ComputeU) m_matrixU = MatrixUType::Identity(size,size);
  if(ComputeV) m_matrixV = MatrixUType::Identity(size,size);
  m_singularValues.resize(size);

sweep_again:
  for(int p = 1; p < size; ++p)
  {
    for(int q = 0; q < p; ++q)
    {
      Scalar c, s;
      while(std::max(ei_abs(work_matrix.coeff(p,q)),ei_abs(work_matrix.coeff(q,p)))
          > std::max(ei_abs(work_matrix.coeff(p,p)),ei_abs(work_matrix.coeff(q,q)))*machine_epsilon<Scalar>())
      {  
        if(work_matrix.makeJacobiForAtA(p,q,&c,&s))
        {
          work_matrix.applyJacobiOnTheRight(p,q,c,s);
          if(ComputeV) m_matrixV.applyJacobiOnTheRight(p,q,c,s);
        }
        if(work_matrix.makeJacobiForAAt(p,q,&c,&s))
        {
          Scalar x = ei_abs2(work_matrix.coeff(p,p)) + ei_abs2(work_matrix.coeff(p,q));
          Scalar y = ei_conj(work_matrix.coeff(q,p))*work_matrix.coeff(p,p) + ei_conj(work_matrix.coeff(q,q))*work_matrix.coeff(p,q);
          Scalar z = ei_abs2(work_matrix.coeff(q,p)) + ei_abs2(work_matrix.coeff(q,q));
          work_matrix.applyJacobiOnTheLeft(p,q,c,s);
          if(ComputeU) m_matrixU.applyJacobiOnTheRight(p,q,c,s);
          if(std::max(ei_abs(work_matrix.coeff(p,q)),ei_abs(work_matrix.coeff(q,p)))
           > std::max(ei_abs(work_matrix.coeff(p,p)),ei_abs(work_matrix.coeff(q,q))) )
          {
            work_matrix.row(p).swap(work_matrix.row(q));
            if(ComputeU) m_matrixU.col(p).swap(m_matrixU.col(q));
          }
        }
      }
    }
  }
  
  RealScalar biggestOnDiag = work_matrix.diagonal().cwise().abs().maxCoeff();
  RealScalar maxAllowedOffDiag = biggestOnDiag * machine_epsilon<Scalar>();
  for(int p = 0; p < size; ++p)
  {
    for(int q = 0; q < p; ++q)
      if(ei_abs(work_matrix.coeff(p,q)) > maxAllowedOffDiag)
        goto sweep_again;
    for(int q = p+1; q < size; ++q)
      if(ei_abs(work_matrix.coeff(p,q)) > maxAllowedOffDiag)
        goto sweep_again;
  }

  m_singularValues = work_matrix.diagonal().cwise().abs();
  RealScalar biggestSingularValue = m_singularValues.maxCoeff();

  for(int i = 0; i < size; ++i)
  {
    RealScalar a = ei_abs(work_matrix.coeff(i,i));
    m_singularValues.coeffRef(i) = a;
    if(ComputeU && !ei_isMuchSmallerThan(a, biggestSingularValue)) m_matrixU.col(i) *= work_matrix.coeff(i,i)/a;
  }

  for(int i = 0; i < size; i++)
  {
    int pos;
    m_singularValues.end(size-i).maxCoeff(&pos);
    if(pos)
    {
      pos += i;
      std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
      if(ComputeU) m_matrixU.col(pos).swap(m_matrixU.col(i));
      if(ComputeV) m_matrixV.col(pos).swap(m_matrixV.col(i));
    }
  }

  m_isInitialized = true;
}
#endif // EIGEN_JACOBISQUARESVD_H