Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ > Class Template Reference

Modified Incomplete Cholesky with dual threshold. More...

#include <IncompleteCholesky.h>

+ Inheritance 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)
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
void setInitialShift (RealScalar shift)
 Set the initial shift parameter \( \sigma \). More...
 
template<typename MatrixType >
void analyzePattern (const MatrixType &mat)
 Computes the fill reducing permutation vector using the sparsity pattern of mat. More...
 
template<typename MatrixType >
void factorize (const MatrixType &mat)
 Performs the numerical factorization of the input matrix mat. More...
 
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
 
RealScalar shift () const
 
template<typename MatrixType_ >
void factorize (const MatrixType_ &mat)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > & derived ()
 
const IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > & derived () const
 
const Solve< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
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
 
RealScalar m_shift
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > >
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
UpLo_The triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
OrderingType_The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<int>.

\implsparsesolverconcept

It performs the following incomplete factorization: \( S P A P' S + \sigma I \approx L L' \) where L is a lower triangular factor, S is a diagonal scaling matrix, P is a fill-in reducing permutation as computed by the ordering method, and \( \sigma \) is a shift for ensuring the decomposed matrix is positive definite.

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, and \sigma = 0. Otherwise, the factorization is performed on the shifted matrix \( B + \sigma I \) for a shifting factor \( \sigma \). We start with \( \sigma = \sigma_0 - \beta \), where \( \sigma_0 \) is the initial shift value as returned and set by setInitialShift() method. The default value is \( \sigma_0 = 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 
@ ColsAtCompileTime
Definition: IncompleteCholesky.h:67
@ MaxColsAtCompileTime
Definition: IncompleteCholesky.h:67
const int Dynamic
Definition: Constants.h:25

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&)
76 : m_initialShift(1e-3), m_analysisIsOk(false), m_factorizationIsOk(false) {}
Array< double, 1, 3 > e(1./3., 0.5, 2.)
bool m_factorizationIsOk
Definition: IncompleteCholesky.h:189
bool m_analysisIsOk
Definition: IncompleteCholesky.h:188
RealScalar m_initialShift
Definition: IncompleteCholesky.h:187

◆ 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.

82  : m_initialShift(1e-3), m_analysisIsOk(false), m_factorizationIsOk(false) {
83  compute(matrix);
84  }
void compute(const MatrixType &mat)
Definition: IncompleteCholesky.h:143
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::compute(), and matrix().

Member Function Documentation

◆ _solve_impl()

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
150  {
151  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
152  if (m_perm.rows() == b.rows())
153  x = m_perm * b;
154  else
155  x = b;
156  x = m_scale.asDiagonal() * x;
157  x = m_L.template triangularView<Lower>().solve(x);
158  x = m_L.adjoint().template triangularView<Upper>().solve(x);
159  x = m_scale.asDiagonal() * x;
160  if (m_perm.rows() == b.rows()) x = m_perm.inverse() * x;
161  }
#define eigen_assert(x)
Definition: Macros.h:910
Scalar * b
Definition: benchVecAdd.cpp:17
PermutationType m_perm
Definition: IncompleteCholesky.h:191
VectorRx m_scale
Definition: IncompleteCholesky.h:186
FactorType m_L
Definition: IncompleteCholesky.h:185
const AdjointReturnType adjoint() const
Definition: SparseMatrixBase.h:360
list x
Definition: plotDoE.py:28

References Eigen::SparseMatrixBase< Derived >::adjoint(), b, eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_L, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_perm, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_scale, and plotDoE::x.

◆ 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.

112  {
113  OrderingType ord;
115  ord(mat.template selfadjointView<UpLo>(), pinv);
116  if (pinv.size() > 0)
117  m_perm = pinv.inverse();
118  else
119  m_perm.resize(0);
120  m_L.resize(mat.rows(), mat.cols());
121  m_analysisIsOk = true;
122  m_isInitialized = true;
123  m_info = Success;
124  }
ComputationInfo m_info
Definition: IncompleteCholesky.h:190
OrderingType_ OrderingType
Definition: IncompleteCholesky.h:58
OrderingType::PermutationType PermutationType
Definition: IncompleteCholesky.h:59
bool m_isInitialized
Definition: SparseSolverBase.h:110
Index cols() const
Definition: SparseMatrix.h:161
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:734
Index rows() const
Definition: SparseMatrix.h:159
@ Success
Definition: Constants.h:440
bool pinv
Definition: jeffery_hamel.cc:86

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), 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, Global_Physical_Variables::pinv, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), and Eigen::Success.

Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::compute().

