Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::SimplicialCholeskyBase< Derived > Class Template Reference

A base class for direct sparse Cholesky factorizations. More...

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

+ Inheritance diagram for Eigen::SimplicialCholeskyBase< Derived >:
+ Collaboration diagram for Eigen::SimplicialCholeskyBase< Derived >:

Classes

struct  keep_diag
 

Public Types

enum  { UpLo = internal::traits<Derived>::UpLo }
 
enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef internal::traits< Derived >::MatrixType MatrixType
 
typedef internal::traits< Derived >::OrderingType OrderingType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexCholMatrixType
 
typedef CholMatrixType constConstCholMatrixPtr
 
typedef Matrix< Scalar, Dynamic, 1 > VectorType
 
typedef Matrix< StorageIndex, Dynamic, 1 > VectorI
 

Public Member Functions

 SimplicialCholeskyBase ()
 
 SimplicialCholeskyBase (const MatrixType &matrix)
 
 ~SimplicialCholeskyBase ()
 
Derived & derived ()
 
const Derived & 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
 
Derived & setShift (const RealScalar &offset, const RealScalar &scale=1)
 
template<typename Stream >
void dumpMemory (Stream &s)
 
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
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
 

Protected Member Functions

template<bool DoLDLT>
void compute (const MatrixType &matrix)
 
template<bool DoLDLT>
void factorize (const MatrixType &a)
 
template<bool DoLDLT>
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 Types

typedef SparseSolverBase< Derived > Base
 

Private Attributes

bool m_isInitialized
 

Detailed Description

template<typename Derived>
class Eigen::SimplicialCholeskyBase< Derived >

A base class for direct sparse Cholesky factorizations.

This is a base class for LL^T and LDL^T Cholesky factorizations of sparse matrices that are selfadjoint and positive definite. These factorizations allow 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
Derivedthe type of the derived class, that is the actual factorization type.

Member Typedef Documentation

◆ Base

template<typename Derived >
typedef SparseSolverBase<Derived> Eigen::SimplicialCholeskyBase< Derived >::Base
private

◆ CholMatrixType

template<typename Derived >
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SimplicialCholeskyBase< Derived >::CholMatrixType

◆ ConstCholMatrixPtr

template<typename Derived >
typedef CholMatrixType const* Eigen::SimplicialCholeskyBase< Derived >::ConstCholMatrixPtr

◆ MatrixType

template<typename Derived >
typedef internal::traits<Derived>::MatrixType Eigen::SimplicialCholeskyBase< Derived >::MatrixType

◆ OrderingType

template<typename Derived >
typedef internal::traits<Derived>::OrderingType Eigen::SimplicialCholeskyBase< Derived >::OrderingType

◆ RealScalar

template<typename Derived >
typedef MatrixType::RealScalar Eigen::SimplicialCholeskyBase< Derived >::RealScalar

◆ Scalar

template<typename Derived >
typedef MatrixType::Scalar Eigen::SimplicialCholeskyBase< Derived >::Scalar

◆ StorageIndex

template<typename Derived >
typedef MatrixType::StorageIndex Eigen::SimplicialCholeskyBase< Derived >::StorageIndex

◆ VectorI

template<typename Derived >
typedef Matrix<StorageIndex,Dynamic,1> Eigen::SimplicialCholeskyBase< Derived >::VectorI

◆ VectorType

template<typename Derived >
typedef Matrix<Scalar,Dynamic,1> Eigen::SimplicialCholeskyBase< Derived >::VectorType

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
UpLo 
63{ UpLo = internal::traits<Derived>::UpLo };
@ UpLo
Definition SimplicialCholesky.h:63

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
72 {
73 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
74 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
75 };
@ ColsAtCompileTime
Definition SimplicialCholesky.h:73
@ MaxColsAtCompileTime
Definition SimplicialCholesky.h:74

Constructor & Destructor Documentation

◆ SimplicialCholeskyBase() [1/2]

template<typename Derived >
Eigen::SimplicialCholeskyBase< Derived >::SimplicialCholeskyBase ( )
inline

Default constructor

84 {}
RealScalar m_shiftOffset
Definition SimplicialCholesky.h:262
RealScalar m_shiftScale
Definition SimplicialCholesky.h:263
ComputationInfo m_info
Definition SimplicialCholesky.h:251
@ Success
Definition Constants.h:432

◆ SimplicialCholeskyBase() [2/2]

template<typename Derived >
Eigen::SimplicialCholeskyBase< Derived >::SimplicialCholeskyBase ( const MatrixType matrix)
inlineexplicit
88 {
89 derived().compute(matrix);
90 }
Derived & derived()
Definition SimplicialCholesky.h:96

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

+ Here is the call graph for this function:

◆ ~SimplicialCholeskyBase()

