Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType > Class Template Reference

Modified Incomplete Cholesky with dual threshold. More...

#include <src/eigen/Eigen/src/IterativeLinearSolvers/IncompleteCholesky.h>

+ Inheritance diagram for Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >:
+ Collaboration diagram for Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >:

Public Types

enum  { UpLo = _UpLo }
 
enum  { ColsAtCompileTime = Dynamic , MaxColsAtCompileTime = Dynamic }
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef _OrderingType OrderingType
 
typedef OrderingType::PermutationType PermutationType
 
typedef PermutationType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexFactorType
 
typedef Matrix< Scalar, Dynamic, 1 > VectorSx
 
typedef Matrix< RealScalar, Dynamic, 1 > VectorRx
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorIx
 
typedef std::vector< std::list< StorageIndex > > VectorList
 

Public Member Functions

 IncompleteCholesky ()
 
template<typename MatrixType >
 IncompleteCholesky (const MatrixType &matrix)
 
Index rows () const
 
Index cols () const
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
void setInitialShift (RealScalar shift)
 Set the initial shift parameter $ \sigma $.
 
template<typename MatrixType >
void analyzePattern (const MatrixType &mat)
 Computes the fill reducing permutation vector using the sparsity pattern of mat.
 
template<typename MatrixType >
void factorize (const MatrixType &mat)
 Performs the numerical factorization of the input matrix mat.
 
template<typename MatrixType >
void compute (const MatrixType &mat)
 
template<typename Rhs , typename Dest >
void _solve_impl (const Rhs &b, Dest &x) const
 
const FactorTypematrixL () const
 
const VectorRxscalingS () const
 
const PermutationTypepermutationP () const
 
template<typename _MatrixType >
void factorize (const _MatrixType &mat)
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< IncompleteCholesky< Scalar, _UpLo, _OrderingType > > Base
 

Protected Attributes

FactorType m_L
 
VectorRx m_scale
 
RealScalar m_initialShift
 
bool m_analysisIsOk
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
PermutationType m_perm
 
bool m_isInitialized
 

Private Member Functions

void updateList (Ref< const VectorIx > colPtr, Ref< VectorIx > rowIdx, Ref< VectorSx > vals, const Index &col, const Index &jk, VectorIx &firstElt, VectorList &listCol)
 

Detailed Description

template<typename Scalar, int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
class Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >

Modified Incomplete Cholesky with dual threshold.

References : C-J. Lin and J. J. Moré, Incomplete Cholesky Factorizations with Limited memory, SIAM J. Sci. Comput. 21(1), pp. 24-45, 1999

Template Parameters
Scalarthe scalar type of the input matrices
_UpLoThe triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
_OrderingTypeThe ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<int>, unless EIGEN_MPL2_ONLY is defined, in which case the default is NaturalOrdering<int>.

\implsparsesolverconcept

It performs the following incomplete factorization: $ S P A P' S \approx L L' $ where L is a lower triangular factor, S is a diagonal scaling matrix, and P is a fill-in reducing permutation as computed by the ordering method.

Shifting strategy: Let $ B = S P A P' S $ be the scaled matrix on which the factorization is carried out, and $ \beta $ be the minimum value of the diagonal. If $ \beta > 0 $ then, the factorization is directly performed on the matrix B. Otherwise, the factorization is performed on the shifted matrix $ B + (\sigma+|\beta| I $ where $ \sigma $ is the initial shift value as returned and set by setInitialShift() method. The default value is $ \sigma = 10^{-3} $. If the factorization fails, then the shift in doubled until it succeed or a maximum of ten attempts. If it still fails, as returned by the info() method, then you can either increase the initial shift, or better use another preconditioning technique.

Member Typedef Documentation

◆ Base

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef SparseSolverBase<IncompleteCholesky<Scalar,_UpLo,_OrderingType> > Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::Base
protected

◆ FactorType

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::FactorType

◆ OrderingType

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef _OrderingType Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::OrderingType

◆ PermutationType

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef OrderingType::PermutationType Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::PermutationType

◆ RealScalar

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef NumTraits<Scalar>::Real Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::RealScalar

◆ StorageIndex

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef PermutationType::StorageIndex Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::StorageIndex

