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
176
177
178
179
180
181
182
183
184
185
186
|
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// 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/>.
#ifndef EIGEN_CHOLESKY_H
#define EIGEN_CHOLESKY_H
/** \ingroup Cholesky_Module
*
* \class Cholesky
*
* \deprecated this class has been renamed LLT
*/
template<typename MatrixType> class Cholesky
{
private:
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
enum {
PacketSize = ei_packet_traits<Scalar>::size,
AlignmentMask = int(PacketSize)-1
};
public:
Cholesky(const MatrixType& matrix)
: m_matrix(matrix.rows(), matrix.cols())
{
compute(matrix);
}
inline Part<MatrixType, Lower> matrixL(void) const { return m_matrix; }
/** \returns true if the matrix is positive definite */
inline bool isPositiveDefinite(void) const { return m_isPositiveDefinite; }
template<typename Derived>
typename Derived::Eval solve(const MatrixBase<Derived> &b) const EIGEN_DEPRECATED;
template<typename RhsDerived, typename ResDerived>
bool solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const;
template<typename Derived>
bool solveInPlace(MatrixBase<Derived> &bAndX) const;
void compute(const MatrixType& matrix);
protected:
/** \internal
* Used to compute and store L
* The strict upper part is not used and even not initialized.
*/
MatrixType m_matrix;
bool m_isPositiveDefinite;
};
/** Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of \a matrix
*/
template<typename MatrixType>
void Cholesky<MatrixType>::compute(const MatrixType& a)
{
assert(a.rows()==a.cols());
const int size = a.rows();
m_matrix.resize(size, size);
const RealScalar eps = ei_sqrt(precision<Scalar>());
RealScalar x;
x = ei_real(a.coeff(0,0));
m_isPositiveDefinite = x > eps && ei_isMuchSmallerThan(ei_imag(a.coeff(0,0)), RealScalar(1));
m_matrix.coeffRef(0,0) = ei_sqrt(x);
m_matrix.col(0).end(size-1) = a.row(0).end(size-1).adjoint() / ei_real(m_matrix.coeff(0,0));
for (int j = 1; j < size; ++j)
{
Scalar tmp = ei_real(a.coeff(j,j)) - m_matrix.row(j).start(j).norm2();
x = ei_real(tmp);
if (x < eps || (!ei_isMuchSmallerThan(ei_imag(tmp), RealScalar(1))))
{
m_isPositiveDefinite = false;
return;
}
m_matrix.coeffRef(j,j) = x = ei_sqrt(x);
int endSize = size-j-1;
if (endSize>0) {
// Note that when all matrix columns have good alignment, then the following
// product is guaranteed to be optimal with respect to alignment.
m_matrix.col(j).end(endSize) =
(m_matrix.block(j+1, 0, endSize, j) * m_matrix.row(j).start(j).adjoint()).lazy();
// FIXME could use a.col instead of a.row
m_matrix.col(j).end(endSize) = (a.row(j).end(endSize).adjoint()
- m_matrix.col(j).end(endSize) ) / x;
}
}
}
/** \deprecated */
template<typename MatrixType>
template<typename Derived>
typename Derived::Eval Cholesky<MatrixType>::solve(const MatrixBase<Derived> &b) const
{
const int size = m_matrix.rows();
ei_assert(size==b.rows());
typename ei_eval_to_column_major<Derived>::type x(b);
solveInPlace(x);
return x;
}
/** Computes the solution x of \f$ A x = b \f$ using the current decomposition of A.
* The result is stored in \a bAndx
*
* \returns true in case of success, false otherwise.
*
* In other words, it computes \f$ b = A^{-1} b \f$ with
* \f$ {L^{*}}^{-1} L^{-1} b \f$ from right to left.
* \param bAndX stores both the matrix \f$ b \f$ and the result \f$ x \f$
*
* Example: \include Cholesky_solve.cpp
* Output: \verbinclude Cholesky_solve.out
*
* \sa MatrixBase::cholesky(), Cholesky::solveInPlace()
*/
template<typename MatrixType>
template<typename RhsDerived, typename ResDerived>
bool Cholesky<MatrixType>::solve(const MatrixBase<RhsDerived> &b, MatrixBase<ResDerived> *result) const
{
const int size = m_matrix.rows();
ei_assert(size==b.rows() && "Cholesky::solve(): invalid number of rows of the right hand side matrix b");
return solveInPlace((*result) = b);
}
/** This is the \em in-place version of solve().
*
* \param bAndX represents both the right-hand side matrix b and result x.
*
* This version avoids a copy when the right hand side matrix b is not
* needed anymore.
*
* \sa Cholesky::solve(), MatrixBase::cholesky()
*/
template<typename MatrixType>
template<typename Derived>
bool Cholesky<MatrixType>::solveInPlace(MatrixBase<Derived> &bAndX) const
{
const int size = m_matrix.rows();
ei_assert(size==bAndX.rows());
if (!m_isPositiveDefinite)
return false;
matrixL().solveTriangularInPlace(bAndX);
m_matrix.adjoint().template part<Upper>().solveTriangularInPlace(bAndX);
return true;
}
/** \cholesky_module
* \deprecated has been renamed llt()
*/
template<typename Derived>
inline const Cholesky<typename MatrixBase<Derived>::EvalType>
MatrixBase<Derived>::cholesky() const
{
return Cholesky<typename ei_eval<Derived>::type>(derived());
}
#endif // EIGEN_CHOLESKY_H
|