IncompleteLUT.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_INCOMPLETE_LUT_H
00011 #define EIGEN_INCOMPLETE_LUT_H
00012 
00013 namespace Eigen { 
00014 
00043 template <typename _Scalar>
00044 class IncompleteLUT : internal::noncopyable
00045 {
00046     typedef _Scalar Scalar;
00047     typedef typename NumTraits<Scalar>::Real RealScalar;
00048     typedef Matrix<Scalar,Dynamic,1> Vector;
00049     typedef SparseMatrix<Scalar,RowMajor> FactorType;
00050     typedef SparseMatrix<Scalar,ColMajor> PermutType;
00051     typedef typename FactorType::Index Index;
00052 
00053   public:
00054     typedef Matrix<Scalar,Dynamic,Dynamic> MatrixType;
00055     
00056     IncompleteLUT()
00057       : m_droptol(NumTraits<Scalar>::dummy_precision()), m_fillfactor(10),
00058         m_analysisIsOk(false), m_factorizationIsOk(false), m_isInitialized(false)
00059     {}
00060     
00061     template<typename MatrixType>
00062     IncompleteLUT(const MatrixType& mat, RealScalar droptol=NumTraits<Scalar>::dummy_precision(), int fillfactor = 10)
00063       : m_droptol(droptol),m_fillfactor(fillfactor),
00064         m_analysisIsOk(false),m_factorizationIsOk(false),m_isInitialized(false)
00065     {
00066       eigen_assert(fillfactor != 0);
00067       compute(mat); 
00068     }
00069     
00070     Index rows() const { return m_lu.rows(); }
00071     
00072     Index cols() const { return m_lu.cols(); }
00073 
00079     ComputationInfo info() const
00080     {
00081       eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
00082       return m_info;
00083     }
00084     
00085     template<typename MatrixType>
00086     void analyzePattern(const MatrixType& amat);
00087     
00088     template<typename MatrixType>
00089     void factorize(const MatrixType& amat);
00090     
00096     template<typename MatrixType>
00097     IncompleteLUT<Scalar>& compute(const MatrixType& amat)
00098     {
00099       analyzePattern(amat); 
00100       factorize(amat);
00101       eigen_assert(m_factorizationIsOk == true); 
00102       m_isInitialized = true;
00103       return *this;
00104     }
00105 
00106     void setDroptol(RealScalar droptol); 
00107     void setFillfactor(int fillfactor); 
00108     
00109     template<typename Rhs, typename Dest>
00110     void _solve(const Rhs& b, Dest& x) const
00111     {
00112       x = m_Pinv * b;  
00113       x = m_lu.template triangularView<UnitLower>().solve(x);
00114       x = m_lu.template triangularView<Upper>().solve(x);
00115       x = m_P * x; 
00116     }
00117 
00118     template<typename Rhs> inline const internal::solve_retval<IncompleteLUT, Rhs>
00119      solve(const MatrixBase<Rhs>& b) const
00120     {
00121       eigen_assert(m_isInitialized && "IncompleteLUT is not initialized.");
00122       eigen_assert(cols()==b.rows()
00123                 && "IncompleteLUT::solve(): invalid number of rows of the right hand side matrix b");
00124       return internal::solve_retval<IncompleteLUT, Rhs>(*this, b.derived());
00125     }
00126 
00127 protected:
00128 
00129     template <typename VectorV, typename VectorI>
00130     int QuickSplit(VectorV &row, VectorI &ind, int ncut);
00131 
00132 
00134     struct keep_diag {
00135       inline bool operator() (const Index& row, const Index& col, const Scalar&) const
00136       {
00137         return row!=col;
00138       }
00139     };
00140 
00141 protected:
00142 
00143     FactorType m_lu;
00144     RealScalar m_droptol;
00145     int m_fillfactor;
00146     bool m_analysisIsOk;
00147     bool m_factorizationIsOk;
00148     bool m_isInitialized;
00149     ComputationInfo m_info;
00150     PermutationMatrix<Dynamic,Dynamic,Index> m_P;     // Fill-reducing permutation
00151     PermutationMatrix<Dynamic,Dynamic,Index> m_Pinv;  // Inverse permutation
00152 };
00153 
00158 template<typename Scalar>
00159 void IncompleteLUT<Scalar>::setDroptol(RealScalar droptol)
00160 {
00161   this->m_droptol = droptol;   
00162 }
00163 
00168 template<typename Scalar>
00169 void IncompleteLUT<Scalar>::setFillfactor(int fillfactor)
00170 {
00171   this->m_fillfactor = fillfactor;   
00172 }
00173 
00174 
00184 template <typename Scalar>
00185 template <typename VectorV, typename VectorI>
00186 int IncompleteLUT<Scalar>::QuickSplit(VectorV &row, VectorI &ind, int ncut)
00187 {
00188   using std::swap;
00189   int mid;
00190   int n = row.size(); /* length of the vector */
00191   int first, last ; 
00192   
00193   ncut--; /* to fit the zero-based indices */
00194   first = 0; 
00195   last = n-1; 
00196   if (ncut < first || ncut > last ) return 0;
00197   
00198   do {
00199     mid = first; 
00200     RealScalar abskey = std::abs(row(mid)); 
00201     for (int j = first + 1; j <= last; j++) {
00202       if ( std::abs(row(j)) > abskey) {
00203         ++mid;
00204         swap(row(mid), row(j));
00205         swap(ind(mid), ind(j));
00206       }
00207     }
00208     /* Interchange for the pivot element */
00209     swap(row(mid), row(first));
00210     swap(ind(mid), ind(first));
00211     
00212     if (mid > ncut) last = mid - 1;
00213     else if (mid < ncut ) first = mid + 1; 
00214   } while (mid != ncut );
00215   
00216   return 0; /* mid is equal to ncut */ 
00217 }
00218 
00219 template <typename Scalar>
00220 template<typename _MatrixType>
00221 void IncompleteLUT<Scalar>::analyzePattern(const _MatrixType& amat)
00222 {
00223   // Compute the Fill-reducing permutation
00224   SparseMatrix<Scalar,ColMajor, Index> mat1 = amat;
00225   SparseMatrix<Scalar,ColMajor, Index> mat2 = amat.transpose();
00226   // Symmetrize the pattern
00227   // FIXME for a matrix with nearly symmetric pattern, mat2+mat1 is the appropriate choice.
00228   //       on the other hand for a really non-symmetric pattern, mat2*mat1 should be prefered...
00229   SparseMatrix<Scalar,ColMajor, Index> AtA = mat2 + mat1;
00230   AtA.prune(keep_diag());
00231   internal::minimum_degree_ordering<Scalar, Index>(AtA, m_P);  // Then compute the AMD ordering...
00232 
00233   m_Pinv  = m_P.inverse(); // ... and the inverse permutation
00234 
00235   m_analysisIsOk = true;
00236 }
00237 
00238 template <typename Scalar>
00239 template<typename _MatrixType>
00240 void IncompleteLUT<Scalar>::factorize(const _MatrixType& amat)
00241 {
00242   using std::sqrt;
00243   using std::swap;
00244   using std::abs;
00245 
00246   eigen_assert((amat.rows() == amat.cols()) && "The factorization should be done on a square matrix");
00247   int n = amat.cols();  // Size of the matrix
00248   m_lu.resize(n,n);
00249   // Declare Working vectors and variables
00250   Vector u(n) ;     // real values of the row -- maximum size is n --
00251   VectorXi ju(n);   // column position of the values in u -- maximum size  is n
00252   VectorXi jr(n);   // Indicate the position of the nonzero elements in the vector u -- A zero location is indicated by -1
00253 
00254   // Apply the fill-reducing permutation
00255   eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
00256   SparseMatrix<Scalar,RowMajor, Index> mat;
00257   mat = amat.twistedBy(m_Pinv);
00258 
00259   // Initialization
00260   jr.fill(-1);
00261   ju.fill(0);
00262   u.fill(0);
00263 
00264   // number of largest elements to keep in each row:
00265   int fill_in =   static_cast<int> (amat.nonZeros()*m_fillfactor)/n+1;
00266   if (fill_in > n) fill_in = n;
00267 
00268   // number of largest nonzero elements to keep in the L and the U part of the current row:
00269   int nnzL = fill_in/2;
00270   int nnzU = nnzL;
00271   m_lu.reserve(n * (nnzL + nnzU + 1));
00272 
00273   // global loop over the rows of the sparse matrix
00274   for (int ii = 0; ii < n; ii++)
00275   {
00276     // 1 - copy the lower and the upper part of the row i of mat in the working vector u
00277 
00278     int sizeu = 1; // number of nonzero elements in the upper part of the current row
00279     int sizel = 0; // number of nonzero elements in the lower part of the current row
00280     ju(ii)    = ii;
00281     u(ii)     = 0;
00282     jr(ii)    = ii;
00283     RealScalar rownorm = 0;
00284 
00285     typename FactorType::InnerIterator j_it(mat, ii); // Iterate through the current row ii
00286     for (; j_it; ++j_it)
00287     {
00288       int k = j_it.index();
00289       if (k < ii)
00290       {
00291         // copy the lower part
00292         ju(sizel) = k;
00293         u(sizel) = j_it.value();
00294         jr(k) = sizel;
00295         ++sizel;
00296       }
00297       else if (k == ii)
00298       {
00299         u(ii) = j_it.value();
00300       }
00301       else
00302       {
00303         // copy the upper part
00304         int jpos = ii + sizeu;
00305         ju(jpos) = k;
00306         u(jpos) = j_it.value();
00307         jr(k) = jpos;
00308         ++sizeu;
00309       }
00310       rownorm += internal::abs2(j_it.value());
00311     }
00312 
00313     // 2 - detect possible zero row
00314     if(rownorm==0)
00315     {
00316       m_info = NumericalIssue;
00317       return;
00318     }
00319     // Take the 2-norm of the current row as a relative tolerance
00320     rownorm = sqrt(rownorm);
00321 
00322     // 3 - eliminate the previous nonzero rows
00323     int jj = 0;
00324     int len = 0;
00325     while (jj < sizel)
00326     {
00327       // In order to eliminate in the correct order,
00328       // we must select first the smallest column index among  ju(jj:sizel)
00329       int k;
00330       int minrow = ju.segment(jj,sizel-jj).minCoeff(&k); // k is relative to the segment
00331       k += jj;
00332       if (minrow != ju(jj))
00333       {
00334         // swap the two locations
00335         int j = ju(jj);
00336         swap(ju(jj), ju(k));
00337         jr(minrow) = jj;   jr(j) = k;
00338         swap(u(jj), u(k));
00339       }
00340       // Reset this location
00341       jr(minrow) = -1;
00342 
00343       // Start elimination
00344       typename FactorType::InnerIterator ki_it(m_lu, minrow);
00345       while (ki_it && ki_it.index() < minrow) ++ki_it;
00346       eigen_internal_assert(ki_it && ki_it.col()==minrow);
00347       Scalar fact = u(jj) / ki_it.value();
00348 
00349       // drop too small elements
00350       if(abs(fact) <= m_droptol)
00351       {
00352         jj++;
00353         continue;
00354       }
00355 
00356       // linear combination of the current row ii and the row minrow
00357       ++ki_it;
00358       for (; ki_it; ++ki_it)
00359       {
00360         Scalar prod = fact * ki_it.value();
00361         int j       = ki_it.index();
00362         int jpos    = jr(j);
00363         if (jpos == -1) // fill-in element
00364         {
00365           int newpos;
00366           if (j >= ii) // dealing with the upper part
00367           {
00368             newpos = ii + sizeu;
00369             sizeu++;
00370             eigen_internal_assert(sizeu<=n);
00371           }
00372           else // dealing with the lower part
00373           {
00374             newpos = sizel;
00375             sizel++;
00376             eigen_internal_assert(sizel<=ii);
00377           }
00378           ju(newpos) = j;
00379           u(newpos) = -prod;
00380           jr(j) = newpos;
00381         }
00382         else
00383           u(jpos) -= prod;
00384       }
00385       // store the pivot element
00386       u(len) = fact;
00387       ju(len) = minrow;
00388       ++len;
00389 
00390       jj++;
00391     } // end of the elimination on the row ii
00392 
00393     // reset the upper part of the pointer jr to zero
00394     for(int k = 0; k <sizeu; k++) jr(ju(ii+k)) = -1;
00395 
00396     // 4 - partially sort and insert the elements in the m_lu matrix
00397 
00398     // sort the L-part of the row
00399     sizel = len;
00400     len = (std::min)(sizel, nnzL);
00401     typename Vector::SegmentReturnType ul(u.segment(0, sizel));
00402     typename VectorXi::SegmentReturnType jul(ju.segment(0, sizel));
00403     QuickSplit(ul, jul, len);
00404 
00405     // store the largest m_fill elements of the L part
00406     m_lu.startVec(ii);
00407     for(int k = 0; k < len; k++)
00408       m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
00409 
00410     // store the diagonal element
00411     // apply a shifting rule to avoid zero pivots (we are doing an incomplete factorization)
00412     if (u(ii) == Scalar(0))
00413       u(ii) = sqrt(m_droptol) * rownorm;
00414     m_lu.insertBackByOuterInnerUnordered(ii, ii) = u(ii);
00415 
00416     // sort the U-part of the row
00417     // apply the dropping rule first
00418     len = 0;
00419     for(int k = 1; k < sizeu; k++)
00420     {
00421       if(abs(u(ii+k)) > m_droptol * rownorm )
00422       {
00423         ++len;
00424         u(ii + len)  = u(ii + k);
00425         ju(ii + len) = ju(ii + k);
00426       }
00427     }
00428     sizeu = len + 1; // +1 to take into account the diagonal element
00429     len = (std::min)(sizeu, nnzU);
00430     typename Vector::SegmentReturnType uu(u.segment(ii+1, sizeu-1));
00431     typename VectorXi::SegmentReturnType juu(ju.segment(ii+1, sizeu-1));
00432     QuickSplit(uu, juu, len);
00433 
00434     // store the largest elements of the U part
00435     for(int k = ii + 1; k < ii + len; k++)
00436       m_lu.insertBackByOuterInnerUnordered(ii,ju(k)) = u(k);
00437   }
00438 
00439   m_lu.finalize();
00440   m_lu.makeCompressed();
00441 
00442   m_factorizationIsOk = true;
00443   m_info = Success;
00444 }
00445 
00446 namespace internal {
00447 
00448 template<typename _MatrixType, typename Rhs>
00449 struct solve_retval<IncompleteLUT<_MatrixType>, Rhs>
00450   : solve_retval_base<IncompleteLUT<_MatrixType>, Rhs>
00451 {
00452   typedef IncompleteLUT<_MatrixType> Dec;
00453   EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
00454 
00455   template<typename Dest> void evalTo(Dest& dst) const
00456   {
00457     dec()._solve(rhs(),dst);
00458   }
00459 };
00460 
00461 } // end namespace internal
00462 
00463 } // end namespace Eigen
00464 
00465 #endif // EIGEN_INCOMPLETE_LUT_H
00466 


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:10:49