◆ VectorIx

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef Matrix<StorageIndex,Dynamic, 1> Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::VectorIx

◆ VectorList

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef std::vector<std::list<StorageIndex> > Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::VectorList

◆ VectorRx

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef Matrix<RealScalar,Dynamic,1> Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::VectorRx

◆ VectorSx

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
typedef Matrix<Scalar,Dynamic,1> Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::VectorSx

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
anonymous enum
Enumerator
UpLo 
66{ UpLo = _UpLo };
@ UpLo
Definition IncompleteCholesky.h:66

◆ anonymous enum

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
67 {
70 };
@ MaxColsAtCompileTime
Definition IncompleteCholesky.h:69
@ ColsAtCompileTime
Definition IncompleteCholesky.h:68
const int Dynamic
Definition Constants.h:21

Constructor & Destructor Documentation

◆ IncompleteCholesky() [1/2]

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::IncompleteCholesky ( )
inline

Default constructor leaving the object in a partly non-initialized stage.

You must call compute() or the pair analyzePattern()/factorize() to make it valid.

See also
IncompleteCholesky(const MatrixType&)
79: m_initialShift(1e-3),m_factorizationIsOk(false) {}
RealScalar m_initialShift
Definition IncompleteCholesky.h:180
bool m_factorizationIsOk
Definition IncompleteCholesky.h:182

◆ IncompleteCholesky() [2/2]

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename MatrixType >
Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::IncompleteCholesky ( const MatrixType &  matrix)
inline

Constructor computing the incomplete factorization for the given matrix matrix.

85 {
86 compute(matrix);
87 }
void compute(const MatrixType &mat)
Definition IncompleteCholesky.h:147

References Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::compute().

+ Here is the call graph for this function:

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename Rhs , typename Dest >
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::_solve_impl ( const Rhs &  b,
Dest &  x 
) const
inline
156 {
157 eigen_assert(m_factorizationIsOk && "factorize() should be called first");
158 if (m_perm.rows() == b.rows()) x = m_perm * b;
159 else x = b;
160 x = m_scale.asDiagonal() * x;
161 x = m_L.template triangularView<Lower>().solve(x);
162 x = m_L.adjoint().template triangularView<Upper>().solve(x);
163 x = m_scale.asDiagonal() * x;
164 if (m_perm.rows() == b.rows())
165 x = m_perm.inverse() * x;
166 }
#define eigen_assert(x)
Definition Macros.h:579
VectorRx m_scale
Definition IncompleteCholesky.h:179
FactorType m_L
Definition IncompleteCholesky.h:178
PermutationType m_perm
Definition IncompleteCholesky.h:184
const AdjointReturnType adjoint() const
Definition SparseMatrixBase.h:351
TCoord< P > x(const P &p)
Definition geometry_traits.hpp:297

References Eigen::SparseMatrixBase< Derived >::adjoint(), eigen_assert, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_factorizationIsOk, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_L, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_perm, and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_scale.

+ Here is the call graph for this function:

◆ _solve_impl() [2/2]

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inlineinherited
112 {
113 internal::solve_sparse_through_dense_panels(derived(), b.derived(), dest.derived());
114 }
Derived & derived()
Definition SparseSolverBase.h:79
enable_if< Rhs::ColsAtCompileTime!=1 &&Dest::ColsAtCompileTime!=1 >::type solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition SparseSolverBase.h:23

References Eigen::SparseSolverBase< Derived >::derived(), Eigen::SparseMatrixBase< Derived >::derived(), and Eigen::internal::solve_sparse_through_dense_panels().

+ Here is the call graph for this function:

◆ analyzePattern()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename MatrixType >
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::analyzePattern ( const MatrixType &  mat)
inline

Computes the fill reducing permutation vector using the sparsity pattern of mat.

