Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering > Class Template Reference

A direct sparse LLT Cholesky factorizations. More...

#include <src/eigen/Eigen/src/SparseCholesky/SimplicialCholesky.h>

+ Inheritance diagram for Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >:
+ Collaboration diagram for Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >:

Public Types

enum  { UpLo = _UpLo }
 
typedef _MatrixType MatrixType
 
typedef SimplicialCholeskyBase< SimplicialLLTBase
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, IndexCholMatrixType
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
typedef internal::traits< SimplicialLLTTraits
 
typedef Traits::MatrixL MatrixL
 
typedef Traits::MatrixU MatrixU
 
enum  
 
enum  
 
typedef internal::traits< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::OrderingType OrderingType
 
typedef CholMatrixType constConstCholMatrixPtr
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorI
 

Public Member Functions

 SimplicialLLT ()
 
 SimplicialLLT (const MatrixType &matrix)
 
const MatrixL matrixL () const
 
const MatrixU matrixU () const
 
SimplicialLLTcompute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &a)
 
void factorize (const MatrixType &a)
 
Scalar determinant () const
 
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived ()
 
const SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived () const
 
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived ()
 
const SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived () const
 
Index cols () const
 
Index rows () const
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationP () const
 
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & permutationPinv () const
 
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
void dumpMemory (Stream &s)
 
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) 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
 

Protected Member Functions

void compute (const MatrixType &matrix)
 
void factorize (const MatrixType &a)
 
void factorize_preordered (const CholMatrixType &a)
 
void analyzePattern (const MatrixType &a, bool doLDLT)
 
void analyzePattern_preordered (const CholMatrixType &a, bool doLDLT)
 
void ordering (const MatrixType &a, ConstCholMatrixPtr &pmat, CholMatrixType &ap)
 

Protected Attributes

ComputationInfo m_info
 
bool m_factorizationIsOk
 
bool m_analysisIsOk
 
CholMatrixType m_matrix
 
VectorType m_diag
 
VectorI m_parent
 
VectorI m_nonZerosPerCol
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_P
 
PermutationMatrix< Dynamic, Dynamic, StorageIndexm_Pinv
 
RealScalar m_shiftOffset
 
RealScalar m_shiftScale
 

Private Attributes

bool m_isInitialized
 

Detailed Description

template<typename _MatrixType, int _UpLo, typename _Ordering>
class Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >

A direct sparse LLT Cholesky factorizations.

This class provides a LL^T Cholesky factorizations of sparse matrices that are selfadjoint and positive definite. The factorization allows for solving A.X = B where X and B can be either dense or sparse.

In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization such that the factorized matrix is P A P^-1.

Template Parameters
_MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
_UpLothe triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
_OrderingThe ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>

\implsparsesolverconcept

See also
class SimplicialLDLT, class AMDOrdering, class NaturalOrdering

Member Typedef Documentation

◆ Base

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef SimplicialCholeskyBase<SimplicialLLT> Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::Base

◆ CholMatrixType

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef SparseMatrix<Scalar,ColMajor,Index> Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::CholMatrixType

◆ ConstCholMatrixPtr

typedef CholMatrixType const* Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::ConstCholMatrixPtr
inherited

◆ MatrixL

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef Traits::MatrixL Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::MatrixL

◆ MatrixType

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef _MatrixType Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::MatrixType

◆ MatrixU

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef Traits::MatrixU Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::MatrixU

◆ OrderingType

typedef internal::traits<SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::OrderingType Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::OrderingType
inherited

◆ RealScalar

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef MatrixType::RealScalar Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::RealScalar

◆ Scalar

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef MatrixType::Scalar Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::Scalar

◆ StorageIndex

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef MatrixType::StorageIndex Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::StorageIndex

◆ Traits

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef internal::traits<SimplicialLLT> Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::Traits

◆ VectorI

typedef Matrix<StorageIndex,Dynamic,1> Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::VectorI
inherited

◆ VectorType

template<typename _MatrixType , int _UpLo, typename _Ordering >
typedef Matrix<Scalar,Dynamic,1> Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::VectorType

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
inherited
63{ UpLo = internal::traits<Derived>::UpLo };

◆ anonymous enum

anonymous enum
inherited
72 {
73 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
74 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
75 };

◆ anonymous enum

template<typename _MatrixType , int _UpLo, typename _Ordering >
anonymous enum
Enumerator
UpLo 
334{ UpLo = _UpLo };
@ UpLo
Definition SimplicialCholesky.h:334

Constructor & Destructor Documentation

