$darkmode
Eigen  5.0.1-dev
SimplicialCholesky_impl.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2012 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 /*
11 NOTE: these functions have been adapted from the LDL library:
12 
13 LDL Copyright (c) 2005 by Timothy A. Davis. All Rights Reserved.
14 
15 The author of LDL, Timothy A. Davis., has executed a license with Google LLC
16 to permit distribution of this code and derivative works as part of Eigen under
17 the Mozilla Public License v. 2.0, as stated at the top of this file.
18  */
19 
20 #ifndef EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
21 #define EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
22 
23 // IWYU pragma: private
24 #include "./InternalHeaderCheck.h"
25 
26 namespace Eigen {
27 
28 namespace internal {
29 
30 template <typename Scalar, typename StorageIndex>
31 struct simpl_chol_helper {
32  using CholMatrixType = SparseMatrix<Scalar, ColMajor, StorageIndex>;
33  using InnerIterator = typename CholMatrixType::InnerIterator;
34  using VectorI = Matrix<StorageIndex, Dynamic, 1>;
35  static constexpr StorageIndex kEmpty = -1;
36 
37  // Implementation of a stack or last-in first-out structure with some debugging machinery.
38  struct Stack {
39  StorageIndex* m_data;
40  Index m_size;
41 #ifndef EIGEN_NO_DEBUG
42  const Index m_maxSize;
43  Stack(StorageIndex* data, StorageIndex size, StorageIndex maxSize)
44  : m_data(data), m_size(size), m_maxSize(maxSize) {
45  eigen_assert(size >= 0);
46  eigen_assert(maxSize >= size);
47  }
48 #else
49  Stack(StorageIndex* data, StorageIndex size, StorageIndex /*maxSize*/) : m_data(data), m_size(size) {}
50 #endif
51  bool empty() const { return m_size == 0; }
52  Index size() const { return m_size; }
53  StorageIndex back() const {
54  eigen_assert(m_size > 0);
55  return m_data[m_size - 1];
56  }
57  void push(const StorageIndex& value) {
58 #ifndef EIGEN_NO_DEBUG
59  eigen_assert(m_size < m_maxSize);
60 #endif
61  m_data[m_size] = value;
62  m_size++;
63  }
64  void pop() {
65  eigen_assert(m_size > 0);
66  m_size--;
67  }
68  };
69 
70  // Implementation of a disjoint-set or union-find structure with path compression.
71  struct DisjointSet {
72  StorageIndex* m_set;
73  DisjointSet(StorageIndex* set, StorageIndex size) : m_set(set) { std::iota(set, set + size, 0); }
74  // Find the set representative or root of `u`.
75  StorageIndex find(StorageIndex u) const {
76  eigen_assert(u != kEmpty);
77  while (m_set[u] != u) {
78  // manually unroll the loop by a factor of 2 to improve performance
79  u = m_set[m_set[u]];
80  }
81  return u;
82  }
83  // Perform full path compression such that each node from `u` to `v` points to `v`.
84  void compress(StorageIndex u, StorageIndex v) {
85  eigen_assert(u != kEmpty);
86  eigen_assert(v != kEmpty);
87  while (m_set[u] != v) {
88  StorageIndex next = m_set[u];
89  m_set[u] = v;
90  u = next;
91  }
92  };
93  };
94 
95  // Computes the higher adjacency pattern by transposing the input lower adjacency matrix.
96  // Only the index arrays are calculated, as the values are not needed for the symbolic factorization.
97  // The outer index array provides the size requirements of the inner index array.
98 
99  // Computes the outer index array of the higher adjacency matrix.
100  static void calc_hadj_outer(const StorageIndex size, const CholMatrixType& ap, StorageIndex* outerIndex) {
101  for (StorageIndex j = 1; j < size; j++) {
102  for (InnerIterator it(ap, j); it; ++it) {
103  StorageIndex i = it.index();
104  if (i < j) outerIndex[i + 1]++;
105  }
106  }
107  std::partial_sum(outerIndex, outerIndex + size + 1, outerIndex);
108  }
109 
110  // inner index array
111  static void calc_hadj_inner(const StorageIndex size, const CholMatrixType& ap, const StorageIndex* outerIndex,
112  StorageIndex* innerIndex, StorageIndex* tmp) {
113  std::fill_n(tmp, size, 0);
114 
115  for (StorageIndex j = 1; j < size; j++) {
116  for (InnerIterator it(ap, j); it; ++it) {
117  StorageIndex i = it.index();
118  if (i < j) {
119  StorageIndex b = outerIndex[i] + tmp[i];
120  innerIndex[b] = j;
121  tmp[i]++;
122  }
123  }
124  }
125  }
126 
127  // Adapted from:
128  // Joseph W. Liu. (1986).
129  // A compact row storage scheme for Cholesky factors using elimination trees.
130  // ACM Trans. Math. Softw. 12, 2 (June 1986), 127-148. https://doi.org/10.1145/6497.6499
131 
132  // Computes the elimination forest of the lower adjacency matrix, a compact representation of the sparse L factor.
133  // The L factor may contain multiple elimination trees if a column contains only its diagonal element.
134  // Each elimination tree is an n-ary tree in which each node points to its parent.
135  static void calc_etree(const StorageIndex size, const CholMatrixType& ap, StorageIndex* parent, StorageIndex* tmp) {
136  std::fill_n(parent, size, kEmpty);
137 
138  DisjointSet ancestor(tmp, size);
139 
140  for (StorageIndex j = 1; j < size; j++) {
141  for (InnerIterator it(ap, j); it; ++it) {
142  StorageIndex i = it.index();
143  if (i < j) {
144  StorageIndex r = ancestor.find(i);
145  if (r != j) parent[r] = j;
146  ancestor.compress(i, j);
147  }
148  }
149  }
150  }
151 
152  // Computes the child pointers of the parent tree to facilitate a depth-first search traversal.
153  static void calc_lineage(const StorageIndex size, const StorageIndex* parent, StorageIndex* firstChild,
154  StorageIndex* firstSibling) {
155  std::fill_n(firstChild, size, kEmpty);
156  std::fill_n(firstSibling, size, kEmpty);
157 
158  for (StorageIndex j = 0; j < size; j++) {
159  StorageIndex p = parent[j];
160  if (p == kEmpty) continue;
161  StorageIndex c = firstChild[p];
162  if (c == kEmpty)
163  firstChild[p] = j;
164  else {
165  while (firstSibling[c] != kEmpty) c = firstSibling[c];
166  firstSibling[c] = j;
167  }
168  }
169  }
170 
171  // Computes a post-ordered traversal of the elimination tree.
172  static void calc_post(const StorageIndex size, const StorageIndex* parent, StorageIndex* firstChild,
173  const StorageIndex* firstSibling, StorageIndex* post, StorageIndex* dfs) {
174  Stack post_stack(post, 0, size);
175  for (StorageIndex j = 0; j < size; j++) {
176  if (parent[j] != kEmpty) continue;
177  // Begin at a root
178  Stack dfs_stack(dfs, 0, size);
179  dfs_stack.push(j);
180  while (!dfs_stack.empty()) {
181  StorageIndex i = dfs_stack.back();
182  StorageIndex c = firstChild[i];
183  if (c == kEmpty) {
184  post_stack.push(i);
185  dfs_stack.pop();
186  } else {
187  dfs_stack.push(c);
188  // Remove the path from `i` to `c` for future traversals.
189  firstChild[i] = firstSibling[c];
190  }
191  }
192  }
193  eigen_assert(post_stack.size() == size);
194  eigen_assert(std::all_of(firstChild, firstChild + size, [](StorageIndex a) { return a == kEmpty; }));
195  }
196 
197  // Adapted from:
198  // Gilbert, J. R., Ng, E., & Peyton, B. W. (1994).
199  // An efficient algorithm to compute row and column counts for sparse Cholesky factorization.
200  // SIAM Journal on Matrix Analysis and Applications, 15(4), 1075-1091.
201 
202  // Computes the non-zero pattern of the L factor.
203  static void calc_colcount(const StorageIndex size, const StorageIndex* hadjOuter, const StorageIndex* hadjInner,
204  const StorageIndex* parent, StorageIndex* prevLeaf, StorageIndex* tmp,
205  const StorageIndex* post, StorageIndex* nonZerosPerCol, bool doLDLT) {
206  // initialize nonZerosPerCol with 1 for leaves, 0 for non-leaves
207  std::fill_n(nonZerosPerCol, size, 1);
208  for (StorageIndex j = 0; j < size; j++) {
209  StorageIndex p = parent[j];
210  // p is not a leaf
211  if (p != kEmpty) nonZerosPerCol[p] = 0;
212  }
213 
214  DisjointSet parentSet(tmp, size);
215  // prevLeaf is already initialized
216  eigen_assert(std::all_of(prevLeaf, prevLeaf + size, [](StorageIndex a) { return a == kEmpty; }));
217 
218  for (StorageIndex j_ = 0; j_ < size; j_++) {
219  StorageIndex j = post[j_];
220  nonZerosPerCol[j] += hadjOuter[j + 1] - hadjOuter[j];
221  for (StorageIndex k = hadjOuter[j]; k < hadjOuter[j + 1]; k++) {
222  StorageIndex i = hadjInner[k];
223  eigen_assert(i > j);
224  StorageIndex prev = prevLeaf[i];
225  if (prev != kEmpty) {
226  StorageIndex q = parentSet.find(prev);
227  parentSet.compress(prev, q);
228  nonZerosPerCol[q]--;
229  }
230  prevLeaf[i] = j;
231  }
232  StorageIndex p = parent[j];
233  if (p != kEmpty) parentSet.compress(j, p);
234  }
235 
236  for (StorageIndex j = 0; j < size; j++) {
237  StorageIndex p = parent[j];
238  if (p != kEmpty) nonZerosPerCol[p] += nonZerosPerCol[j] - 1;
239  if (doLDLT) nonZerosPerCol[j]--;
240  }
241  }
242 
243  // Finalizes the non zero pattern of the L factor and allocates the memory for the factorization.
244  static void init_matrix(const StorageIndex size, const StorageIndex* nonZerosPerCol, CholMatrixType& L) {
245  eigen_assert(L.outerIndexPtr()[0] == 0);
246  std::partial_sum(nonZerosPerCol, nonZerosPerCol + size, L.outerIndexPtr() + 1);
247  L.resizeNonZeros(L.outerIndexPtr()[size]);
248  }
249 
250  // Driver routine for the symbolic sparse Cholesky factorization.
251  static void run(const StorageIndex size, const CholMatrixType& ap, CholMatrixType& L, VectorI& parent,
252  VectorI& workSpace, bool doLDLT) {
253  parent.resize(size);
254  workSpace.resize(4 * size);
255  L.resize(size, size);
256 
257  StorageIndex* tmp1 = workSpace.data();
258  StorageIndex* tmp2 = workSpace.data() + size;
259  StorageIndex* tmp3 = workSpace.data() + 2 * size;
260  StorageIndex* tmp4 = workSpace.data() + 3 * size;
261 
262  // Borrow L's outer index array for the higher adjacency pattern.
263  StorageIndex* hadj_outer = L.outerIndexPtr();
264  calc_hadj_outer(size, ap, hadj_outer);
265  // Request additional temporary storage for the inner indices of the higher adjacency pattern.
266  ei_declare_aligned_stack_constructed_variable(StorageIndex, hadj_inner, hadj_outer[size], nullptr);
267  calc_hadj_inner(size, ap, hadj_outer, hadj_inner, tmp1);
268 
269  calc_etree(size, ap, parent.data(), tmp1);
270  calc_lineage(size, parent.data(), tmp1, tmp2);
271  calc_post(size, parent.data(), tmp1, tmp2, tmp3, tmp4);
272  calc_colcount(size, hadj_outer, hadj_inner, parent.data(), tmp1, tmp2, tmp3, tmp4, doLDLT);
273  init_matrix(size, tmp4, L);
274  }
275 };
276 
277 // Symbol is ODR-used, so we need a definition.
278 template <typename Scalar, typename StorageIndex>
279 constexpr StorageIndex simpl_chol_helper<Scalar, StorageIndex>::kEmpty;
280 
281 } // namespace internal
282 
283 template <typename Derived>
284 void SimplicialCholeskyBase<Derived>::analyzePattern_preordered(const CholMatrixType& ap, bool doLDLT) {
285  using Helper = internal::simpl_chol_helper<Scalar, StorageIndex>;
286 
287  eigen_assert(ap.innerSize() == ap.outerSize());
288  const StorageIndex size = internal::convert_index<StorageIndex>(ap.outerSize());
289 
290  Helper::run(size, ap, m_matrix, m_parent, m_workSpace, doLDLT);
291 
292  m_isInitialized = true;
293  m_info = Success;
294  m_analysisIsOk = true;
295  m_factorizationIsOk = false;
296 }
297 
298 template <typename Derived>
299 template <bool DoLDLT, bool NonHermitian>
300 void SimplicialCholeskyBase<Derived>::factorize_preordered(const CholMatrixType& ap) {
301  using std::sqrt;
302  const StorageIndex size = StorageIndex(ap.rows());
303 
304  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
305  eigen_assert(ap.rows() == ap.cols());
306  eigen_assert(m_parent.size() == size);
307  eigen_assert(m_workSpace.size() >= 3 * size);
308 
309  const StorageIndex* Lp = m_matrix.outerIndexPtr();
310  StorageIndex* Li = m_matrix.innerIndexPtr();
311  Scalar* Lx = m_matrix.valuePtr();
312 
313  ei_declare_aligned_stack_constructed_variable(Scalar, y, size, 0);
314  StorageIndex* nonZerosPerCol = m_workSpace.data();
315  StorageIndex* pattern = m_workSpace.data() + size;
316  StorageIndex* tags = m_workSpace.data() + 2 * size;
317 
318  bool ok = true;
319  m_diag.resize(DoLDLT ? size : 0);
320 
321  for (StorageIndex k = 0; k < size; ++k) {
322  // compute nonzero pattern of kth row of L, in topological order
323  y[k] = Scalar(0); // Y(0:k) is now all zero
324  StorageIndex top = size; // stack for pattern is empty
325  tags[k] = k; // mark node k as visited
326  nonZerosPerCol[k] = 0; // count of nonzeros in column k of L
327  for (typename CholMatrixType::InnerIterator it(ap, k); it; ++it) {
328  StorageIndex i = it.index();
329  if (i <= k) {
330  y[i] += getSymm(it.value()); /* scatter A(i,k) into Y (sum duplicates) */
331  Index len;
332  for (len = 0; tags[i] != k; i = m_parent[i]) {
333  pattern[len++] = i; /* L(k,i) is nonzero */
334  tags[i] = k; /* mark i as visited */
335  }
336  while (len > 0) pattern[--top] = pattern[--len];
337  }
338  }
339 
340  /* compute numerical values kth row of L (a sparse triangular solve) */
341 
342  DiagonalScalar d =
343  getDiag(y[k]) * m_shiftScale + m_shiftOffset; // get D(k,k), apply the shift function, and clear Y(k)
344  y[k] = Scalar(0);
345  for (; top < size; ++top) {
346  Index i = pattern[top]; /* pattern[top:n-1] is pattern of L(:,k) */
347  Scalar yi = y[i]; /* get and clear Y(i) */
348  y[i] = Scalar(0);
349 
350  /* the nonzero entry L(k,i) */
351  Scalar l_ki;
352  if (DoLDLT)
353  l_ki = yi / getDiag(m_diag[i]);
354  else
355  yi = l_ki = yi / Lx[Lp[i]];
356 
357  Index p2 = Lp[i] + nonZerosPerCol[i];
358  Index p;
359  for (p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p) y[Li[p]] -= getSymm(Lx[p]) * yi;
360  d -= getDiag(l_ki * getSymm(yi));
361  Li[p] = k; /* store L(k,i) in column form of L */
362  Lx[p] = l_ki;
363  ++nonZerosPerCol[i]; /* increment count of nonzeros in col i */
364  }
365  if (DoLDLT) {
366  m_diag[k] = d;
367  if (d == RealScalar(0)) {
368  ok = false; /* failure, D(k,k) is zero */
369  break;
370  }
371  } else {
372  Index p = Lp[k] + nonZerosPerCol[k]++;
373  Li[p] = k; /* store L(k,k) = sqrt (d) in column k */
374  if (NonHermitian ? d == RealScalar(0) : numext::real(d) <= RealScalar(0)) {
375  ok = false; /* failure, matrix is not positive definite */
376  break;
377  }
378  Lx[p] = sqrt(d);
379  }
380  }
381 
382  m_info = ok ? Success : NumericalIssue;
383  m_factorizationIsOk = true;
384 }
385 
386 } // end namespace Eigen
387 
388 #endif // EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Definition: Constants.h:442
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: Constants.h:440