template<typename Derived >
Eigen::SimplicialCholeskyBase< Derived >::~SimplicialCholeskyBase ( )
inline
93 {
94 }

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SimplicialCholeskyBase< Derived >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
inline
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
CholMatrixType m_matrix
Definition SimplicialCholesky.h:255
VectorType m_diag
Definition SimplicialCholesky.h:256
PermutationMatrix< Dynamic, Dynamic, StorageIndex > m_Pinv
Definition SimplicialCholesky.h:260
bool m_factorizationIsOk
Definition SimplicialCholesky.h:252
Index nonZeros() const
Definition SparseCompressedBase.h:56
Index rows() const
Definition SparseMatrix.h:136

References Eigen::SimplicialCholeskyBase< Derived >::derived(), eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_diag, Eigen::SimplicialCholeskyBase< Derived >::m_factorizationIsOk, Eigen::SimplicialCholeskyBase< Derived >::m_info, Eigen::SimplicialCholeskyBase< Derived >::m_matrix, Eigen::SimplicialCholeskyBase< Derived >::m_P, Eigen::SimplicialCholeskyBase< Derived >::m_Pinv, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::nonZeros(), Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows(), Eigen::PermutationBase< Derived >::size(), and Eigen::Success.

+ Here is the call graph for this function:

◆ _solve_impl() [2/2]

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SimplicialCholeskyBase< Derived >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inline
184 {
186 }
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::SimplicialCholeskyBase< Derived >::derived(), and Eigen::internal::solve_sparse_through_dense_panels().

+ Here is the call graph for this function:

◆ analyzePattern()

template<typename Derived >
void Eigen::SimplicialCholeskyBase< Derived >::analyzePattern ( const MatrixType a,
bool  doLDLT 
)
inlineprotected
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

References Eigen::SimplicialCholeskyBase< Derived >::analyzePattern_preordered(), eigen_assert, and ordering.

Referenced by Eigen::SimplicialLLT< _MatrixType, _UpLo, _Ordering >::analyzePattern(), Eigen::SimplicialLDLT< _MatrixType, _UpLo, _Ordering >::analyzePattern(), and Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >::analyzePattern().

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

◆ analyzePattern_preordered()

template<typename Derived >
void Eigen::SimplicialCholeskyBase< Derived >::analyzePattern_preordered ( const CholMatrixType a,
bool  doLDLT 
)
protected
52{
53 const StorageIndex size = StorageIndex(ap.rows());
54 m_matrix.resize(size, size);
55 m_parent.resize(size);
57
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
VectorI m_parent
Definition SimplicialCholesky.h:257
bool m_analysisIsOk
Definition SimplicialCholesky.h:253
MatrixType::StorageIndex StorageIndex
Definition SimplicialCholesky.h:66
bool m_isInitialized
Definition SparseSolverBase.h:119
VectorI m_nonZerosPerCol
Definition SimplicialCholesky.h:258
Base::InnerIterator InnerIterator
Definition SparseMatrix.h:112
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

References ei_declare_aligned_stack_constructed_variable, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows(), and Eigen::Success.

Referenced by Eigen::SimplicialCholeskyBase< Derived >::analyzePattern(), and Eigen::SimplicialCholeskyBase< Derived >::compute().

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

◆ cols()

template<typename Derived >
Index Eigen::SimplicialCholeskyBase< Derived >::cols ( ) const
inline
99{ return m_matrix.cols(); }
Index cols() const
Definition SparseMatrix.h:138

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

+ Here is the call graph for this function:

◆ compute()

template<typename Derived >
template<bool DoLDLT>
void Eigen::SimplicialCholeskyBase< Derived >::compute ( const MatrixType matrix)
inlineprotected

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 }

References Eigen::SimplicialCholeskyBase< Derived >::analyzePattern_preordered(), eigen_assert, and ordering.

+ Here is the call graph for this function:

◆ derived() [1/4]

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

◆ derived() [2/4]

template<typename Derived >
Derived & Eigen::SimplicialCholeskyBase< Derived >::derived ( )
inline
96{ return *static_cast<Derived*>(this); }

Referenced by Eigen::SimplicialCholeskyBase< Derived >::SimplicialCholeskyBase(), Eigen::SimplicialCholeskyBase< Derived >::_solve_impl(), Eigen::SimplicialCholeskyBase< Derived >::_solve_impl(), and Eigen::SimplicialCholeskyBase< Derived >::setShift().

+ Here is the caller graph for this function:

◆ derived() [3/4]

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

◆ derived() [4/4]

template<typename Derived >
const Derived & Eigen::SimplicialCholeskyBase< Derived >::derived ( ) const
inline
97{ return *static_cast<const Derived*>(this); }

◆ dumpMemory()

template<typename Derived >
template<typename Stream >
void Eigen::SimplicialCholeskyBase< Derived >::dumpMemory ( Stream &  s)
inline
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

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::cols(), Eigen::SimplicialCholeskyBase< Derived >::m_diag, Eigen::SimplicialCholeskyBase< Derived >::m_matrix, Eigen::SimplicialCholeskyBase< Derived >::m_nonZerosPerCol, Eigen::SimplicialCholeskyBase< Derived >::m_P, Eigen::SimplicialCholeskyBase< Derived >::m_parent, Eigen::SimplicialCholeskyBase< Derived >::m_Pinv, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::nonZeros(), and Eigen::PermutationBase< Derived >::size().