◆ SimplicialLLT() [1/2]

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::SimplicialLLT ( )
inline

Default constructor

346: Base() {}
SimplicialCholeskyBase< SimplicialLLT > Base
Definition SimplicialCholesky.h:335

◆ SimplicialLLT() [2/2]

template<typename _MatrixType , int _UpLo, typename _Ordering >
Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::SimplicialLLT ( const MatrixType matrix)
inlineexplicit

Constructs and performs the LLT factorization of matrix

349 : Base(matrix) {}

Member Function Documentation

◆ _solve_impl() [1/2]

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
inlineinherited
157 {
158 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
159 eigen_assert(m_matrix.rows()==b.rows());
160
161 if(m_info!=Success)
162 return;
163
164 if(m_P.size()>0)
165 dest = m_P * b;
166 else
167 dest = b;
168
169 if(m_matrix.nonZeros()>0) // otherwise L==I
170 derived().matrixL().solveInPlace(dest);
171
172 if(m_diag.size()>0)
173 dest = m_diag.asDiagonal().inverse() * dest;
174
175 if (m_matrix.nonZeros()>0) // otherwise U==I
176 derived().matrixU().solveInPlace(dest);
177
178 if(m_P.size()>0)
179 dest = m_Pinv * dest;
180 }
#define eigen_assert(x)
Definition Macros.h:579
Index size() const
Definition PermutationMatrix.h:108
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_P
Definition SimplicialCholesky.h:259
SimplicialLLT< _MatrixType, _UpLo, _Ordering > & derived()
Definition SimplicialCholesky.h:96
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition SimplicialCholesky.h:260
Index nonZeros() const
Definition SparseCompressedBase.h:56
Index rows() const
Definition SparseMatrix.h:136

◆ _solve_impl() [2/2]

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inlineinherited
184 {
185 internal::solve_sparse_through_dense_panels(derived(), b, dest);
186 }

◆ analyzePattern() [1/2]

template<typename _MatrixType , int _UpLo, typename _Ordering >
void Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::analyzePattern ( const MatrixType a)
inline

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
377 {
378 Base::analyzePattern(a, false);
379 }
void analyzePattern(const MatrixType &a, bool doLDLT)
Definition SimplicialCholesky.h:230

References Eigen::SimplicialCholeskyBase< Derived >::analyzePattern().

+ Here is the call graph for this function:

◆ analyzePattern() [2/2]

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::analyzePattern ( const MatrixType a,
bool  doLDLT 
)
inlineprotectedinherited
231 {
232 eigen_assert(a.rows()==a.cols());
233 Index size = a.cols();
234 CholMatrixType tmp(size,size);
236 ordering(a, pmat, tmp);
237 analyzePattern_preordered(*pmat,doLDLT);
238 }
CholMatrixType const * ConstCholMatrixPtr
Definition SimplicialCholesky.h:68
void analyzePattern_preordered(const CholMatrixType &a, bool doLDLT)
Definition SimplicialCholesky_impl.h:51
SparseMatrix< Scalar, ColMajor, StorageIndex > CholMatrixType
Definition SimplicialCholesky.h:67
static enum @17 ordering
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

◆ analyzePattern_preordered()

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::analyzePattern_preordered ( const CholMatrixType a,
bool  doLDLT 
)
protectedinherited
52{
53 const StorageIndex size = StorageIndex(ap.rows());
54 m_matrix.resize(size, size);
55 m_parent.resize(size);
57
58 ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
59
60 for(StorageIndex k = 0; k < size; ++k)
61 {
62 /* L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) */
63 m_parent[k] = -1; /* parent of k is not yet known */
64 tags[k] = k; /* mark node k as visited */
65 m_nonZerosPerCol[k] = 0; /* count of nonzeros in column k of L */
66 for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
67 {
68 StorageIndex i = it.index();
69 if(i < k)
70 {
71 /* follow path from i to root of etree, stop at flagged node */
72 for(; tags[i] != k; i = m_parent[i])
73 {
74 /* find parent of i if not yet determined */
75 if (m_parent[i] == -1)
76 m_parent[i] = k;
77 m_nonZerosPerCol[i]++; /* L (k,i) is nonzero */
78 tags[i] = k; /* mark i as visited */
79 }
80 }
81 }
82 }
83
84 /* construct Lp index array from m_nonZerosPerCol column counts */
86 Lp[0] = 0;
87 for(StorageIndex k = 0; k < size; ++k)
88 Lp[k+1] = Lp[k] + m_nonZerosPerCol[k] + (doLDLT ? 0 : 1);
89
90 m_matrix.resizeNonZeros(Lp[size]);
91
92 m_isInitialized = true;
94 m_analysisIsOk = true;
95 m_factorizationIsOk = false;
96}
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition Memory.h:644
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
MatrixType::StorageIndex StorageIndex
Definition SimplicialCholesky.h:66
void resizeNonZeros(Index size)
Definition SparseMatrix.h:644
const StorageIndex * outerIndexPtr() const
Definition SparseMatrix.h:166
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
@ Success
Definition Constants.h:432

