Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
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 <[email protected]>
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/*
11NOTE: these functions have been adapted from the LDL library:
12
13LDL Copyright (c) 2005 by Timothy A. Davis. All Rights Reserved.
14
15The author of LDL, Timothy A. Davis., has executed a license with Google LLC
16to permit distribution of this code and derivative works as part of Eigen under
17the 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
26namespace Eigen {
27
28template <typename Derived>
29void SimplicialCholeskyBase<Derived>::analyzePattern_preordered(const CholMatrixType& ap, bool doLDLT) {
30 const StorageIndex size = StorageIndex(ap.rows());
31 m_matrix.resize(size, size);
32 m_parent.resize(size);
33 m_nonZerosPerCol.resize(size);
34
35 ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
36
37 for (StorageIndex k = 0; k < size; ++k) {
38 /* L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) */
39 m_parent[k] = -1; /* parent of k is not yet known */
40 tags[k] = k; /* mark node k as visited */
41 m_nonZerosPerCol[k] = 0; /* count of nonzeros in column k of L */
42 for (typename CholMatrixType::InnerIterator it(ap, k); it; ++it) {
43 StorageIndex i = it.index();
44 if (i < k) {
45 /* follow path from i to root of etree, stop at flagged node */
46 for (; tags[i] != k; i = m_parent[i]) {
47 /* find parent of i if not yet determined */
48 if (m_parent[i] == -1) m_parent[i] = k;
49 m_nonZerosPerCol[i]++; /* L (k,i) is nonzero */
50 tags[i] = k; /* mark i as visited */
51 }
52 }
53 }
54 }
55
56 /* construct Lp index array from m_nonZerosPerCol column counts */
57 StorageIndex* Lp = m_matrix.outerIndexPtr();
58 Lp[0] = 0;
59 for (StorageIndex k = 0; k < size; ++k) Lp[k + 1] = Lp[k] + m_nonZerosPerCol[k] + (doLDLT ? 0 : 1);
60
61 m_matrix.resizeNonZeros(Lp[size]);
62
63 m_isInitialized = true;
64 m_info = Success;
65 m_analysisIsOk = true;
66 m_factorizationIsOk = false;
67}
68
69template <typename Derived>
70template <bool DoLDLT, bool NonHermitian>
71void SimplicialCholeskyBase<Derived>::factorize_preordered(const CholMatrixType& ap) {
72 using std::sqrt;
73
74 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
75 eigen_assert(ap.rows() == ap.cols());
76 eigen_assert(m_parent.size() == ap.rows());
77 eigen_assert(m_nonZerosPerCol.size() == ap.rows());
78
79 const StorageIndex size = StorageIndex(ap.rows());
80 const StorageIndex* Lp = m_matrix.outerIndexPtr();
81 StorageIndex* Li = m_matrix.innerIndexPtr();
82 Scalar* Lx = m_matrix.valuePtr();
83
84 ei_declare_aligned_stack_constructed_variable(Scalar, y, size, 0);
85 ei_declare_aligned_stack_constructed_variable(StorageIndex, pattern, size, 0);
86 ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
87
88 bool ok = true;
89 m_diag.resize(DoLDLT ? size : 0);
90
91 for (StorageIndex k = 0; k < size; ++k) {
92 // compute nonzero pattern of kth row of L, in topological order
93 y[k] = Scalar(0); // Y(0:k) is now all zero
94 StorageIndex top = size; // stack for pattern is empty
95 tags[k] = k; // mark node k as visited
96 m_nonZerosPerCol[k] = 0; // count of nonzeros in column k of L
97 for (typename CholMatrixType::InnerIterator it(ap, k); it; ++it) {
98 StorageIndex i = it.index();
99 if (i <= k) {
100 y[i] += getSymm(it.value()); /* scatter A(i,k) into Y (sum duplicates) */
101 Index len;
102 for (len = 0; tags[i] != k; i = m_parent[i]) {
103 pattern[len++] = i; /* L(k,i) is nonzero */
104 tags[i] = k; /* mark i as visited */
105 }
106 while (len > 0) pattern[--top] = pattern[--len];
107 }
108 }
109
110 /* compute numerical values kth row of L (a sparse triangular solve) */
111
112 DiagonalScalar d =
113 getDiag(y[k]) * m_shiftScale + m_shiftOffset; // get D(k,k), apply the shift function, and clear Y(k)
114 y[k] = Scalar(0);
115 for (; top < size; ++top) {
116 Index i = pattern[top]; /* pattern[top:n-1] is pattern of L(:,k) */
117 Scalar yi = y[i]; /* get and clear Y(i) */
118 y[i] = Scalar(0);
119
120 /* the nonzero entry L(k,i) */
121 Scalar l_ki;
122 if (DoLDLT)
123 l_ki = yi / getDiag(m_diag[i]);
124 else
125 yi = l_ki = yi / Lx[Lp[i]];
126
127 Index p2 = Lp[i] + m_nonZerosPerCol[i];
128 Index p;
129 for (p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p) y[Li[p]] -= getSymm(Lx[p]) * yi;
130 d -= getDiag(l_ki * getSymm(yi));
131 Li[p] = k; /* store L(k,i) in column form of L */
132 Lx[p] = l_ki;
133 ++m_nonZerosPerCol[i]; /* increment count of nonzeros in col i */
134 }
135 if (DoLDLT) {
136 m_diag[k] = d;
137 if (d == RealScalar(0)) {
138 ok = false; /* failure, D(k,k) is zero */
139 break;
140 }
141 } else {
142 Index p = Lp[k] + m_nonZerosPerCol[k]++;
143 Li[p] = k; /* store L(k,k) = sqrt (d) in column k */
144 if (NonHermitian ? d == RealScalar(0) : numext::real(d) <= RealScalar(0)) {
145 ok = false; /* failure, matrix is not positive definite */
146 break;
147 }
148 Lx[p] = sqrt(d);
149 }
150 }
151
152 m_info = ok ? Success : NumericalIssue;
153 m_factorizationIsOk = true;
154}
155
156} // end namespace Eigen
157
158#endif // EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
@ NumericalIssue
Definition Constants.h:442
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)