◆ cols()

◆ 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()
143  {
145  factorize(mat);
146  }
void factorize(const MatrixType &mat)
Performs the numerical factorization of the input matrix mat.
void analyzePattern(const MatrixType &mat)
Computes the fill reducing permutation vector using the sparsity pattern of mat.
Definition: IncompleteCholesky.h:112

References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::analyzePattern(), and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::factorize().

Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::IncompleteCholesky().

◆ 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)

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().

◆ 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)
205  {
206  using std::sqrt;
207  eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
208 
209  // Dropping strategy : Keep only the p largest elements per column, where p is the number of elements in the column of
210  // the original matrix. Other strategies will be added
211 
212  // Apply the fill-reducing permutation computed in analyzePattern()
213  if (m_perm.rows() == mat.rows()) // To detect the null permutation
214  {
215  // The temporary is needed to make sure that the diagonal entry is properly sorted
216  FactorType tmp(mat.rows(), mat.cols());
217  tmp = mat.template selfadjointView<UpLo_>().twistedBy(m_perm);
218  m_L.template selfadjointView<Lower>() = tmp.template selfadjointView<Lower>();
219  } else {
220  m_L.template selfadjointView<Lower>() = mat.template selfadjointView<UpLo_>();
221  }
222 
223  // The algorithm will insert increasingly large shifts on the diagonal until
224  // factorization succeeds. Therefore we have to make sure that there is a
225  // space in the datastructure to store such values, even if the original
226  // matrix has a zero on the diagonal.
227  bool modified = false;
228  for (Index i = 0; i < mat.cols(); ++i) {
229  bool inserted = false;
230  m_L.findOrInsertCoeff(i, i, &inserted);
231  if (inserted) {
232  modified = true;
233  }
234  }
235  if (modified) m_L.makeCompressed();
236 
237  Index n = m_L.cols();
238  Index nnz = m_L.nonZeros();
239  Map<VectorSx> vals(m_L.valuePtr(), nnz); // values
240  Map<VectorIx> rowIdx(m_L.innerIndexPtr(), nnz); // Row indices
241  Map<VectorIx> colPtr(m_L.outerIndexPtr(), n + 1); // Pointer to the beginning of each row
242  VectorIx firstElt(n - 1); // for each j, points to the next entry in vals that will be used in the factorization
243  VectorList listCol(n); // listCol(j) is a linked list of columns to update column j
244  VectorSx col_vals(n); // Store a nonzero values in each column
245  VectorIx col_irow(n); // Row indices of nonzero elements in each column
246  VectorIx col_pattern(n);
247  col_pattern.fill(-1);
248  StorageIndex col_nnz;
249 
250  // Computes the scaling factors
251  m_scale.resize(n);
252  m_scale.setZero();
253  for (Index j = 0; j < n; j++)
254  for (Index k = colPtr[j]; k < colPtr[j + 1]; k++) {
255  m_scale(j) += numext::abs2(vals(k));
256  if (rowIdx[k] != j) m_scale(rowIdx[k]) += numext::abs2(vals(k));
257  }
258 
259  m_scale = m_scale.cwiseSqrt().cwiseSqrt();
260 
261  for (Index j = 0; j < n; ++j)
263  m_scale(j) = RealScalar(1) / m_scale(j);
264  else
265  m_scale(j) = 1;
266 
267  // TODO disable scaling if not needed, i.e., if it is roughly uniform? (this will make solve() faster)
268 
269  // Scale and compute the shift for the matrix
270  RealScalar mindiag = NumTraits<RealScalar>::highest();
271  for (Index j = 0; j < n; j++) {
272  for (Index k = colPtr[j]; k < colPtr[j + 1]; k++) vals[k] *= (m_scale(j) * m_scale(rowIdx[k]));
273  eigen_internal_assert(rowIdx[colPtr[j]] == j &&
274  "IncompleteCholesky: only the lower triangular part must be stored");
275  mindiag = numext::mini(numext::real(vals[colPtr[j]]), mindiag);
276  }
277 
278  FactorType L_save = m_L;
279 
280  m_shift = RealScalar(0);
281  if (mindiag <= RealScalar(0.)) m_shift = m_initialShift - mindiag;
282 
284 
285  // Try to perform the incomplete factorization using the current shift
286  int iter = 0;
287  do {
288  // Apply the shift to the diagonal elements of the matrix
289  for (Index j = 0; j < n; j++) vals[colPtr[j]] += m_shift;
290 
291  // jki version of the Cholesky factorization
292  Index j = 0;
293  for (; j < n; ++j) {
294  // Left-looking factorization of the j-th column
295  // First, load the j-th column into col_vals
296  Scalar diag = vals[colPtr[j]]; // It is assumed that only the lower part is stored
297  col_nnz = 0;
298  for (Index i = colPtr[j] + 1; i < colPtr[j + 1]; i++) {
299  StorageIndex l = rowIdx[i];
300  col_vals(col_nnz) = vals[i];
301  col_irow(col_nnz) = l;
302  col_pattern(l) = col_nnz;
303  col_nnz++;
304  }
305  {
306  typename std::list<StorageIndex>::iterator k;
307  // Browse all previous columns that will update column j
308  for (k = listCol[j].begin(); k != listCol[j].end(); k++) {
309  Index jk = firstElt(*k); // First element to use in the column
310  eigen_internal_assert(rowIdx[jk] == j);
311  Scalar v_j_jk = numext::conj(vals[jk]);
312 
313  jk += 1;
314  for (Index i = jk; i < colPtr[*k + 1]; i++) {
315  StorageIndex l = rowIdx[i];
316  if (col_pattern[l] < 0) {
317  col_vals(col_nnz) = vals[i] * v_j_jk;
318  col_irow[col_nnz] = l;
319  col_pattern(l) = col_nnz;
320  col_nnz++;
321  } else
322  col_vals(col_pattern[l]) -= vals[i] * v_j_jk;
323  }
324  updateList(colPtr, rowIdx, vals, *k, jk, firstElt, listCol);
325  }
326  }
327 
328  // Scale the current column
329  if (numext::real(diag) <= 0) {
330  if (++iter >= 10) return;
331 
332  // increase shift
334  // restore m_L, col_pattern, and listCol
335  vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
336  rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
337  colPtr = Map<const VectorIx>(L_save.outerIndexPtr(), n + 1);
338  col_pattern.fill(-1);
339  for (Index i = 0; i < n; ++i) listCol[i].clear();
340 
341  break;
342  }
343 
344  RealScalar rdiag = sqrt(numext::real(diag));
345  vals[colPtr[j]] = rdiag;
346  for (Index k = 0; k < col_nnz; ++k) {
347  Index i = col_irow[k];
348  // Scale
349  col_vals(k) /= rdiag;
350  // Update the remaining diagonals with col_vals
351  vals[colPtr[i]] -= numext::abs2(col_vals(k));
352  }
353  // Select the largest p elements
354  // p is the original number of elements in the column (without the diagonal)
355  Index p = colPtr[j + 1] - colPtr[j] - 1;
356  Ref<VectorSx> cvals = col_vals.head(col_nnz);
357  Ref<VectorIx> cirow = col_irow.head(col_nnz);
358  internal::QuickSplit(cvals, cirow, p);
359  // Insert the largest p elements in the matrix
360  Index cpt = 0;
361  for (Index i = colPtr[j] + 1; i < colPtr[j + 1]; i++) {
362  vals[i] = col_vals(cpt);
363  rowIdx[i] = col_irow(cpt);
364  // restore col_pattern:
365  col_pattern(col_irow(cpt)) = -1;
366  cpt++;
367  }
368  // Get the first smallest row index and put it after the diagonal element
369  Index jk = colPtr(j) + 1;
370  updateList(colPtr, rowIdx, vals, j, jk, firstElt, listCol);
371  }
372 
373  if (j == n) {
374  m_factorizationIsOk = true;
375  m_info = Success;
376  }
377  } while (m_info != Success);
378 }
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define eigen_internal_assert(x)
Definition: Macros.h:916
float * p
Definition: Tutorial_Map_using.cpp:9
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
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:381
Matrix< StorageIndex, Dynamic, 1 > VectorIx
Definition: IncompleteCholesky.h:64
std::vector< std::list< StorageIndex > > VectorList
Definition: IncompleteCholesky.h:65
Matrix< Scalar, Dynamic, 1 > VectorSx
Definition: IncompleteCholesky.h:62
RealScalar m_shift
Definition: IncompleteCholesky.h:192
PermutationType::StorageIndex StorageIndex
Definition: IncompleteCholesky.h:60
SparseMatrix< Scalar, ColMajor, StorageIndex > FactorType
Definition: IncompleteCholesky.h:61
NumTraits< Scalar >::Real RealScalar
Definition: IncompleteCholesky.h:57
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
SparseSymmetricPermutationProduct< Derived, Upper|Lower > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition: SparseMatrixBase.h:325
Index nonZeros() const
Definition: SparseCompressedBase.h:64
void makeCompressed()
Definition: SparseMatrix.h:589
const Scalar * valuePtr() const
Definition: SparseMatrix.h:171
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:189
Scalar & findOrInsertCoeff(Index row, Index col, bool *inserted)
Definition: SparseMatrix.h:231
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:180
float real
Definition: datatypes.h:10
#define min(a, b)
Definition: datatypes.h:22
@ NumericalIssue
Definition: Constants.h:442
const char const char const char * diag
Definition: level2_impl.h:86
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
Definition: IncompleteLUT.h:31
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:926
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Eigen::numext::abs2(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), conj(), diag, eigen_assert, eigen_internal_assert, i, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::innerIndexPtr(), j, k, Eigen::numext::maxi(), min, Eigen::numext::mini(), n, Eigen::NumericalIssue, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::outerIndexPtr(), p, Eigen::internal::QuickSplit(), Eigen::PlainObjectBase< Derived >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), sqrt(), Eigen::Success, tmp, Eigen::SparseMatrixBase< Derived >::twistedBy(), and Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::valuePtr().