◆ cols()

Index Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::cols ( ) const
inlineinherited
99{ return m_matrix.cols(); }
Index cols() const
Definition SparseMatrix.h:138

◆ compute() [1/2]

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::compute ( const MatrixType matrix)
inlineprotectedinherited

Computes the sparse Cholesky decomposition of matrix

195 {
196 eigen_assert(matrix.rows()==matrix.cols());
197 Index size = matrix.cols();
198 CholMatrixType tmp(size,size);
200 ordering(matrix, pmat, tmp);
201 analyzePattern_preordered(*pmat, DoLDLT);
202 factorize_preordered<DoLDLT>(*pmat);
203 }

◆ compute() [2/2]

template<typename _MatrixType , int _UpLo, typename _Ordering >
SimplicialLLT & Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::compute ( const MatrixType matrix)
inline

Computes the sparse Cholesky decomposition of matrix

365 {
366 Base::template compute<false>(matrix);
367 return *this;
368 }

Referenced by igl::embree::bone_heat().

+ Here is the caller graph for this function:

◆ derived() [1/4]

SimplicialLLT< _MatrixType, _UpLo, _Ordering > & Eigen::SparseSolverBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::derived ( )
inlineinherited
79{ return *static_cast<Derived*>(this); }

◆ derived() [2/4]

SimplicialLLT< _MatrixType, _UpLo, _Ordering > & Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::derived ( )
inlineinherited
96{ return *static_cast<Derived*>(this); }

◆ derived() [3/4]

const SimplicialLLT< _MatrixType, _UpLo, _Ordering > & Eigen::SparseSolverBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::derived ( ) const
inlineinherited
80{ return *static_cast<const Derived*>(this); }

◆ derived() [4/4]

const SimplicialLLT< _MatrixType, _UpLo, _Ordering > & Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::derived ( ) const
inlineinherited
97{ return *static_cast<const Derived*>(this); }

◆ determinant()

template<typename _MatrixType , int _UpLo, typename _Ordering >
Scalar Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization
394 {
395 Scalar detL = Base::m_matrix.diagonal().prod();
396 return numext::abs2(detL);
397 }
MatrixType::Scalar Scalar
Definition SimplicialCholesky.h:336
const ConstDiagonalReturnType diagonal() const
Definition SparseMatrix.h:650

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::diagonal(), and Eigen::SimplicialCholeskyBase< Derived >::m_matrix.

+ Here is the call graph for this function:

◆ dumpMemory()

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::dumpMemory ( Stream &  s)
inlineinherited
143 {
144 int total = 0;
145 s << " L: " << ((total+=(m_matrix.cols()+1) * sizeof(int) + m_matrix.nonZeros()*(sizeof(int)+sizeof(Scalar))) >> 20) << "Mb" << "\n";
146 s << " diag: " << ((total+=m_diag.size() * sizeof(Scalar)) >> 20) << "Mb" << "\n";
147 s << " tree: " << ((total+=m_parent.size() * sizeof(int)) >> 20) << "Mb" << "\n";
148 s << " nonzeros: " << ((total+=m_nonZerosPerCol.size() * sizeof(int)) >> 20) << "Mb" << "\n";
149 s << " perm: " << ((total+=m_P.size() * sizeof(int)) >> 20) << "Mb" << "\n";
150 s << " perm^-1: " << ((total+=m_Pinv.size() * sizeof(int)) >> 20) << "Mb" << "\n";
151 s << " TOTAL: " << (total>> 20) << "Mb" << "\n";
152 }
MatrixType::Scalar Scalar
Definition SimplicialCholesky.h:64

◆ factorize() [1/2]

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::factorize ( const MatrixType a)
inlineprotectedinherited
207 {
208 eigen_assert(a.rows()==a.cols());
209 Index size = a.cols();
210 CholMatrixType tmp(size,size);
212
213 if(m_P.size()==0 && (UpLo&Upper)==Upper)
214 {
215 // If there is no ordering, try to directly use the input matrix without any copy
216 internal::simplicial_cholesky_grab_input<CholMatrixType,MatrixType>::run(a, pmat, tmp);
217 }
218 else
219 {
220 tmp.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>().twistedBy(m_P);
221 pmat = &tmp;
222 }
223
224 factorize_preordered<DoLDLT>(*pmat);
225 }

