ilu.hh

Go to the documentation of this file.
00001 #ifndef DUNE_ILU_HH
00002 #define DUNE_ILU_HH
00003 
00004 #include<cmath>
00005 #include<complex>
00006 #include<iostream>
00007 #include<iomanip>
00008 #include<string>
00009 #include<set>
00010 #include<map>
00011 
00012 #include "istlexception.hh"
00013 #include "io.hh"
00014 
00019 namespace Dune {
00020    
00025 
00026   template<class M>
00027   void bilu0_decomposition (M& A)
00028   {
00029         // iterator types
00030         typedef typename M::RowIterator rowiterator;
00031         typedef typename M::ColIterator coliterator;
00032         typedef typename M::block_type block;
00033 
00034         // implement left looking variant with stored inverse
00035         rowiterator endi=A.end();
00036         for (rowiterator i=A.begin(); i!=endi; ++i)
00037           {
00038                 // coliterator is diagonal after the following loop
00039                 coliterator endij=(*i).end(); // end of row i 
00040                 coliterator ij;
00041 
00042                 // eliminate entries left of diagonal; store L factor
00043                 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
00044                   {
00045                         // find A_jj which eliminates A_ij
00046                         coliterator jj = A[ij.index()].find(ij.index());
00047 
00048                         // compute L_ij = A_jj^-1 * A_ij
00049                         (*ij).rightmultiply(*jj);
00050 
00051                         // modify row
00052                         coliterator endjk=A[ij.index()].end(); // end of row j
00053                         coliterator jk=jj; ++jk;
00054                         coliterator ik=ij; ++ik;
00055                         while (ik!=endij && jk!=endjk)
00056                           if (ik.index()==jk.index())
00057                                 {
00058                                   block B(*jk);
00059                                   B.leftmultiply(*ij);
00060                                   *ik -= B;
00061                                   ++ik; ++jk;
00062                                 }
00063                           else
00064                                 {
00065                                   if (ik.index()<jk.index())
00066                                         ++ik;
00067                                   else
00068                                         ++jk;
00069                                 }
00070                   }
00071 
00072                 // invert pivot and store it in A
00073                 if (ij.index()!=i.index())
00074                   DUNE_THROW(ISTLError,"diagonal entry missing");
00075                 (*ij).invert(); // compute inverse of diagonal block
00076           }
00077   }
00078 
00080   template<class M, class X, class Y>
00081   void bilu_backsolve (const M& A, X& v, const Y& d)
00082   {
00083         // iterator types
00084         typedef typename M::ConstRowIterator rowiterator;
00085         typedef typename M::ConstColIterator coliterator;
00086         typedef typename Y::block_type dblock;
00087         typedef typename X::block_type vblock;
00088 
00089         // lower triangular solve
00090         rowiterator endi=A.end();
00091         for (rowiterator i=A.begin(); i!=endi; ++i)
00092           {
00093                 dblock rhs(d[i.index()]);
00094                 for (coliterator j=(*i).begin(); j.index()<i.index(); ++j)
00095                   (*j).mmv(v[j.index()],rhs);
00096                 v[i.index()] = rhs; // Lii = I
00097           }
00098 
00099         // upper triangular solve
00100         rowiterator rendi=A.rend();
00101         for (rowiterator i=A.rbegin(); i!=rendi; --i)
00102           {
00103                 vblock rhs(v[i.index()]);
00104                 coliterator j;
00105                 for (j=(*i).rbegin(); j.index()>i.index(); --j)
00106                         (*j).mmv(v[j.index()],rhs);
00107                 v[i.index()] = 0;
00108                 (*j).umv(rhs,v[i.index()]); // diagonal stores inverse!
00109           }
00110   }
00111 
00112 
00113 
00114   // recursive function template to access first entry of a matrix
00115   template<class M>
00116   typename M::field_type& firstmatrixelement (M& A)
00117   {
00118         return firstmatrixelement(*(A.begin()->begin()));
00119   }
00120 
00121   template<class K, int n, int m>
00122   K& firstmatrixelement (FieldMatrix<K,n,m>& A)
00123   {
00124         return A[0][0];
00125   }
00126 
00127   template<class K>
00128   K& firstmatrixelement (FieldMatrix<K,1,1>& A)
00129   {
00130         return A[0][0];
00131   }
00132 
00133 
00140   template<class M>
00141   void bilu_decomposition (const M& A, int n, M& ILU)
00142   {
00143         // iterator types
00144         typedef typename M::RowIterator rowiterator;
00145         typedef typename M::ColIterator coliterator;
00146         typedef typename M::ConstRowIterator crowiterator;
00147         typedef typename M::ConstColIterator ccoliterator;
00148         typedef typename M::CreateIterator createiterator;
00149         typedef typename M::block_type block;
00150         typedef typename M::field_type K;
00151         typedef typename std::map<int,int> map;
00152         typedef typename std::map<int,int>::iterator mapiterator;
00153 
00154         // symbolic factorization phase, store generation number in first matrix element
00155         crowiterator endi=A.end();
00156         createiterator ci=ILU.createbegin();
00157         for (crowiterator i=A.begin(); i!=endi; ++i)
00158           {
00159                 //              std::cout << "in row " << i.index() << std::endl;
00160                 map rowpattern; // maps column index to generation
00161 
00162                 // initialize pattern with row of A
00163                 for (ccoliterator j=(*i).begin(); j!=(*i).end(); ++j)
00164                   rowpattern[j.index()] = 0;
00165 
00166                 // eliminate entries in row which are to the left of the diagonal
00167                 for (mapiterator ik=rowpattern.begin(); (*ik).first<i.index(); ++ik)
00168                   {
00169                         if ((*ik).second<n)
00170                           {
00171 //                              std::cout << "  eliminating " << i.index() << "," << (*ik).first
00172 //                                                << " level " << (*ik).second << std::endl;
00173 
00174                                 coliterator endk = ILU[(*ik).first].end(); // end of row k
00175                                 coliterator kj = ILU[(*ik).first].find((*ik).first); // diagonal in k
00176                                 for (++kj; kj!=endk; ++kj) // row k eliminates in row i
00177                                   {
00178                                         int generation = (int) firstmatrixelement(*kj);
00179                                         if (generation<n)
00180                                           {
00181                                                 mapiterator ij = rowpattern.find(kj.index());
00182                                                 if (ij==rowpattern.end())
00183                                                   {
00184                                                         //std::cout << "    new entry " << i.index() << "," << kj.index()
00185 //                                                << " level " << (*ik).second+1 << std::endl;
00186 
00187                                                         rowpattern[kj.index()] = generation+1;
00188                                                   }
00189                                           }
00190                                   }
00191                           }
00192                   }
00193 
00194                 // create row
00195                 for (mapiterator ik=rowpattern.begin(); ik!=rowpattern.end(); ++ik)
00196                   ci.insert((*ik).first);
00197                 ++ci; // now row i exist
00198 
00199                 // write generation index into entries
00200                 coliterator endILUij = ILU[i.index()].end();;
00201                 for (coliterator ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
00202                   firstmatrixelement(*ILUij) = (K) rowpattern[ILUij.index()]; 
00203           }
00204 
00205         //      printmatrix(std::cout,ILU,"ilu pattern","row",10,2);
00206 
00207         // copy entries of A
00208         for (crowiterator i=A.begin(); i!=endi; ++i)
00209           {
00210                 coliterator ILUij;
00211                 coliterator endILUij = ILU[i.index()].end();;
00212                 for (ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
00213                   (*ILUij) = 0; // clear row
00214                 ccoliterator Aij = (*i).begin();
00215                 ccoliterator endAij = (*i).end();
00216                 ILUij = ILU[i.index()].begin();
00217                 while (Aij!=endAij && ILUij!=endILUij)
00218                   {
00219                         if (Aij.index()==ILUij.index())
00220                           {
00221                                 *ILUij = *Aij;
00222                                 ++Aij; ++ILUij;
00223                           }
00224                         else
00225                           {
00226                                 if (Aij.index()<ILUij.index())
00227                                   ++Aij;
00228                                 else
00229                                   ++ILUij;
00230                           }
00231                   }
00232           }
00233 
00234         // call decomposition on pattern
00235         bilu0_decomposition(ILU);
00236   }
00237 
00238 
00241 } // end namespace
00242 
00243 #endif

Generated on 9 Apr 2008 with Doxygen (ver 1.5.2) [logfile].