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