◆ 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.
100  {
101  eigen_assert(m_isInitialized && "IncompleteCholesky is not initialized.");
102  return m_info;
103  }

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
164  {
165  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
166  return m_L;
167  }

References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, 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)
176  {
177  eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
178  return m_perm;
179  }

References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_analysisIsOk, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_perm.

◆ rows()

◆ 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
170  {
171  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
172  return m_scale;
173  }

References eigen_assert, Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_factorizationIsOk, 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 \).

107 { m_initialShift = shift; }
RealScalar shift() const
Definition: IncompleteCholesky.h:182

References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_initialShift, and Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::shift().

◆ shift()

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
RealScalar Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::shift ( ) const
inline
Returns
the final shift parameter from the computation
182 { return m_shift; }

References Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_shift.

Referenced by Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::setInitialShift().

◆ 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
384  {
385  if (jk < colPtr(col + 1)) {
386  Index p = colPtr(col + 1) - jk;
387  Index minpos;
388  rowIdx.segment(jk, p).minCoeff(&minpos);
389  minpos += jk;
390  if (rowIdx(minpos) != rowIdx(jk)) {
391  // Swap
392  std::swap(rowIdx(jk), rowIdx(minpos));
393  std::swap(vals(jk), vals(minpos));
394  }
395  firstElt(col) = internal::convert_index<StorageIndex, Index>(jk);
396  listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex, Index>(col));
397  }
398 }
m col(1)
EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:117

References col(), p, and swap().

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

◆ m_shift

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
RealScalar Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::m_shift
protected

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