118 {
119 OrderingType ord;
121 ord(mat.template selfadjointView<UpLo>(), pinv);
122 if(pinv.size()>0) m_perm = pinv.inverse();
123 else m_perm.resize(0);
124 m_L.resize(mat.rows(), mat.cols());
125 m_analysisIsOk = true;
126 m_isInitialized = true;
127 m_info = Success;
128 }
bool m_analysisIsOk
Definition IncompleteCholesky.h:181
OrderingType::PermutationType PermutationType
Definition IncompleteCholesky.h:59
ComputationInfo m_info
Definition IncompleteCholesky.h:183
bool m_isInitialized
Definition SparseSolverBase.h:119
_OrderingType OrderingType
Definition IncompleteCholesky.h:58
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
@ Success
Definition Constants.h:432
void pinv(const Eigen::MatrixBase< DerivedA > &A, typename DerivedA::Scalar tol, Eigen::PlainObjectBase< DerivedX > &X)
Definition pinv.cpp:6

References Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_analysisIsOk, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_info, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_isInitialized, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_L, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_perm, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::resize(), and Eigen::Success.

Referenced by Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::compute().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ cols()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
Index Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::cols ( ) const
inline
Returns
number of columns of the factored matrix
93{ return m_L.cols(); }
Index cols() const
Definition SparseMatrix.h:138

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::cols(), and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_L.

+ Here is the call graph for this function:

◆ compute()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename MatrixType >
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::compute ( const MatrixType &  mat)
inline

Computes or re-computes the incomplete Cholesky factorization of the input matrix mat

It is a shortcut for a sequential call to the analyzePattern() and factorize() methods.

See also
analyzePattern(), factorize()
148 {
149 analyzePattern(mat);
150 factorize(mat);
151 }
void analyzePattern(const MatrixType &mat)
Computes the fill reducing permutation vector using the sparsity pattern of mat.
Definition IncompleteCholesky.h:117
void factorize(const MatrixType &mat)
Performs the numerical factorization of the input matrix mat.

References Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::analyzePattern(), and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::factorize().

Referenced by Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::IncompleteCholesky().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ derived() [1/2]

template<typename Derived >
Derived & Eigen::SparseSolverBase< Derived >::derived ( )
inlineinherited
79{ return *static_cast<Derived*>(this); }

Referenced by Eigen::SparseSolverBase< Derived >::_solve_impl(), and Eigen::SparseSolverBase< Derived >::solve().

+ Here is the caller graph for this function:

◆ derived() [2/2]

template<typename Derived >
const Derived & Eigen::SparseSolverBase< Derived >::derived ( ) const
inlineinherited
80{ return *static_cast<const Derived*>(this); }