◆ factorize() [2/2]

template<typename _MatrixType , int _UpLo, typename _Ordering >
void Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::factorize ( const MatrixType a)
inline

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
388 {
389 Base::template factorize<false>(a);
390 }

◆ factorize_preordered()

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::factorize_preordered ( const CholMatrixType a)
protectedinherited
102{
103 using std::sqrt;
104
105 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
106 eigen_assert(ap.rows()==ap.cols());
107 eigen_assert(m_parent.size()==ap.rows());
108 eigen_assert(m_nonZerosPerCol.size()==ap.rows());
109
110 const StorageIndex size = StorageIndex(ap.rows());
111 const StorageIndex* Lp = m_matrix.outerIndexPtr();
113 Scalar* Lx = m_matrix.valuePtr();
114
116 ei_declare_aligned_stack_constructed_variable(StorageIndex, pattern, size, 0);
117 ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
118
119 bool ok = true;
120 m_diag.resize(DoLDLT ? size : 0);
121
122 for(StorageIndex k = 0; k < size; ++k)
123 {
124 // compute nonzero pattern of kth row of L, in topological order
125 y[k] = 0.0; // Y(0:k) is now all zero
126 StorageIndex top = size; // stack for pattern is empty
127 tags[k] = k; // mark node k as visited
128 m_nonZerosPerCol[k] = 0; // count of nonzeros in column k of L
129 for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
130 {
131 StorageIndex i = it.index();
132 if(i <= k)
133 {
134 y[i] += numext::conj(it.value()); /* scatter A(i,k) into Y (sum duplicates) */
135 Index len;
136 for(len = 0; tags[i] != k; i = m_parent[i])
137 {
138 pattern[len++] = i; /* L(k,i) is nonzero */
139 tags[i] = k; /* mark i as visited */
140 }
141 while(len > 0)
142 pattern[--top] = pattern[--len];
143 }
144 }
145
146 /* compute numerical values kth row of L (a sparse triangular solve) */
147
148 RealScalar d = numext::real(y[k]) * m_shiftScale + m_shiftOffset; // get D(k,k), apply the shift function, and clear Y(k)
149 y[k] = 0.0;
150 for(; top < size; ++top)
151 {
152 Index i = pattern[top]; /* pattern[top:n-1] is pattern of L(:,k) */
153 Scalar yi = y[i]; /* get and clear Y(i) */
154 y[i] = 0.0;
155
156 /* the nonzero entry L(k,i) */
157 Scalar l_ki;
158 if(DoLDLT)
159 l_ki = yi / m_diag[i];
160 else
161 yi = l_ki = yi / Lx[Lp[i]];
162
163 Index p2 = Lp[i] + m_nonZerosPerCol[i];
164 Index p;
165 for(p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p)
166 y[Li[p]] -= numext::conj(Lx[p]) * yi;
167 d -= numext::real(l_ki * numext::conj(yi));
168 Li[p] = k; /* store L(k,i) in column form of L */
169 Lx[p] = l_ki;
170 ++m_nonZerosPerCol[i]; /* increment count of nonzeros in col i */
171 }
172 if(DoLDLT)
173 {
174 m_diag[k] = d;
175 if(d == RealScalar(0))
176 {
177 ok = false; /* failure, D(k,k) is zero */
178 break;
179 }
180 }
181 else
182 {
183 Index p = Lp[k] + m_nonZerosPerCol[k]++;
184 Li[p] = k ; /* store L(k,k) = sqrt (d) in column k */
185 if(d <= RealScalar(0)) {
186 ok = false; /* failure, matrix is not positive definite */
187 break;
188 }
189 Lx[p] = sqrt(d) ;
190 }
191 }
192
194 m_factorizationIsOk = true;
195}
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition ArrayCwiseUnaryOps.h:152
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
Definition PlainObjectBase.h:151
MatrixType::RealScalar RealScalar
Definition SimplicialCholesky.h:65
const StorageIndex * innerIndexPtr() const
Definition SparseMatrix.h:157
const Scalar * valuePtr() const
Definition SparseMatrix.h:148
@ NumericalIssue
Definition Constants.h:434
const Scalar & y
Definition MathFunctions.h:552

◆ info()

