Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/SparseCholesky/SimplicialCholesky_impl.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008-2012 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 /*
0011 NOTE: these functions have been adapted from the LDL library:
0012 
0013 LDL Copyright (c) 2005 by Timothy A. Davis.  All Rights Reserved.
0014 
0015 The author of LDL, Timothy A. Davis., has executed a license with Google LLC
0016 to permit distribution of this code and derivative works as part of Eigen under
0017 the Mozilla Public License v. 2.0, as stated at the top of this file.
0018  */
0019 
0020 #ifndef EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
0021 #define EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
0022 
0023 namespace Eigen {
0024 
0025 template<typename Derived>
0026 void SimplicialCholeskyBase<Derived>::analyzePattern_preordered(const CholMatrixType& ap, bool doLDLT)
0027 {
0028   const StorageIndex size = StorageIndex(ap.rows());
0029   m_matrix.resize(size, size);
0030   m_parent.resize(size);
0031   m_nonZerosPerCol.resize(size);
0032 
0033   ei_declare_aligned_stack_constructed_variable(StorageIndex, tags, size, 0);
0034 
0035   for(StorageIndex k = 0; k < size; ++k)
0036   {
0037     /* L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) */
0038     m_parent[k] = -1;             /* parent of k is not yet known */
0039     tags[k] = k;                  /* mark node k as visited */
0040     m_nonZerosPerCol[k] = 0;      /* count of nonzeros in column k of L */
0041     for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
0042     {
0043       StorageIndex i = it.index();
0044       if(i < k)
0045       {
0046         /* follow path from i to root of etree, stop at flagged node */
0047         for(; tags[i] != k; i = m_parent[i])
0048         {
0049           /* find parent of i if not yet determined */
0050           if (m_parent[i] == -1)
0051             m_parent[i] = k;
0052           m_nonZerosPerCol[i]++;        /* L (k,i) is nonzero */
0053           tags[i] = k;                  /* mark i as visited */
0054         }
0055       }
0056     }
0057   }
0058 
0059   /* construct Lp index array from m_nonZerosPerCol column counts */
0060   StorageIndex* Lp = m_matrix.outerIndexPtr();
0061   Lp[0] = 0;
0062   for(StorageIndex k = 0; k < size; ++k)
0063     Lp[k+1] = Lp[k] + m_nonZerosPerCol[k] + (doLDLT ? 0 : 1);
0064 
0065   m_matrix.resizeNonZeros(Lp[size]);
0066 
0067   m_isInitialized     = true;
0068   m_info              = Success;
0069   m_analysisIsOk      = true;
0070   m_factorizationIsOk = false;
0071 }
0072 
0073 
0074 template<typename Derived>
0075 template<bool DoLDLT>
0076 void SimplicialCholeskyBase<Derived>::factorize_preordered(const CholMatrixType& ap)
0077 {
0078   using std::sqrt;
0079 
0080   eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
0081   eigen_assert(ap.rows()==ap.cols());
0082   eigen_assert(m_parent.size()==ap.rows());
0083   eigen_assert(m_nonZerosPerCol.size()==ap.rows());
0084 
0085   const StorageIndex size = StorageIndex(ap.rows());
0086   const StorageIndex* Lp = m_matrix.outerIndexPtr();
0087   StorageIndex* Li = m_matrix.innerIndexPtr();
0088   Scalar* Lx = m_matrix.valuePtr();
0089 
0090   ei_declare_aligned_stack_constructed_variable(Scalar, y, size, 0);
0091   ei_declare_aligned_stack_constructed_variable(StorageIndex,  pattern, size, 0);
0092   ei_declare_aligned_stack_constructed_variable(StorageIndex,  tags, size, 0);
0093 
0094   bool ok = true;
0095   m_diag.resize(DoLDLT ? size : 0);
0096 
0097   for(StorageIndex k = 0; k < size; ++k)
0098   {
0099     // compute nonzero pattern of kth row of L, in topological order
0100     y[k] = Scalar(0);                     // Y(0:k) is now all zero
0101     StorageIndex top = size;               // stack for pattern is empty
0102     tags[k] = k;                    // mark node k as visited
0103     m_nonZerosPerCol[k] = 0;        // count of nonzeros in column k of L
0104     for(typename CholMatrixType::InnerIterator it(ap,k); it; ++it)
0105     {
0106       StorageIndex i = it.index();
0107       if(i <= k)
0108       {
0109         y[i] += numext::conj(it.value());            /* scatter A(i,k) into Y (sum duplicates) */
0110         Index len;
0111         for(len = 0; tags[i] != k; i = m_parent[i])
0112         {
0113           pattern[len++] = i;     /* L(k,i) is nonzero */
0114           tags[i] = k;            /* mark i as visited */
0115         }
0116         while(len > 0)
0117           pattern[--top] = pattern[--len];
0118       }
0119     }
0120 
0121     /* compute numerical values kth row of L (a sparse triangular solve) */
0122 
0123     RealScalar d = numext::real(y[k]) * m_shiftScale + m_shiftOffset;    // get D(k,k), apply the shift function, and clear Y(k)
0124     y[k] = Scalar(0);
0125     for(; top < size; ++top)
0126     {
0127       Index i = pattern[top];       /* pattern[top:n-1] is pattern of L(:,k) */
0128       Scalar yi = y[i];             /* get and clear Y(i) */
0129       y[i] = Scalar(0);
0130 
0131       /* the nonzero entry L(k,i) */
0132       Scalar l_ki;
0133       if(DoLDLT)
0134         l_ki = yi / numext::real(m_diag[i]);
0135       else
0136         yi = l_ki = yi / Lx[Lp[i]];
0137 
0138       Index p2 = Lp[i] + m_nonZerosPerCol[i];
0139       Index p;
0140       for(p = Lp[i] + (DoLDLT ? 0 : 1); p < p2; ++p)
0141         y[Li[p]] -= numext::conj(Lx[p]) * yi;
0142       d -= numext::real(l_ki * numext::conj(yi));
0143       Li[p] = k;                          /* store L(k,i) in column form of L */
0144       Lx[p] = l_ki;
0145       ++m_nonZerosPerCol[i];              /* increment count of nonzeros in col i */
0146     }
0147     if(DoLDLT)
0148     {
0149       m_diag[k] = d;
0150       if(d == RealScalar(0))
0151       {
0152         ok = false;                         /* failure, D(k,k) is zero */
0153         break;
0154       }
0155     }
0156     else
0157     {
0158       Index p = Lp[k] + m_nonZerosPerCol[k]++;
0159       Li[p] = k ;                /* store L(k,k) = sqrt (d) in column k */
0160       if(d <= RealScalar(0)) {
0161         ok = false;              /* failure, matrix is not positive definite */
0162         break;
0163       }
0164       Lx[p] = sqrt(d) ;
0165     }
0166   }
0167 
0168   m_info = ok ? Success : NumericalIssue;
0169   m_factorizationIsOk = true;
0170 }
0171 
0172 } // end namespace Eigen
0173 
0174 #endif // EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H