◆ factorize() [1/2]

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename _MatrixType >
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::factorize ( const _MatrixType &  mat)
197{
198 using std::sqrt;
199 eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
200
201 // Dropping strategy : Keep only the p largest elements per column, where p is the number of elements in the column of the original matrix. Other strategies will be added
202
203 // Apply the fill-reducing permutation computed in analyzePattern()
204 if (m_perm.rows() == mat.rows() ) // To detect the null permutation
205 {
206 // The temporary is needed to make sure that the diagonal entry is properly sorted
207 FactorType tmp(mat.rows(), mat.cols());
208 tmp = mat.template selfadjointView<_UpLo>().twistedBy(m_perm);
209 m_L.template selfadjointView<Lower>() = tmp.template selfadjointView<Lower>();
210 }
211 else
212 {
213 m_L.template selfadjointView<Lower>() = mat.template selfadjointView<_UpLo>();
214 }
215
216 Index n = m_L.cols();
217 Index nnz = m_L.nonZeros();
218 Map<VectorSx> vals(m_L.valuePtr(), nnz); //values
219 Map<VectorIx> rowIdx(m_L.innerIndexPtr(), nnz); //Row indices
220 Map<VectorIx> colPtr( m_L.outerIndexPtr(), n+1); // Pointer to the beginning of each row
221 VectorIx firstElt(n-1); // for each j, points to the next entry in vals that will be used in the factorization
222 VectorList listCol(n); // listCol(j) is a linked list of columns to update column j
223 VectorSx col_vals(n); // Store a nonzero values in each column
224 VectorIx col_irow(n); // Row indices of nonzero elements in each column
225 VectorIx col_pattern(n);
226 col_pattern.fill(-1);
227 StorageIndex col_nnz;
228
229
230 // Computes the scaling factors
231 m_scale.resize(n);
233 for (Index j = 0; j < n; j++)
234 for (Index k = colPtr[j]; k < colPtr[j+1]; k++)
235 {
236 m_scale(j) += numext::abs2(vals(k));
237 if(rowIdx[k]!=j)
238 m_scale(rowIdx[k]) += numext::abs2(vals(k));
239 }
240
241 m_scale = m_scale.cwiseSqrt().cwiseSqrt();
242
243 for (Index j = 0; j < n; ++j)
244 if(m_scale(j)>(std::numeric_limits<RealScalar>::min)())
245 m_scale(j) = RealScalar(1)/m_scale(j);
246 else
247 m_scale(j) = 1;
248
249 // TODO disable scaling if not needed, i.e., if it is roughly uniform? (this will make solve() faster)
250
251 // Scale and compute the shift for the matrix
252 RealScalar mindiag = NumTraits<RealScalar>::highest();
253 for (Index j = 0; j < n; j++)
254 {
255 for (Index k = colPtr[j]; k < colPtr[j+1]; k++)
256 vals[k] *= (m_scale(j)*m_scale(rowIdx[k]));
257 eigen_internal_assert(rowIdx[colPtr[j]]==j && "IncompleteCholesky: only the lower triangular part must be stored");
258 mindiag = numext::mini(numext::real(vals[colPtr[j]]), mindiag);
259 }
260
261 FactorType L_save = m_L;
262
263 RealScalar shift = 0;
264 if(mindiag <= RealScalar(0.))
265 shift = m_initialShift - mindiag;
266
268
269 // Try to perform the incomplete factorization using the current shift
270 int iter = 0;
271 do
272 {
273 // Apply the shift to the diagonal elements of the matrix
274 for (Index j = 0; j < n; j++)
275 vals[colPtr[j]] += shift;
276
277 // jki version of the Cholesky factorization
278 Index j=0;
279 for (; j < n; ++j)
280 {
281 // Left-looking factorization of the j-th column
282 // First, load the j-th column into col_vals
283 Scalar diag = vals[colPtr[j]]; // It is assumed that only the lower part is stored
284 col_nnz = 0;
285 for (Index i = colPtr[j] + 1; i < colPtr[j+1]; i++)
286 {
287 StorageIndex l = rowIdx[i];
288 col_vals(col_nnz) = vals[i];
289 col_irow(col_nnz) = l;
290 col_pattern(l) = col_nnz;
291 col_nnz++;
292 }
293 {
294 typename std::list<StorageIndex>::iterator k;
295 // Browse all previous columns that will update column j
296 for(k = listCol[j].begin(); k != listCol[j].end(); k++)
297 {
298 Index jk = firstElt(*k); // First element to use in the column
299 eigen_internal_assert(rowIdx[jk]==j);
300 Scalar v_j_jk = numext::conj(vals[jk]);
301
302 jk += 1;
303 for (Index i = jk; i < colPtr[*k+1]; i++)
304 {
305 StorageIndex l = rowIdx[i];
306 if(col_pattern[l]<0)
307 {
308 col_vals(col_nnz) = vals[i] * v_j_jk;
309 col_irow[col_nnz] = l;
310 col_pattern(l) = col_nnz;
311 col_nnz++;
312 }
313 else
314 col_vals(col_pattern[l]) -= vals[i] * v_j_jk;
315 }
316 updateList(colPtr,rowIdx,vals, *k, jk, firstElt, listCol);
317 }
318 }
319
320 // Scale the current column
321 if(numext::real(diag) <= 0)
322 {
323 if(++iter>=10)
324 return;
325
326 // increase shift
327 shift = numext::maxi(m_initialShift,RealScalar(2)*shift);
328 // restore m_L, col_pattern, and listCol
329 vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
330 rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
331 colPtr = Map<const VectorIx>(L_save.outerIndexPtr(), n+1);
332 col_pattern.fill(-1);
333 for(Index i=0; i<n; ++i)
334 listCol[i].clear();
335
336 break;
337 }
338
339 RealScalar rdiag = sqrt(numext::real(diag));
340 vals[colPtr[j]] = rdiag;
341 for (Index k = 0; k<col_nnz; ++k)
342 {
343 Index i = col_irow[k];
344 //Scale
345 col_vals(k) /= rdiag;
346 //Update the remaining diagonals with col_vals
347 vals[colPtr[i]] -= numext::abs2(col_vals(k));
348 }
349 // Select the largest p elements
350 // p is the original number of elements in the column (without the diagonal)
351 Index p = colPtr[j+1] - colPtr[j] - 1 ;
352 Ref<VectorSx> cvals = col_vals.head(col_nnz);
353 Ref<VectorIx> cirow = col_irow.head(col_nnz);
354 internal::QuickSplit(cvals,cirow, p);
355 // Insert the largest p elements in the matrix
356 Index cpt = 0;
357 for (Index i = colPtr[j]+1; i < colPtr[j+1]; i++)
358 {
359 vals[i] = col_vals(cpt);
360 rowIdx[i] = col_irow(cpt);
361 // restore col_pattern:
362 col_pattern(col_irow(cpt)) = -1;
363 cpt++;
364 }
365 // Get the first smallest row index and put it after the diagonal element
366 Index jk = colPtr(j)+1;
367 updateList(colPtr,rowIdx,vals,j,jk,firstElt,listCol);
368 }
369
370 if(j==n)
371 {
372 m_factorizationIsOk = true;
373 m_info = Success;
374 }
375 } while(m_info!=Success);
376}
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition ArrayCwiseUnaryOps.h:152
#define eigen_internal_assert(x)
Definition Macros.h:585
PermutationType::StorageIndex StorageIndex
Definition IncompleteCholesky.h:60
Matrix< StorageIndex, Dynamic, 1 > VectorIx
Definition IncompleteCholesky.h:64
SparseMatrix< Scalar, ColMajor, StorageIndex > FactorType
Definition IncompleteCholesky.h:61
Matrix< Scalar, Dynamic, 1 > VectorSx
Definition IncompleteCholesky.h:62
std::vector< std::list< StorageIndex > > VectorList
Definition IncompleteCholesky.h:65
NumTraits< Scalar >::Real RealScalar
Definition IncompleteCholesky.h:57
void updateList(Ref< const VectorIx > colPtr, Ref< VectorIx > rowIdx, Ref< VectorSx > vals, const Index &col, const Index &jk, VectorIx &firstElt, VectorList &listCol)
Definition IncompleteCholesky.h:379
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition CwiseNullaryOp.h:515
Index nonZeros() const
Definition SparseCompressedBase.h:56
const StorageIndex * innerIndexPtr() const
Definition SparseMatrix.h:157
const StorageIndex * outerIndexPtr() const
Definition SparseMatrix.h:166
const Scalar * valuePtr() const
Definition SparseMatrix.h:148
@ NumericalIssue
Definition Constants.h:434
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition IncompleteLUT.h:29
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition MathFunctions.h:823
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition MathFunctions.h:815
T * begin(Slic3r::Mat< N, M, T > &mat)
Definition Point.hpp:605
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33
typename Traits< remove_cvref_t< L > >::Scalar Scalar
Definition Line.hpp:36
IGL_INLINE void diag(const Eigen::SparseMatrix< T > &X, Eigen::SparseVector< T > &V)
Definition diag.cpp:17

