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

Modified Incomplete Cholesky with dual threshold. More...

+ Inheritance diagram for Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >:

Public Types

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

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_impl (const Rhs &b, Dest &x) const
 
template<typename MatrixType >
void analyzePattern (const MatrixType &mat)
 Computes the fill reducing permutation vector using the sparsity pattern of mat. More...
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename MatrixType >
void compute (const MatrixType &mat)
 
template<typename MatrixType >
void factorize (const MatrixType &mat)
 Performs the numerical factorization of the input matrix mat. More...
 
template<typename MatrixType_ >
void factorize (const MatrixType_ &mat)
 
 IncompleteCholesky ()
 
template<typename MatrixType >
 IncompleteCholesky (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const FactorTypematrixL () const
 
const PermutationTypepermutationP () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
const VectorRxscalingS () const
 
void setInitialShift (RealScalar shift)
 Set the initial shift parameter \( \sigma \). More...
 
- Public Member Functions inherited from Eigen::SparseSolverBase< IncompleteCholesky< Scalar, Lower, AMDOrdering< int > > >
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
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Protected Types

typedef SparseSolverBase< IncompleteCholesky< Scalar, UpLo_, OrderingType_ > > Base
 

Protected Attributes

bool m_analysisIsOk
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
RealScalar m_initialShift
 
bool m_isInitialized
 
FactorType m_L
 
PermutationType m_perm
 
VectorRx m_scale
 
- 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>.

This class follows the sparse solver concept .

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.

Definition at line 46 of file IncompleteCholesky.h.

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

Definition at line 49 of file IncompleteCholesky.h.

◆ FactorType

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

Definition at line 56 of file IncompleteCholesky.h.

◆ OrderingType

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
typedef OrderingType_ Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::OrderingType

Definition at line 53 of file IncompleteCholesky.h.

◆ PermutationType

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

Definition at line 54 of file IncompleteCholesky.h.

◆ RealScalar

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

Definition at line 52 of file IncompleteCholesky.h.

◆ StorageIndex

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

Definition at line 55 of file IncompleteCholesky.h.

◆ VectorIx

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

Definition at line 59 of file IncompleteCholesky.h.

◆ VectorList

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

Definition at line 60 of file IncompleteCholesky.h.

◆ VectorRx

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

Definition at line 58 of file IncompleteCholesky.h.

◆ VectorSx

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

Definition at line 57 of file IncompleteCholesky.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
anonymous enum
Enumerator
UpLo 

Definition at line 61 of file IncompleteCholesky.h.

◆ anonymous enum

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 62 of file IncompleteCholesky.h.

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

Definition at line 74 of file IncompleteCholesky.h.

Array< double, 1, 3 > e(1./3., 0.5, 2.)

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

Definition at line 79 of file IncompleteCholesky.h.

80  {
81  compute(matrix);
82  }
void compute(const MatrixType &mat)

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

Definition at line 150 of file IncompleteCholesky.h.

151  {
152  eigen_assert(m_factorizationIsOk && "factorize() should be called first");
153  if (m_perm.rows() == b.rows()) x = m_perm * b;
154  else x = b;
155  x = m_scale.asDiagonal() * x;
156  x = m_L.template triangularView<Lower>().solve(x);
157  x = m_L.adjoint().template triangularView<Upper>().solve(x);
158  x = m_scale.asDiagonal() * x;
159  if (m_perm.rows() == b.rows())
160  x = m_perm.inverse() * x;
161  }
Array< int, 3, 1 > b
#define eigen_assert(x)
Definition: Macros.h:902
const AdjointReturnType adjoint() const

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

Definition at line 112 of file IncompleteCholesky.h.

113  {
114  OrderingType ord;
115  PermutationType pinv;
116  ord(mat.template selfadjointView<UpLo>(), pinv);
117  if(pinv.size()>0) m_perm = pinv.inverse();
118  else m_perm.resize(0);
119  m_L.resize(mat.rows(), mat.cols());
120  m_analysisIsOk = true;
121  m_isInitialized = true;
122  m_info = Success;
123  }
OrderingType::PermutationType PermutationType
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:739
@ Success
Definition: Constants.h:446

◆ cols()

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
EIGEN_CONSTEXPR Index Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::cols ( void  ) const
inline
Returns
number of columns of the factored matrix

Definition at line 88 of file IncompleteCholesky.h.

88 { return m_L.cols(); }
Index cols() const
Definition: SparseMatrix.h:167

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

Definition at line 142 of file IncompleteCholesky.h.

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.

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

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

Definition at line 191 of file IncompleteCholesky.h.

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

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

Definition at line 99 of file IncompleteCholesky.h.

100  {
101  eigen_assert(m_isInitialized && "IncompleteCholesky is not initialized.");
102  return m_info;
103  }

◆ 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

Definition at line 164 of file IncompleteCholesky.h.

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

Definition at line 170 of file IncompleteCholesky.h.

170 { eigen_assert(m_analysisIsOk && "analyzePattern() should be called first"); return m_perm; }

◆ rows()

template<typename Scalar , int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int>>
EIGEN_CONSTEXPR Index Eigen::IncompleteCholesky< Scalar, UpLo_, OrderingType_ >::rows ( void  ) const
inline
Returns
number of rows of the factored matrix

Definition at line 85 of file IncompleteCholesky.h.

85 { return m_L.rows(); }
Index rows() const
Definition: SparseMatrix.h:165

◆ 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

Definition at line 167 of file IncompleteCholesky.h.

167 { eigen_assert(m_factorizationIsOk && "factorize() should be called first"); return 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 \).

Definition at line 107 of file IncompleteCholesky.h.

107 { m_initialShift = shift; }

◆ 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

Definition at line 374 of file IncompleteCholesky.h.

375 {
376  if (jk < colPtr(col+1) )
377  {
378  Index p = colPtr(col+1) - jk;
379  Index minpos;
380  rowIdx.segment(jk,p).minCoeff(&minpos);
381  minpos += jk;
382  if (rowIdx(minpos) != rowIdx(jk))
383  {
384  //Swap
385  std::swap(rowIdx(jk),rowIdx(minpos));
386  std::swap(vals(jk),vals(minpos));
387  }
388  firstElt(col) = internal::convert_index<StorageIndex,Index>(jk);
389  listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex,Index>(col));
390  }
391 }
ColXpr col(Index i)
This is the const version of col().
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:788

Member Data Documentation

◆ m_analysisIsOk

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

Definition at line 176 of file IncompleteCholesky.h.

◆ m_factorizationIsOk

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

Definition at line 177 of file IncompleteCholesky.h.

◆ m_info

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

Definition at line 178 of file IncompleteCholesky.h.

◆ m_initialShift

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

Definition at line 175 of file IncompleteCholesky.h.

◆ m_isInitialized

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

Definition at line 123 of file SparseSolverBase.h.

◆ m_L

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

Definition at line 173 of file IncompleteCholesky.h.

◆ m_perm

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

Definition at line 179 of file IncompleteCholesky.h.

◆ m_scale

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

Definition at line 174 of file IncompleteCholesky.h.


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