00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #ifndef EIGEN_SPARSELLT_H
00026 #define EIGEN_SPARSELLT_H
00027
00038 template<typename MatrixType, int Backend = DefaultBackend>
00039 class SparseLLT
00040 {
00041 protected:
00042 typedef typename MatrixType::Scalar Scalar;
00043 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00044 typedef SparseMatrix<Scalar,LowerTriangular> CholMatrixType;
00045
00046 enum {
00047 SupernodalFactorIsDirty = 0x10000,
00048 MatrixLIsDirty = 0x20000
00049 };
00050
00051 public:
00052
00054 SparseLLT(int flags = 0)
00055 : m_flags(flags), m_status(0)
00056 {
00057 m_precision = RealScalar(0.1) * Eigen::precision<RealScalar>();
00058 }
00059
00062 SparseLLT(const MatrixType& matrix, int flags = 0)
00063 : m_matrix(matrix.rows(), matrix.cols()), m_flags(flags), m_status(0)
00064 {
00065 m_precision = RealScalar(0.1) * Eigen::precision<RealScalar>();
00066 compute(matrix);
00067 }
00068
00082 void setPrecision(RealScalar v) { m_precision = v; }
00083
00087 RealScalar precision() const { return m_precision; }
00088
00099 void setFlags(int f) { m_flags = f; }
00101 int flags() const { return m_flags; }
00102
00104 void compute(const MatrixType& matrix);
00105
00107 inline const CholMatrixType& matrixL(void) const { return m_matrix; }
00108
00109 template<typename Derived>
00110 bool solveInPlace(MatrixBase<Derived> &b) const;
00111
00113 inline bool succeeded(void) const { return m_succeeded; }
00114
00115 protected:
00116 CholMatrixType m_matrix;
00117 RealScalar m_precision;
00118 int m_flags;
00119 mutable int m_status;
00120 bool m_succeeded;
00121 };
00122
00126 template<typename MatrixType, int Backend>
00127 void SparseLLT<MatrixType,Backend>::compute(const MatrixType& a)
00128 {
00129 assert(a.rows()==a.cols());
00130 const int size = a.rows();
00131 m_matrix.resize(size, size);
00132
00133
00134 AmbiVector<Scalar> tempVector(size);
00135 RealScalar density = a.nonZeros()/RealScalar(size*size);
00136
00137
00138 m_matrix.startFill(a.nonZeros()*2);
00139 for (int j = 0; j < size; ++j)
00140 {
00141 Scalar x = ei_real(a.coeff(j,j));
00142
00143
00144 tempVector.init(density>0.001? IsDense : IsSparse);
00145 tempVector.setBounds(j+1,size);
00146 tempVector.setZero();
00147
00148 {
00149 typename MatrixType::InnerIterator it(a,j);
00150 ++it;
00151 for (; it; ++it)
00152 tempVector.coeffRef(it.index()) = it.value();
00153 }
00154 for (int k=0; k<j+1; ++k)
00155 {
00156 typename CholMatrixType::InnerIterator it(m_matrix, k);
00157 while (it && it.index()<j)
00158 ++it;
00159 if (it && it.index()==j)
00160 {
00161 Scalar y = it.value();
00162 x -= ei_abs2(y);
00163 ++it;
00164 tempVector.restart();
00165 for (; it; ++it)
00166 {
00167 tempVector.coeffRef(it.index()) -= it.value() * y;
00168 }
00169 }
00170 }
00171
00172
00173 RealScalar rx = ei_sqrt(ei_real(x));
00174 m_matrix.fill(j,j) = rx;
00175 Scalar y = Scalar(1)/rx;
00176 for (typename AmbiVector<Scalar>::Iterator it(tempVector, m_precision*rx); it; ++it)
00177 {
00178 m_matrix.fill(it.index(), j) = it.value() * y;
00179 }
00180 }
00181 m_matrix.endFill();
00182 }
00183
00185 template<typename MatrixType, int Backend>
00186 template<typename Derived>
00187 bool SparseLLT<MatrixType, Backend>::solveInPlace(MatrixBase<Derived> &b) const
00188 {
00189 const int size = m_matrix.rows();
00190 ei_assert(size==b.rows());
00191
00192 m_matrix.solveTriangularInPlace(b);
00193
00194 if (NumTraits<Scalar>::IsComplex)
00195 {
00196 CholMatrixType aux = m_matrix.conjugate();
00197 aux.transpose().solveTriangularInPlace(b);
00198 }
00199 else
00200 m_matrix.transpose().solveTriangularInPlace(b);
00201
00202 return true;
00203 }
00204
00205 #endif // EIGEN_SPARSELLT_H