xref: /aosp_15_r20/external/eigen/Eigen/src/SparseCholesky/SimplicialCholesky_impl.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008-2012 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li /*
11*bf2c3715SXin Li NOTE: these functions have been adapted from the LDL library:
12*bf2c3715SXin Li 
13*bf2c3715SXin Li LDL Copyright (c) 2005 by Timothy A. Davis.  All Rights Reserved.
14*bf2c3715SXin Li 
15*bf2c3715SXin Li The author of LDL, Timothy A. Davis., has executed a license with Google LLC
16*bf2c3715SXin Li to permit distribution of this code and derivative works as part of Eigen under
17*bf2c3715SXin Li the Mozilla Public License v. 2.0, as stated at the top of this file.
18*bf2c3715SXin Li  */
19*bf2c3715SXin Li 
20*bf2c3715SXin Li #ifndef EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
21*bf2c3715SXin Li #define EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
22*bf2c3715SXin Li 
23*bf2c3715SXin Li namespace Eigen {
24*bf2c3715SXin Li 
25*bf2c3715SXin Li template<typename Derived>
analyzePattern_preordered(const CholMatrixType & ap,bool doLDLT)26*bf2c3715SXin Li void SimplicialCholeskyBase<Derived>::analyzePattern_preordered(const CholMatrixType& ap, bool doLDLT)
27*bf2c3715SXin Li {
28*bf2c3715SXin Li   const StorageIndex size = StorageIndex(ap.rows());
29*bf2c3715SXin Li   m_matrix.resize(size, size);
30*bf2c3715SXin Li   m_parent.resize(size);
31*bf2c3715SXin Li   m_nonZerosPerCol.resize(size);
32*bf2c3715SXin Li 
33*bf2c3715SXin Li   ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
34*bf2c3715SXin Li 
35*bf2c3715SXin Li   for(StorageIndex k = 0; k < size; ++k)
36*bf2c3715SXin Li   {
37*bf2c3715SXin Li     /* L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) */
38*bf2c3715SXin Li     m_parent[k] = -1;             /* parent of k is not yet known */
39*bf2c3715SXin Li     tags[k] = k;                  /* mark node k as visited */
40*bf2c3715SXin Li     m_nonZerosPerCol[k] = 0;      /* count of nonzeros in column k of L */
41*bf2c3715SXin Li     for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
42*bf2c3715SXin Li     {
43*bf2c3715SXin Li       StorageIndex i = it.index();
44*bf2c3715SXin Li       if(i < k)
45*bf2c3715SXin Li       {
46*bf2c3715SXin Li         /* follow path from i to root of etree, stop at flagged node */
47*bf2c3715SXin Li         for(; tags[i] != k; i = m_parent[i])
48*bf2c3715SXin Li         {
49*bf2c3715SXin Li           /* find parent of i if not yet determined */
50*bf2c3715SXin Li           if (m_parent[i] == -1)
51*bf2c3715SXin Li             m_parent[i] = k;
52*bf2c3715SXin Li           m_nonZerosPerCol[i]++;        /* L (k,i) is nonzero */
53*bf2c3715SXin Li           tags[i] = k;                  /* mark i as visited */
54*bf2c3715SXin Li         }
55*bf2c3715SXin Li       }
56*bf2c3715SXin Li     }
57*bf2c3715SXin Li   }
58*bf2c3715SXin Li 
59*bf2c3715SXin Li   /* construct Lp index array from m_nonZerosPerCol column counts */
60*bf2c3715SXin Li   StorageIndex* Lp = m_matrix.outerIndexPtr();
61*bf2c3715SXin Li   Lp[0] = 0;
62*bf2c3715SXin Li   for(StorageIndex k = 0; k < size; ++k)
63*bf2c3715SXin Li     Lp[k+1] = Lp[k] + m_nonZerosPerCol[k] + (doLDLT ? 0 : 1);
64*bf2c3715SXin Li 
65*bf2c3715SXin Li   m_matrix.resizeNonZeros(Lp[size]);
66*bf2c3715SXin Li 
67*bf2c3715SXin Li   m_isInitialized     = true;
68*bf2c3715SXin Li   m_info              = Success;
69*bf2c3715SXin Li   m_analysisIsOk      = true;
70*bf2c3715SXin Li   m_factorizationIsOk = false;
71*bf2c3715SXin Li }
72*bf2c3715SXin Li 
73*bf2c3715SXin Li 
74*bf2c3715SXin Li template<typename Derived>
75*bf2c3715SXin Li template<bool DoLDLT>
factorize_preordered(const CholMatrixType & ap)76*bf2c3715SXin Li void SimplicialCholeskyBase<Derived>::factorize_preordered(const CholMatrixType& ap)
77*bf2c3715SXin Li {
78*bf2c3715SXin Li   using std::sqrt;
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
81*bf2c3715SXin Li   eigen_assert(ap.rows()==ap.cols());
82*bf2c3715SXin Li   eigen_assert(m_parent.size()==ap.rows());
83*bf2c3715SXin Li   eigen_assert(m_nonZerosPerCol.size()==ap.rows());
84*bf2c3715SXin Li 
85*bf2c3715SXin Li   const StorageIndex size = StorageIndex(ap.rows());
86*bf2c3715SXin Li   const StorageIndex* Lp = m_matrix.outerIndexPtr();
87*bf2c3715SXin Li   StorageIndex* Li = m_matrix.innerIndexPtr();
88*bf2c3715SXin Li   Scalar* Lx = m_matrix.valuePtr();
89*bf2c3715SXin Li 
90*bf2c3715SXin Li   ei_declare_aligned_stack_constructed_variable(Scalar, y, size, 0);
91*bf2c3715SXin Li   ei_declare_aligned_stack_constructed_variable(StorageIndex,  pattern, size, 0);
92*bf2c3715SXin Li   ei_declare_aligned_stack_constructed_variable(StorageIndex,  tags, size, 0);
93*bf2c3715SXin Li 
94*bf2c3715SXin Li   bool ok = true;
95*bf2c3715SXin Li   m_diag.resize(DoLDLT ? size : 0);
96*bf2c3715SXin Li 
97*bf2c3715SXin Li   for(StorageIndex k = 0; k < size; ++k)
98*bf2c3715SXin Li   {
99*bf2c3715SXin Li     // compute nonzero pattern of kth row of L, in topological order
100*bf2c3715SXin Li     y[k] = Scalar(0);                     // Y(0:k) is now all zero
101*bf2c3715SXin Li     StorageIndex top = size;               // stack for pattern is empty
102*bf2c3715SXin Li     tags[k] = k;                    // mark node k as visited
103*bf2c3715SXin Li     m_nonZerosPerCol[k] = 0;        // count of nonzeros in column k of L
104*bf2c3715SXin Li     for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
105*bf2c3715SXin Li     {
106*bf2c3715SXin Li       StorageIndex i = it.index();
107*bf2c3715SXin Li       if(i <= k)
108*bf2c3715SXin Li       {
109*bf2c3715SXin Li         y[i] += numext::conj(it.value());            /* scatter A(i,k) into Y (sum duplicates) */
110*bf2c3715SXin Li         Index len;
111*bf2c3715SXin Li         for(len = 0; tags[i] != k; i = m_parent[i])
112*bf2c3715SXin Li         {
113*bf2c3715SXin Li           pattern[len++] = i;     /* L(k,i) is nonzero */
114*bf2c3715SXin Li           tags[i] = k;            /* mark i as visited */
115*bf2c3715SXin Li         }
116*bf2c3715SXin Li         while(len > 0)
117*bf2c3715SXin Li           pattern[--top] = pattern[--len];
118*bf2c3715SXin Li       }
119*bf2c3715SXin Li     }
120*bf2c3715SXin Li 
121*bf2c3715SXin Li     /* compute numerical values kth row of L (a sparse triangular solve) */
122*bf2c3715SXin Li 
123*bf2c3715SXin Li     RealScalar d = numext::real(y[k]) * m_shiftScale + m_shiftOffset;    // get D(k,k), apply the shift function, and clear Y(k)
124*bf2c3715SXin Li     y[k] = Scalar(0);
125*bf2c3715SXin Li     for(; top < size; ++top)
126*bf2c3715SXin Li     {
127*bf2c3715SXin Li       Index i = pattern[top];       /* pattern[top:n-1] is pattern of L(:,k) */
128*bf2c3715SXin Li       Scalar yi = y[i];             /* get and clear Y(i) */
129*bf2c3715SXin Li       y[i] = Scalar(0);
130*bf2c3715SXin Li 
131*bf2c3715SXin Li       /* the nonzero entry L(k,i) */
132*bf2c3715SXin Li       Scalar l_ki;
133*bf2c3715SXin Li       if(DoLDLT)
134*bf2c3715SXin Li         l_ki = yi / numext::real(m_diag[i]);
135*bf2c3715SXin Li       else
136*bf2c3715SXin Li         yi = l_ki = yi / Lx[Lp[i]];
137*bf2c3715SXin Li 
138*bf2c3715SXin Li       Index p2 = Lp[i] + m_nonZerosPerCol[i];
139*bf2c3715SXin Li       Index p;
140*bf2c3715SXin Li       for(p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p)
141*bf2c3715SXin Li         y[Li[p]] -= numext::conj(Lx[p]) * yi;
142*bf2c3715SXin Li       d -= numext::real(l_ki * numext::conj(yi));
143*bf2c3715SXin Li       Li[p] = k;                          /* store L(k,i) in column form of L */
144*bf2c3715SXin Li       Lx[p] = l_ki;
145*bf2c3715SXin Li       ++m_nonZerosPerCol[i];              /* increment count of nonzeros in col i */
146*bf2c3715SXin Li     }
147*bf2c3715SXin Li     if(DoLDLT)
148*bf2c3715SXin Li     {
149*bf2c3715SXin Li       m_diag[k] = d;
150*bf2c3715SXin Li       if(d == RealScalar(0))
151*bf2c3715SXin Li       {
152*bf2c3715SXin Li         ok = false;                         /* failure, D(k,k) is zero */
153*bf2c3715SXin Li         break;
154*bf2c3715SXin Li       }
155*bf2c3715SXin Li     }
156*bf2c3715SXin Li     else
157*bf2c3715SXin Li     {
158*bf2c3715SXin Li       Index p = Lp[k] + m_nonZerosPerCol[k]++;
159*bf2c3715SXin Li       Li[p] = k ;                /* store L(k,k) = sqrt (d) in column k */
160*bf2c3715SXin Li       if(d <= RealScalar(0)) {
161*bf2c3715SXin Li         ok = false;              /* failure, matrix is not positive definite */
162*bf2c3715SXin Li         break;
163*bf2c3715SXin Li       }
164*bf2c3715SXin Li       Lx[p] = sqrt(d) ;
165*bf2c3715SXin Li     }
166*bf2c3715SXin Li   }
167*bf2c3715SXin Li 
168*bf2c3715SXin Li   m_info = ok ? Success : NumericalIssue;
169*bf2c3715SXin Li   m_factorizationIsOk = true;
170*bf2c3715SXin Li }
171*bf2c3715SXin Li 
172*bf2c3715SXin Li } // end namespace Eigen
173*bf2c3715SXin Li 
174*bf2c3715SXin Li #endif // EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
175