+ Here is the call graph for this function:

◆ factorize()

template<typename Derived >
template<bool DoLDLT>
void Eigen::SimplicialCholeskyBase< Derived >::factorize ( const MatrixType a)
inlineprotected
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
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 }
@ Upper
Definition Constants.h:206
static void run(const InputMatrixType &input, ConstCholMatrixPtr &pmat, CholMatrixType &tmp)
Definition SimplicialCholesky.h:24

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_P, Eigen::internal::simplicial_cholesky_grab_input< CholMatrixType, InputMatrixType >::run(), Eigen::PermutationBase< Derived >::size(), Eigen::SimplicialCholeskyBase< Derived >::UpLo, and Eigen::Upper.

+ Here is the call graph for this function:

◆ factorize_preordered()

template<typename Derived >
template<bool DoLDLT>
void Eigen::SimplicialCholeskyBase< Derived >::factorize_preordered ( const CholMatrixType a)
protected
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
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

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::cols(), ei_declare_aligned_stack_constructed_variable, eigen_assert, Eigen::NumericalIssue, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::rows(), sqrt(), and Eigen::Success.

+ Here is the call graph for this function:

◆ info()

template<typename Derived >
ComputationInfo Eigen::SimplicialCholeskyBase< Derived >::info ( ) const
inline

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 }

References eigen_assert, Eigen::SimplicialCholeskyBase< Derived >::m_info, and Eigen::SimplicialCholeskyBase< Derived >::m_isInitialized.

Referenced by igl::embree::bone_heat(), igl::Frame_field_deformer::compute_optimal_positions(), igl::eigs(), and igl::Frame_field_deformer::precompute_opt().

+ Here is the caller graph for this function:

◆ ordering()

template<typename Derived >
void Eigen::SimplicialCholeskyBase< Derived >::ordering ( const MatrixType a,
ConstCholMatrixPtr pmat,
CholMatrixType ap 
)
protected
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
684 }
685}
void resize(Index newSize)
Definition PermutationMatrix.h:136
InverseReturnType inverse() const
Definition PermutationMatrix.h:196
internal::traits< Derived >::OrderingType OrderingType
Definition SimplicialCholesky.h:62
@ Lower
Definition Constants.h:204

References eigen_assert, Eigen::Lower, ordering, Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::resize(), and Eigen::internal::simplicial_cholesky_grab_input< CholMatrixType, InputMatrixType >::run().

+ Here is the call graph for this function:

◆ permutationP()

template<typename Derived >
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & Eigen::SimplicialCholeskyBase< Derived >::permutationP ( ) const
inline
Returns
the permutation P
See also
permutationPinv()
116 { return m_P; }

References Eigen::SimplicialCholeskyBase< Derived >::m_P.

◆ permutationPinv()

template<typename Derived >
const PermutationMatrix< Dynamic, Dynamic, StorageIndex > & Eigen::SimplicialCholeskyBase< Derived >::permutationPinv ( ) const
inline
Returns
the inverse P^-1 of the permutation P
See also
permutationP()
121 { return m_Pinv; }

References Eigen::SimplicialCholeskyBase< Derived >::m_Pinv.

◆ rows()

template<typename Derived >
Index Eigen::SimplicialCholeskyBase< Derived >::rows ( ) const
inline
100{ return m_matrix.rows(); }

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

+ Here is the call graph for this function:

◆ setShift()

template<typename Derived >
Derived & Eigen::SimplicialCholeskyBase< Derived >::setShift ( const RealScalar offset,
const RealScalar scale = 1 
)
inline

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

References Eigen::SimplicialCholeskyBase< Derived >::derived(), Eigen::SimplicialCholeskyBase< Derived >::m_shiftOffset, Eigen::SimplicialCholeskyBase< Derived >::m_shiftScale, and scale().

+ Here is the call graph for this function:

◆ 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

template<typename Derived >
bool Eigen::SimplicialCholeskyBase< Derived >::m_analysisIsOk
protected

◆ m_diag

◆ m_factorizationIsOk

◆ m_info

◆ m_isInitialized

template<typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprivate

◆ m_matrix

◆ m_nonZerosPerCol

template<typename Derived >
VectorI Eigen::SimplicialCholeskyBase< Derived >::m_nonZerosPerCol
protected

◆ m_P

◆ m_parent

template<typename Derived >
VectorI Eigen::SimplicialCholeskyBase< Derived >::m_parent
protected

◆ m_Pinv

◆ m_shiftOffset

template<typename Derived >
RealScalar Eigen::SimplicialCholeskyBase< Derived >::m_shiftOffset
protected

◆ m_shiftScale

template<typename Derived >
RealScalar Eigen::SimplicialCholeskyBase< Derived >::m_shiftScale
protected

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