SimplicialCholesky_impl.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
24 
25 namespace Eigen {
26 
27 template<typename Derived>
29 {
30  const StorageIndex size = StorageIndex(ap.rows());
31  m_matrix.resize(size, size);
32  m_parent.resize(size);
33  m_nonZerosPerCol.resize(size);
34 
36 
37  for(StorageIndex k = 0; k < size; ++k)
38  {
39  /* L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) */
40  m_parent[k] = -1; /* parent of k is not yet known */
41  tags[k] = k; /* mark node k as visited */
42  m_nonZerosPerCol[k] = 0; /* count of nonzeros in column k of L */
43  for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
44  {
45  StorageIndex i = it.index();
46  if(i < k)
47  {
48  /* follow path from i to root of etree, stop at flagged node */
49  for(; tags[i] != k; i = m_parent[i])
50  {
51  /* find parent of i if not yet determined */
52  if (m_parent[i] == -1)
53  m_parent[i] = k;
54  m_nonZerosPerCol[i]++; /* L (k,i) is nonzero */
55  tags[i] = k; /* mark i as visited */
56  }
57  }
58  }
59  }
60 
61  /* construct Lp index array from m_nonZerosPerCol column counts */
62  StorageIndex* Lp = m_matrix.outerIndexPtr();
63  Lp[0] = 0;
64  for(StorageIndex k = 0; k < size; ++k)
65  Lp[k+1] = Lp[k] + m_nonZerosPerCol[k] + (doLDLT ? 0 : 1);
66 
67  m_matrix.resizeNonZeros(Lp[size]);
68 
69  m_isInitialized = true;
70  m_info = Success;
71  m_analysisIsOk = true;
72  m_factorizationIsOk = false;
73 }
74 
75 
76 template<typename Derived>
77 template<bool DoLDLT>
79 {
80  using std::sqrt;
81 
82  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
83  eigen_assert(ap.rows()==ap.cols());
84  eigen_assert(m_parent.size()==ap.rows());
85  eigen_assert(m_nonZerosPerCol.size()==ap.rows());
86 
87  const StorageIndex size = StorageIndex(ap.rows());
88  const StorageIndex* Lp = m_matrix.outerIndexPtr();
89  StorageIndex* Li = m_matrix.innerIndexPtr();
90  Scalar* Lx = m_matrix.valuePtr();
91 
95 
96  bool ok = true;
97  m_diag.resize(DoLDLT ? size : 0);
98 
99  for(StorageIndex k = 0; k < size; ++k)
100  {
101  // compute nonzero pattern of kth row of L, in topological order
102  y[k] = Scalar(0); // Y(0:k) is now all zero
103  StorageIndex top = size; // stack for pattern is empty
104  tags[k] = k; // mark node k as visited
105  m_nonZerosPerCol[k] = 0; // count of nonzeros in column k of L
106  for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
107  {
108  StorageIndex i = it.index();
109  if(i <= k)
110  {
111  y[i] += numext::conj(it.value()); /* scatter A(i,k) into Y (sum duplicates) */
112  Index len;
113  for(len = 0; tags[i] != k; i = m_parent[i])
114  {
115  pattern[len++] = i; /* L(k,i) is nonzero */
116  tags[i] = k; /* mark i as visited */
117  }
118  while(len > 0)
119  pattern[--top] = pattern[--len];
120  }
121  }
122 
123  /* compute numerical values kth row of L (a sparse triangular solve) */
124 
125  RealScalar d = numext::real(y[k]) * m_shiftScale + m_shiftOffset; // get D(k,k), apply the shift function, and clear Y(k)
126  y[k] = Scalar(0);
127  for(; top < size; ++top)
128  {
129  Index i = pattern[top]; /* pattern[top:n-1] is pattern of L(:,k) */
130  Scalar yi = y[i]; /* get and clear Y(i) */
131  y[i] = Scalar(0);
132 
133  /* the nonzero entry L(k,i) */
134  Scalar l_ki;
135  if(DoLDLT)
136  l_ki = yi / numext::real(m_diag[i]);
137  else
138  yi = l_ki = yi / Lx[Lp[i]];
139 
140  Index p2 = Lp[i] + m_nonZerosPerCol[i];
141  Index p;
142  for(p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p)
143  y[Li[p]] -= numext::conj(Lx[p]) * yi;
144  d -= numext::real(l_ki * numext::conj(yi));
145  Li[p] = k; /* store L(k,i) in column form of L */
146  Lx[p] = l_ki;
147  ++m_nonZerosPerCol[i]; /* increment count of nonzeros in col i */
148  }
149  if(DoLDLT)
150  {
151  m_diag[k] = d;
152  if(d == RealScalar(0))
153  {
154  ok = false; /* failure, D(k,k) is zero */
155  break;
156  }
157  }
158  else
159  {
160  Index p = Lp[k] + m_nonZerosPerCol[k]++;
161  Li[p] = k ; /* store L(k,k) = sqrt (d) in column k */
162  if(d <= RealScalar(0)) {
163  ok = false; /* failure, matrix is not positive definite */
164  break;
165  }
166  Lx[p] = sqrt(d) ;
167  }
168  }
169 
170  m_info = ok ? Success : NumericalIssue;
171  m_factorizationIsOk = true;
172 }
173 
174 } // end namespace Eigen
175 
176 #endif // EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
const SqrtReturnType sqrt() const
RealReturnType real() const
#define eigen_assert(x)
Definition: Macros.h:902
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:847
float * p
MatrixType::RealScalar RealScalar
void factorize_preordered(const CholMatrixType &a)
MatrixType::StorageIndex StorageIndex
void analyzePattern_preordered(const CholMatrixType &a, bool doLDLT)
Index cols() const
Definition: SparseMatrix.h:167
Index rows() const
Definition: SparseMatrix.h:165
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
const Scalar & y
: InteropHeaders
Definition: Core:139
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)