References Eigen::begin(), eigen_assert, eigen_internal_assert, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::innerIndexPtr(), Eigen::numext::maxi(), Eigen::numext::mini(), Eigen::NumericalIssue, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::outerIndexPtr(), Eigen::internal::QuickSplit(), Eigen::PlainObjectBase< Derived >::resize(), sqrt(), Eigen::Success, Eigen::SparseMatrixBase< Derived >::twistedBy(), and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::valuePtr().

+ Here is the call graph for this function:

◆ factorize() [2/2]

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
template<typename MatrixType >
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::factorize ( const MatrixType &  mat)

Performs the numerical factorization of the input matrix mat.

The method analyzePattern() or compute() must have been called beforehand with a matrix having the same pattern.

See also
compute(), analyzePattern()

Referenced by Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::compute().

+ Here is the caller graph for this function:

◆ info()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
ComputationInfo Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::info ( ) const
inline

Reports whether previous computation was successful.

It triggers an assertion if *this has not been initialized through the respective constructor, or a call to compute() or analyzePattern().

Returns
Success if computation was successful, NumericalIssue if the matrix appears to be negative.
105 {
106 eigen_assert(m_isInitialized && "IncompleteCholesky is not initialized.");
107 return m_info;
108 }

References eigen_assert, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_info, and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_isInitialized.

