gsetc.hh

Go to the documentation of this file.
00001 #ifndef DUNE_GSETC_HH
00002 #define DUNE_GSETC_HH
00003 
00004 #include<cmath>
00005 #include<complex>
00006 #include<iostream>
00007 #include<iomanip>
00008 #include<string>
00009 
00010 #include "istlexception.hh"
00011 
00017 namespace Dune {
00018    
00029   //============================================================
00030   // parameter types
00031   //============================================================
00032 
00034   template<int l>
00035   struct BL {
00036         enum {recursion_level = l};
00037   };
00038 
00039   enum WithDiagType {
00040         withdiag=1,
00041         nodiag=0
00042   };
00043 
00044   enum WithRelaxType {
00045         withrelax=1,
00046         norelax=0
00047   };
00048 
00049   //============================================================
00050   // generic triangular solves
00051   // consider block decomposition A = L + D + U
00052   // we can invert L, L+D, U, U+D
00053   // we can apply relaxation or not
00054   // we can recurse over a fixed number of levels
00055   //============================================================
00056 
00057   // template meta program for triangular solves
00058   template<int I, WithDiagType diag, WithRelaxType relax>
00059   struct algmeta_btsolve {
00060         template<class M, class X, class Y, class K>
00061         static void bltsolve (const M& A, X& v, const Y& d, const K& w)
00062         {
00063           // iterator types
00064           typedef typename M::ConstRowIterator rowiterator;
00065           typedef typename M::ConstColIterator coliterator;
00066           typedef typename Y::block_type bblock;
00067 
00068           // local solve at each block and immediate update
00069           rowiterator endi=A.end();
00070           for (rowiterator i=A.begin(); i!=endi; ++i)
00071                 {
00072                   bblock rhs(d[i.index()]);
00073                   coliterator j;
00074                   for (j=(*i).begin(); j.index()<i.index(); ++j)
00075                         (*j).mmv(v[j.index()],rhs);
00076                   algmeta_btsolve<I-1,diag,relax>::bltsolve(*j,v[i.index()],rhs,w);
00077                 }
00078         }
00079         template<class M, class X, class Y, class K>
00080         static void butsolve (const M& A, X& v, const Y& d, const K& w)
00081         {
00082           // iterator types
00083           typedef typename M::ConstRowIterator rowiterator;
00084           typedef typename M::ConstColIterator coliterator;
00085           typedef typename Y::block_type bblock;
00086 
00087           // local solve at each block and immediate update
00088           rowiterator rendi=A.rend();
00089           for (rowiterator i=A.rbegin(); i!=rendi; --i)
00090                 {
00091                   bblock rhs(d[i.index()]);
00092                   coliterator j;
00093                   for (j=(*i).rbegin(); j.index()>i.index(); --j)
00094                         (*j).mmv(v[j.index()],rhs);
00095                   algmeta_btsolve<I-1,diag,relax>::butsolve(*j,v[i.index()],rhs,w);
00096                 }
00097         }
00098   };
00099 
00100   // recursion end ...
00101   template<>
00102   struct algmeta_btsolve<0,withdiag,withrelax> {
00103         template<class M, class X, class Y, class K>
00104         static void bltsolve (const M& A, X& v, const Y& d, const K& w)
00105         {
00106           A.solve(v,d);
00107           v *= w;
00108         }
00109         template<class M, class X, class Y, class K>
00110         static void butsolve (const M& A, X& v, const Y& d, const K& w)
00111         {
00112           A.solve(v,d);
00113           v *= w;
00114         }
00115   };
00116   template<>
00117   struct algmeta_btsolve<0,withdiag,norelax> {
00118         template<class M, class X, class Y, class K>
00119         static void bltsolve (const M& A, X& v, const Y& d, const K& w)
00120         {
00121           A.solve(v,d);
00122         }
00123         template<class M, class X, class Y, class K>
00124         static void butsolve (const M& A, X& v, const Y& d, const K& w)
00125         {
00126           A.solve(v,d);
00127         }
00128   };
00129   template<>
00130   struct algmeta_btsolve<0,nodiag,withrelax> {
00131         template<class M, class X, class Y, class K>
00132         static void bltsolve (const M& A, X& v, const Y& d, const K& w)
00133         {
00134           v = d;
00135           v *= w;
00136         }
00137         template<class M, class X, class Y, class K>
00138         static void butsolve (const M& A, X& v, const Y& d, const K& w)
00139         {
00140           v = d;
00141           v *= w;
00142         }
00143   };
00144   template<>
00145   struct algmeta_btsolve<0,nodiag,norelax> {
00146         template<class M, class X, class Y, class K>
00147         static void bltsolve (const M& A, X& v, const Y& d, const K& w)
00148         {
00149           v = d;
00150         }
00151         template<class M, class X, class Y, class K>
00152         static void butsolve (const M& A, X& v, const Y& d, const K& w)
00153         {
00154           v = d;
00155         }
00156   };
00157 
00158 
00159   // user calls
00160 
00161   // default block recursion level = 1
00162 
00164   template<class M, class X, class Y>
00165   void bltsolve (const M& A, X& v, const Y& d)
00166   {
00167         typename X::field_type w=1;
00168         algmeta_btsolve<1,withdiag,norelax>::bltsolve(A,v,d,w);
00169   }
00171   template<class M, class X, class Y, class K>
00172   void bltsolve (const M& A, X& v, const Y& d, const K& w)
00173   {
00174         algmeta_btsolve<1,withdiag,withrelax>::bltsolve(A,v,d,w);
00175   }
00177   template<class M, class X, class Y>
00178   void ubltsolve (const M& A, X& v, const Y& d)
00179   {
00180         typename X::field_type w=1;
00181         algmeta_btsolve<1,nodiag,norelax>::bltsolve(A,v,d,w);
00182   }
00184   template<class M, class X, class Y, class K>
00185   void ubltsolve (const M& A, X& v, const Y& d, const K& w)
00186   {
00187         algmeta_btsolve<1,nodiag,withrelax>::bltsolve(A,v,d,w);
00188   }
00189 
00191   template<class M, class X, class Y>
00192   void butsolve (const M& A, X& v, const Y& d)
00193   {
00194         typename X::field_type w=1;
00195         algmeta_btsolve<1,withdiag,norelax>::butsolve(A,v,d,w);
00196   }
00198   template<class M, class X, class Y, class K>
00199   void butsolve (const M& A, X& v, const Y& d, const K& w)
00200   {
00201         algmeta_btsolve<1,withdiag,withrelax>::butsolve(A,v,d,w);
00202   }
00204   template<class M, class X, class Y>
00205   void ubutsolve (const M& A, X& v, const Y& d)
00206   {
00207         typename X::field_type w=1;
00208         algmeta_btsolve<1,nodiag,norelax>::butsolve(A,v,d,w);
00209   }
00211   template<class M, class X, class Y, class K>
00212   void ubutsolve (const M& A, X& v, const Y& d, const K& w)
00213   {
00214         algmeta_btsolve<1,nodiag,withrelax>::butsolve(A,v,d,w);
00215   }
00216 
00217   // general block recursion level >= 0
00218 
00220   template<class M, class X, class Y, int l>
00221   void bltsolve (const M& A, X& v, const Y& d, BL<l> bl)
00222   {
00223         typename X::field_type w=1;
00224         algmeta_btsolve<l,withdiag,norelax>::bltsolve(A,v,d,w);
00225   }
00227   template<class M, class X, class Y, class K, int l>
00228   void bltsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl)
00229   {
00230         algmeta_btsolve<l,withdiag,withrelax>::bltsolve(A,v,d,w);
00231   }
00233   template<class M, class X, class Y, int l>
00234   void ubltsolve (const M& A, X& v, const Y& d, BL<l> bl)
00235   {
00236         typename X::field_type w=1;
00237         algmeta_btsolve<l,nodiag,norelax>::bltsolve(A,v,d,w);
00238   }
00240   template<class M, class X, class Y, class K, int l>
00241   void ubltsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl)
00242   {
00243         algmeta_btsolve<l,nodiag,withrelax>::bltsolve(A,v,d,w);
00244   }
00245 
00247   template<class M, class X, class Y, int l>
00248   void butsolve (const M& A, X& v, const Y& d, BL<l> bl)
00249   {
00250         typename X::field_type w=1;
00251         algmeta_btsolve<l,withdiag,norelax>::butsolve(A,v,d,w);
00252   }
00254   template<class M, class X, class Y, class K, int l>
00255   void butsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl)
00256   {
00257         algmeta_btsolve<l,withdiag,withrelax>::butsolve(A,v,d,w);
00258   }
00260   template<class M, class X, class Y, int l>
00261   void ubutsolve (const M& A, X& v, const Y& d, BL<l> bl)
00262   {
00263         typename X::field_type w=1;
00264         algmeta_btsolve<l,nodiag,norelax>::butsolve(A,v,d,w);
00265   }
00267   template<class M, class X, class Y, class K, int l>
00268   void ubutsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl)
00269   {
00270         algmeta_btsolve<l,nodiag,withrelax>::butsolve(A,v,d,w);
00271   }
00272 
00273 
00274 
00275   //============================================================
00276   // generic block diagonal solves
00277   // consider block decomposition A = L + D + U
00278   // we can apply relaxation or not
00279   // we can recurse over a fixed number of levels
00280   //============================================================
00281 
00282   // template meta program for diagonal solves
00283   template<int I, WithRelaxType relax>
00284   struct algmeta_bdsolve {
00285         template<class M, class X, class Y, class K>
00286         static void bdsolve (const M& A, X& v, const Y& d, const K& w)
00287         {
00288           // iterator types
00289           typedef typename M::ConstRowIterator rowiterator;
00290           typedef typename M::ConstColIterator coliterator;
00291 
00292           // local solve at each block and immediate update
00293           rowiterator rendi=A.rend();
00294           for (rowiterator i=A.rbegin(); i!=rendi; --i)
00295                 {
00296                   coliterator ii=(*i).find(i.index());
00297                   algmeta_bdsolve<I-1,relax>::bdsolve(*ii,v[i.index()],d[i.index()],w);
00298                 }
00299         }
00300   };
00301 
00302   // recursion end ...
00303   template<>
00304   struct algmeta_bdsolve<0,withrelax> {
00305         template<class M, class X, class Y, class K>
00306         static void bdsolve (const M& A, X& v, const Y& d, const K& w)
00307         {
00308           A.solve(v,d);
00309           v *= w;
00310         }
00311   };
00312   template<>
00313   struct algmeta_bdsolve<0,norelax> {
00314         template<class M, class X, class Y, class K>
00315         static void bdsolve (const M& A, X& v, const Y& d, const K& w)
00316         {
00317           A.solve(v,d);
00318         }
00319   };
00320 
00321   // user calls
00322 
00323   // default block recursion level = 1
00324 
00326   template<class M, class X, class Y>
00327   void bdsolve (const M& A, X& v, const Y& d)
00328   {
00329         typename X::field_type w=1;
00330         algmeta_bdsolve<1,norelax>::bdsolve(A,v,d,w);
00331   }
00333   template<class M, class X, class Y, class K>
00334   void bdsolve (const M& A, X& v, const Y& d, const K& w)
00335   {
00336         algmeta_bdsolve<1,withrelax>::bdsolve(A,v,d,w);
00337   }
00338 
00339   // general block recursion level >= 0
00340 
00342   template<class M, class X, class Y, int l>
00343   void bdsolve (const M& A, X& v, const Y& d, BL<l> bl)
00344   {
00345         typename X::field_type w=1;
00346         algmeta_bdsolve<l,norelax>::bdsolve(A,v,d,w);
00347   }
00349   template<class M, class X, class Y, class K, int l>
00350   void bdsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl)
00351   {
00352         algmeta_bdsolve<l,withrelax>::bdsolve(A,v,d,w);
00353   }
00354 
00355 
00356 
00357 
00358 
00359   //============================================================
00360   // generic steps of iteration methods
00361   // Jacobi, Gauss-Seidel, SOR, SSOR
00362   // work directly on Ax=b, ie solve M(x^{i+1}-x^i) = w (b-Ax^i)
00363   // we can recurse over a fixed number of levels
00364   //============================================================
00365 
00366 
00367   // template meta program for iterative solver steps
00368   template<int I>
00369   struct algmeta_itsteps {
00370         template<class M, class X, class Y, class K>
00371         static void dbgs (const M& A, X& x, const Y& b, const K& w)
00372         {
00373           typedef typename M::ConstRowIterator rowiterator;
00374           typedef typename M::ConstColIterator coliterator;
00375           typedef typename Y::block_type bblock;
00376           typedef typename X::block_type xblock;
00377           bblock rhs;
00378 
00379           X xold(x); // remember old x
00380 
00381           rowiterator endi=A.end();
00382           for (rowiterator i=A.begin(); i!=endi; ++i)
00383                 {
00384                   rhs = b[i.index()];
00385                   coliterator endj=(*i).end();
00386                   coliterator j=(*i).begin();
00387                   for (; j.index()<i.index(); ++j)
00388                         (*j).mmv(x[j.index()],rhs);
00389                   coliterator diag=j; 
00390                   for (; j!=endj; ++j)
00391                         (*j).mmv(x[j.index()],rhs);
00392                   algmeta_itsteps<I-1>::dbgs(*diag,x[i.index()],rhs,w);
00393                 }
00394           x *= w;
00395           x.axpy(1-w,xold);
00396         }
00397         template<class M, class X, class Y, class K>
00398         static void bsorf (const M& A, X& x, const Y& b, const K& w)
00399         {
00400           typedef typename M::ConstRowIterator rowiterator;
00401           typedef typename M::ConstColIterator coliterator;
00402           typedef typename Y::block_type bblock;
00403           typedef typename X::block_type xblock;
00404           bblock rhs;
00405           xblock v;
00406 
00407           // Initialize nested data structure if there are entries
00408           if(A.begin()!=A.end())
00409             v=x[0];
00410 
00411           rowiterator endi=A.end();
00412           for (rowiterator i=A.begin(); i!=endi; ++i)
00413                 {
00414                   rhs = b[i.index()];
00415                   coliterator endj=(*i).end();
00416                   coliterator j=(*i).begin();
00417                   for (; j.index()<i.index(); ++j)
00418                         (*j).mmv(x[j.index()],rhs);
00419                   coliterator diag=j; 
00420                   for (; j!=endj; ++j)
00421                         (*j).mmv(x[j.index()],rhs);
00422                   algmeta_itsteps<I-1>::bsorf(*diag,v,rhs,w);
00423                   x[i.index()].axpy(w,v);
00424                 }
00425         }
00426         template<class M, class X, class Y, class K>
00427         static void bsorb (const M& A, X& x, const Y& b, const K& w)
00428         {
00429           typedef typename M::ConstRowIterator rowiterator;
00430           typedef typename M::ConstColIterator coliterator;
00431           typedef typename Y::block_type bblock;
00432           typedef typename X::block_type xblock;
00433           bblock rhs;
00434           xblock v;
00435 
00436           // Initialize nested data structure if there are entries
00437           if(A.begin()!=A.end())
00438             v=x[0];
00439           
00440           rowiterator endi=A.rend();
00441           for (rowiterator i=A.rbegin(); i!=endi; --i)
00442                 {
00443                   rhs = b[i.index()];
00444                   coliterator endj=(*i).end();
00445                   coliterator j=(*i).begin();
00446                   for (; j.index()<i.index(); ++j)
00447                         (*j).mmv(x[j.index()],rhs);
00448                   coliterator diag=j; 
00449                   for (; j!=endj; ++j)
00450                         (*j).mmv(x[j.index()],rhs);
00451                   algmeta_itsteps<I-1>::bsorb(*diag,v,rhs,w);
00452                   x[i.index()].axpy(w,v);
00453                 }
00454         }
00455         template<class M, class X, class Y, class K>
00456         static void dbjac (const M& A, X& x, const Y& b, const K& w)
00457         {
00458           typedef typename M::ConstRowIterator rowiterator;
00459           typedef typename M::ConstColIterator coliterator;
00460           typedef typename Y::block_type bblock;
00461           typedef typename X::block_type xblock;
00462           bblock rhs;
00463 
00464           X v(x); // allocate with same size
00465 
00466           rowiterator endi=A.end();
00467           for (rowiterator i=A.begin(); i!=endi; ++i)
00468                 {
00469                   rhs = b[i.index()];
00470                   coliterator endj=(*i).end();
00471                   coliterator j=(*i).begin();
00472                   for (; j.index()<i.index(); ++j)
00473                         (*j).mmv(x[j.index()],rhs);
00474                   coliterator diag=j; 
00475                   for (; j!=endj; ++j)
00476                         (*j).mmv(x[j.index()],rhs);
00477                   algmeta_itsteps<I-1>::dbjac(*diag,v[i.index()],rhs,w);
00478                 }
00479           x.axpy(w,v);
00480         }
00481   };
00482   // end of recursion
00483   template<>
00484   struct algmeta_itsteps<0> {
00485         template<class M, class X, class Y, class K>
00486         static void dbgs (const M& A, X& x, const Y& b, const K& w)
00487         {
00488           A.solve(x,b);
00489         }
00490         template<class M, class X, class Y, class K>
00491         static void bsorf (const M& A, X& x, const Y& b, const K& w)
00492         {
00493           A.solve(x,b);
00494         }
00495         template<class M, class X, class Y, class K>
00496         static void bsorb (const M& A, X& x, const Y& b, const K& w)
00497         {
00498           A.solve(x,b);
00499         }
00500         template<class M, class X, class Y, class K>
00501         static void dbjac (const M& A, X& x, const Y& b, const K& w)
00502         {
00503           A.solve(x,b);
00504         }
00505   };
00506 
00507 
00508   // user calls
00509 
00511   template<class M, class X, class Y, class K>
00512   void dbgs (const M& A, X& x, const Y& b, const K& w)
00513   {
00514         algmeta_itsteps<1>::dbgs(A,x,b,w);
00515   }
00517   template<class M, class X, class Y, class K, int l>
00518   void dbgs (const M& A, X& x, const Y& b, const K& w, BL<l> bl)
00519   {
00520         algmeta_itsteps<l>::dbgs(A,x,b,w);
00521   }
00523   template<class M, class X, class Y, class K>
00524   void bsorf (const M& A, X& x, const Y& b, const K& w)
00525   {
00526         algmeta_itsteps<1>::bsorf(A,x,b,w);
00527   }
00529   template<class M, class X, class Y, class K, int l>
00530   void bsorf (const M& A, X& x, const Y& b, const K& w, BL<l> bl)
00531   {
00532         algmeta_itsteps<l>::bsorf(A,x,b,w);
00533   }
00535   template<class M, class X, class Y, class K>
00536   void bsorb (const M& A, X& x, const Y& b, const K& w)
00537   {
00538         algmeta_itsteps<1>::bsorb(A,x,b,w);
00539   }
00541   template<class M, class X, class Y, class K, int l>
00542   void bsorb (const M& A, X& x, const Y& b, const K& w, BL<l> bl)
00543   {
00544         algmeta_itsteps<l>::bsorb(A,x,b,w);
00545   }
00547   template<class M, class X, class Y, class K>
00548   void dbjac (const M& A, X& x, const Y& b, const K& w)
00549   {
00550         algmeta_itsteps<1>::dbjac(A,x,b,w);
00551   }
00553   template<class M, class X, class Y, class K, int l>
00554   void dbjac (const M& A, X& x, const Y& b, const K& w, BL<l> bl)
00555   {
00556         algmeta_itsteps<l>::dbjac(A,x,b,w);
00557   }
00558 
00559 
00562 } // end namespace
00563 
00564 #endif

Generated on Sun Nov 15 22:29:35 2009 for dune-istl by  doxygen 1.5.6