ComputationInfo Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::info ( ) const
inlineinherited

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.
108 {
109 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
110 return m_info;
111 }

◆ matrixL()

template<typename _MatrixType , int _UpLo, typename _Ordering >
const MatrixL Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::matrixL ( ) const
inline
Returns
an expression of the factor L
352 {
353 eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
354 return Traits::getL(Base::m_matrix);
355 }

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk, and Eigen::SimplicialCholeskyBase< Derived >::m_matrix.

◆ matrixU()

template<typename _MatrixType , int _UpLo, typename _Ordering >
const MatrixU Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::matrixU ( ) const
inline
Returns
an expression of the factor U (= L^*)
358 {
359 eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
360 return Traits::getU(Base::m_matrix);
361 }

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk, and Eigen::SimplicialCholeskyBase< Derived >::m_matrix.

◆ ordering()

void Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::ordering ( const MatrixType a,
ConstCholMatrixPtr pmat,
CholMatrixType ap 
)
protectedinherited
651{
652 eigen_assert(a.rows()==a.cols());
653 const Index size = a.rows();
654 pmat = &ap;
655 // Note that ordering methods compute the inverse permutation
656 if(!internal::is_same<OrderingType,NaturalOrdering<Index> >::value)
657 {
658 {
660 C = a.template selfadjointView<UpLo>();
661
663 ordering(C,m_Pinv);
664 }
665
666 if(m_Pinv.size()>0) m_P = m_Pinv.inverse();
667 else m_P.resize(0);
668
669 ap.resize(size,size);
670 ap.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>().twistedBy(m_P);
671 }
672 else
673 {
674 m_Pinv.resize(0);
675 m_P.resize(0);
676 if(int(UpLo)==int(Lower) || MatrixType::IsRowMajor)
677 {
678 // we have to transpose the lower part to to the upper one
679 ap.resize(size,size);
680 ap.template selfadjointView<Upper>() = a.template selfadjointView<UpLo>();
681 }
682 else
683 internal::simplicial_cholesky_grab_input<CholMatrixType,MatrixType>::run(a, pmat, ap);
684 }
685}
void resize(Index newSize)
Definition PermutationMatrix.h:136
InverseReturnType inverse() const
Definition PermutationMatrix.h:196
internal::traits< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::OrderingType OrderingType
Definition SimplicialCholesky.h:62

◆ permutationP()

const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::permutationP ( ) const
inlineinherited
Returns
the permutation P
See also
permutationPinv()
116 { return m_P; }

◆ permutationPinv()

const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::permutationPinv ( ) const
inlineinherited
Returns
the inverse P^-1 of the permutation P
See also
permutationP()
121 { return m_Pinv; }

◆ rows()

Index Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::rows ( ) const
inlineinherited
100{ return m_matrix.rows(); }

◆ setShift()

SimplicialLLT< _MatrixType, _UpLo, _Ordering > & Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::setShift ( const RealScalar offset,
const RealScalar scale = 1 
)
inlineinherited

Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical factorization.

During the numerical factorization, the diagonal coefficients are transformed by the following linear model:
d_ii = offset + scale * d_ii

The default is the identity transformation with offset=0, and scale=1.

Returns
a reference to *this.
133 {
136 return derived();
137 }
int scale(const int val)
Definition WipeTowerDialog.cpp:14
void offset(Slic3r::ExPolygon &sh, coord_t distance, const PolygonTag &)
Definition geometries.hpp:132

◆ 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
Derived & derived()
Definition SparseSolverBase.h:79
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 }

Member Data Documentation

◆ m_analysisIsOk

bool Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_analysisIsOk
protectedinherited

◆ m_diag

VectorType Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_diag
protectedinherited

◆ m_factorizationIsOk

bool Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_factorizationIsOk
protectedinherited

◆ m_info

ComputationInfo Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_info
mutableprotectedinherited

◆ m_isInitialized

bool Eigen::SparseSolverBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_isInitialized
mutableprivateinherited

◆ m_matrix

CholMatrixType Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_matrix
protectedinherited

◆ m_nonZerosPerCol

VectorI Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_nonZerosPerCol
protectedinherited

◆ m_P

PermutationMatrix<Dynamic,Dynamic,StorageIndex> Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_P
protectedinherited

◆ m_parent

VectorI Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_parent
protectedinherited

◆ m_Pinv

PermutationMatrix<Dynamic,Dynamic,StorageIndex> Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_Pinv
protectedinherited

◆ m_shiftOffset

RealScalar Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_shiftOffset
protectedinherited

◆ m_shiftScale

RealScalar Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >::m_shiftScale
protectedinherited

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