◆ matrixL()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
const FactorType & Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::matrixL ( ) const
inline
Returns
the sparse lower triangular factor L
169{ eigen_assert("m_factorizationIsOk"); return m_L; }

References eigen_assert, and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_L.

◆ permutationP()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
const PermutationType & Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::permutationP ( ) const
inline
Returns
the fill-in reducing permutation P (can be empty for a natural ordering)
175{ eigen_assert("m_analysisIsOk"); return m_perm; }

References eigen_assert, and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_perm.

◆ rows()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
Index Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::rows ( ) const
inline
Returns
number of rows of the factored matrix
90{ return m_L.rows(); }
Index rows() const
Definition SparseMatrix.h:136

References Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_L, and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows().

+ Here is the call graph for this function:

◆ scalingS()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
const VectorRx & Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::scalingS ( ) const
inline
Returns
a vector representing the scaling factor S
172{ eigen_assert("m_factorizationIsOk"); return m_scale; }

References eigen_assert, and Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_scale.

◆ setInitialShift()

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
void Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::setInitialShift ( RealScalar  shift)
inline

Set the initial shift parameter $ \sigma $.

112{ m_initialShift = shift; }

References Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_initialShift.

◆ solve() [1/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const MatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
89 {
90 eigen_assert(m_isInitialized && "Solver is not initialized.");
91 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
92 return Solve<Derived, Rhs>(derived(), b.derived());
93 }
bool m_isInitialized
Definition SparseSolverBase.h:119
size_t rows(const T &raster)
Definition MarchingSquares.hpp:55

References Eigen::SparseSolverBase< Derived >::derived(), eigen_assert, and Eigen::SparseSolverBase< Derived >::m_isInitialized.

Referenced by igl::embree::bone_heat(), igl::Frame_field_deformer::compute_optimal_positions(), igl::eigs(), igl::copyleft::comiso::FrameInterpolator::interpolateSymmetric(), igl::slim::solve_weighted_arap(), and igl::copyleft::comiso::NRosyField::solveNoRoundings().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [2/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const SparseMatrixBase< Rhs > &  b) const
inlineinherited
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
102 {
103 eigen_assert(m_isInitialized && "Solver is not initialized.");
104 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
105 return Solve<Derived, Rhs>(derived(), b.derived());
106 }

◆ updateList()

template<typename Scalar , int _UpLo, typename OrderingType >
void Eigen::IncompleteCholesky< Scalar, _UpLo, OrderingType >::updateList ( Ref< const VectorIx colPtr,
Ref< VectorIx rowIdx,
Ref< VectorSx vals,
const Index col,
const Index jk,
VectorIx firstElt,
VectorList listCol 
)
inlineprivate
380{
381 if (jk < colPtr(col+1) )
382 {
383 Index p = colPtr(col+1) - jk;
384 Index minpos;
385 rowIdx.segment(jk,p).minCoeff(&minpos);
386 minpos += jk;
387 if (rowIdx(minpos) != rowIdx(jk))
388 {
389 //Swap
390 std::swap(rowIdx(jk),rowIdx(minpos));
391 std::swap(vals(jk),vals(minpos));
392 }
393 firstElt(col) = internal::convert_index<StorageIndex,Index>(jk);
394 listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex,Index>(col));
395 }
396}
EIGEN_DEVICE_FUNC ColXpr col(Index i)
This is the const version of col().
Definition BlockMethods.h:838

References col().

+ Here is the call graph for this function:

Member Data Documentation

◆ m_analysisIsOk

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
bool Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_analysisIsOk
protected

◆ m_factorizationIsOk

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
bool Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_factorizationIsOk
protected

◆ m_info

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
ComputationInfo Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_info
protected

◆ m_initialShift

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
RealScalar Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_initialShift
protected

◆ m_isInitialized

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_L

◆ m_perm

◆ m_scale

template<typename Scalar , int _UpLo = Lower, typename _OrderingType = AMDOrdering<int>>
VectorRx Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >::m_scale
protected

The documentation for this class was generated from the following file: