- Home
- About DUNE
- Download
- Documentation
- Community
- Development
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*- 00002 // vi: set et ts=4 sw=2 sts=2: 00003 #ifndef DUNE_GSETC_HH 00004 #define DUNE_GSETC_HH 00005 00006 #include<cmath> 00007 #include<complex> 00008 #include<iostream> 00009 #include<iomanip> 00010 #include<string> 00011 #include "multitypeblockvector.hh" 00012 #include "multitypeblockmatrix.hh" 00013 00014 #include "istlexception.hh" 00015 00016 00022 namespace Dune { 00023 00034 //============================================================ 00035 // parameter types 00036 //============================================================ 00037 00039 template<int l> 00040 struct BL { 00041 enum {recursion_level = l}; 00042 }; 00043 00044 enum WithDiagType { 00045 withdiag=1, 00046 nodiag=0 00047 }; 00048 00049 enum WithRelaxType { 00050 withrelax=1, 00051 norelax=0 00052 }; 00053 00054 //============================================================ 00055 // generic triangular solves 00056 // consider block decomposition A = L + D + U 00057 // we can invert L, L+D, U, U+D 00058 // we can apply relaxation or not 00059 // we can recurse over a fixed number of levels 00060 //============================================================ 00061 00062 // template meta program for triangular solves 00063 template<int I, WithDiagType diag, WithRelaxType relax> 00064 struct algmeta_btsolve { 00065 template<class M, class X, class Y, class K> 00066 static void bltsolve (const M& A, X& v, const Y& d, const K& w) 00067 { 00068 // iterator types 00069 typedef typename M::ConstRowIterator rowiterator; 00070 typedef typename M::ConstColIterator coliterator; 00071 typedef typename Y::block_type bblock; 00072 00073 // local solve at each block and immediate update 00074 rowiterator endi=A.end(); 00075 for (rowiterator i=A.begin(); i!=endi; ++i) 00076 { 00077 bblock rhs(d[i.index()]); 00078 coliterator j; 00079 for (j=(*i).begin(); j.index()<i.index(); ++j) 00080 (*j).mmv(v[j.index()],rhs); 00081 algmeta_btsolve<I-1,diag,relax>::bltsolve(*j,v[i.index()],rhs,w); 00082 } 00083 } 00084 template<class M, class X, class Y, class K> 00085 static void butsolve (const M& A, X& v, const Y& d, const K& w) 00086 { 00087 // iterator types 00088 typedef typename M::ConstRowIterator rowiterator; 00089 typedef typename M::ConstColIterator coliterator; 00090 typedef typename Y::block_type bblock; 00091 00092 // local solve at each block and immediate update 00093 rowiterator rendi=A.beforeBegin(); 00094 for (rowiterator i=A.beforeEnd(); i!=rendi; --i) 00095 { 00096 bblock rhs(d[i.index()]); 00097 coliterator j; 00098 for (j=(*i).beforeEnd(); j.index()>i.index(); --j) 00099 (*j).mmv(v[j.index()],rhs); 00100 algmeta_btsolve<I-1,diag,relax>::butsolve(*j,v[i.index()],rhs,w); 00101 } 00102 } 00103 }; 00104 00105 // recursion end ... 00106 template<> 00107 struct algmeta_btsolve<0,withdiag,withrelax> { 00108 template<class M, class X, class Y, class K> 00109 static void bltsolve (const M& A, X& v, const Y& d, const K& w) 00110 { 00111 A.solve(v,d); 00112 v *= w; 00113 } 00114 template<class M, class X, class Y, class K> 00115 static void butsolve (const M& A, X& v, const Y& d, const K& w) 00116 { 00117 A.solve(v,d); 00118 v *= w; 00119 } 00120 }; 00121 template<> 00122 struct algmeta_btsolve<0,withdiag,norelax> { 00123 template<class M, class X, class Y, class K> 00124 static void bltsolve (const M& A, X& v, const Y& d, const K& w) 00125 { 00126 A.solve(v,d); 00127 } 00128 template<class M, class X, class Y, class K> 00129 static void butsolve (const M& A, X& v, const Y& d, const K& w) 00130 { 00131 A.solve(v,d); 00132 } 00133 }; 00134 template<> 00135 struct algmeta_btsolve<0,nodiag,withrelax> { 00136 template<class M, class X, class Y, class K> 00137 static void bltsolve (const M& A, X& v, const Y& d, const K& w) 00138 { 00139 v = d; 00140 v *= w; 00141 } 00142 template<class M, class X, class Y, class K> 00143 static void butsolve (const M& A, X& v, const Y& d, const K& w) 00144 { 00145 v = d; 00146 v *= w; 00147 } 00148 }; 00149 template<> 00150 struct algmeta_btsolve<0,nodiag,norelax> { 00151 template<class M, class X, class Y, class K> 00152 static void bltsolve (const M& A, X& v, const Y& d, const K& w) 00153 { 00154 v = d; 00155 } 00156 template<class M, class X, class Y, class K> 00157 static void butsolve (const M& A, X& v, const Y& d, const K& w) 00158 { 00159 v = d; 00160 } 00161 }; 00162 00163 00164 // user calls 00165 00166 // default block recursion level = 1 00167 00169 template<class M, class X, class Y> 00170 void bltsolve (const M& A, X& v, const Y& d) 00171 { 00172 typename X::field_type w=1; 00173 algmeta_btsolve<1,withdiag,norelax>::bltsolve(A,v,d,w); 00174 } 00176 template<class M, class X, class Y, class K> 00177 void bltsolve (const M& A, X& v, const Y& d, const K& w) 00178 { 00179 algmeta_btsolve<1,withdiag,withrelax>::bltsolve(A,v,d,w); 00180 } 00182 template<class M, class X, class Y> 00183 void ubltsolve (const M& A, X& v, const Y& d) 00184 { 00185 typename X::field_type w=1; 00186 algmeta_btsolve<1,nodiag,norelax>::bltsolve(A,v,d,w); 00187 } 00189 template<class M, class X, class Y, class K> 00190 void ubltsolve (const M& A, X& v, const Y& d, const K& w) 00191 { 00192 algmeta_btsolve<1,nodiag,withrelax>::bltsolve(A,v,d,w); 00193 } 00194 00196 template<class M, class X, class Y> 00197 void butsolve (const M& A, X& v, const Y& d) 00198 { 00199 typename X::field_type w=1; 00200 algmeta_btsolve<1,withdiag,norelax>::butsolve(A,v,d,w); 00201 } 00203 template<class M, class X, class Y, class K> 00204 void butsolve (const M& A, X& v, const Y& d, const K& w) 00205 { 00206 algmeta_btsolve<1,withdiag,withrelax>::butsolve(A,v,d,w); 00207 } 00209 template<class M, class X, class Y> 00210 void ubutsolve (const M& A, X& v, const Y& d) 00211 { 00212 typename X::field_type w=1; 00213 algmeta_btsolve<1,nodiag,norelax>::butsolve(A,v,d,w); 00214 } 00216 template<class M, class X, class Y, class K> 00217 void ubutsolve (const M& A, X& v, const Y& d, const K& w) 00218 { 00219 algmeta_btsolve<1,nodiag,withrelax>::butsolve(A,v,d,w); 00220 } 00221 00222 // general block recursion level >= 0 00223 00225 template<class M, class X, class Y, int l> 00226 void bltsolve (const M& A, X& v, const Y& d, BL<l> bl) 00227 { 00228 typename X::field_type w=1; 00229 algmeta_btsolve<l,withdiag,norelax>::bltsolve(A,v,d,w); 00230 } 00232 template<class M, class X, class Y, class K, int l> 00233 void bltsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl) 00234 { 00235 algmeta_btsolve<l,withdiag,withrelax>::bltsolve(A,v,d,w); 00236 } 00238 template<class M, class X, class Y, int l> 00239 void ubltsolve (const M& A, X& v, const Y& d, BL<l> bl) 00240 { 00241 typename X::field_type w=1; 00242 algmeta_btsolve<l,nodiag,norelax>::bltsolve(A,v,d,w); 00243 } 00245 template<class M, class X, class Y, class K, int l> 00246 void ubltsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl) 00247 { 00248 algmeta_btsolve<l,nodiag,withrelax>::bltsolve(A,v,d,w); 00249 } 00250 00252 template<class M, class X, class Y, int l> 00253 void butsolve (const M& A, X& v, const Y& d, BL<l> bl) 00254 { 00255 typename X::field_type w=1; 00256 algmeta_btsolve<l,withdiag,norelax>::butsolve(A,v,d,w); 00257 } 00259 template<class M, class X, class Y, class K, int l> 00260 void butsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl) 00261 { 00262 algmeta_btsolve<l,withdiag,withrelax>::butsolve(A,v,d,w); 00263 } 00265 template<class M, class X, class Y, int l> 00266 void ubutsolve (const M& A, X& v, const Y& d, BL<l> bl) 00267 { 00268 typename X::field_type w=1; 00269 algmeta_btsolve<l,nodiag,norelax>::butsolve(A,v,d,w); 00270 } 00272 template<class M, class X, class Y, class K, int l> 00273 void ubutsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl) 00274 { 00275 algmeta_btsolve<l,nodiag,withrelax>::butsolve(A,v,d,w); 00276 } 00277 00278 00279 00280 //============================================================ 00281 // generic block diagonal solves 00282 // consider block decomposition A = L + D + U 00283 // we can apply relaxation or not 00284 // we can recurse over a fixed number of levels 00285 //============================================================ 00286 00287 // template meta program for diagonal solves 00288 template<int I, WithRelaxType relax> 00289 struct algmeta_bdsolve { 00290 template<class M, class X, class Y, class K> 00291 static void bdsolve (const M& A, X& v, const Y& d, const K& w) 00292 { 00293 // iterator types 00294 typedef typename M::ConstRowIterator rowiterator; 00295 typedef typename M::ConstColIterator coliterator; 00296 00297 // local solve at each block and immediate update 00298 rowiterator rendi=A.beforeBegin(); 00299 for (rowiterator i=A.beforeEnd(); i!=rendi; --i) 00300 { 00301 coliterator ii=(*i).find(i.index()); 00302 algmeta_bdsolve<I-1,relax>::bdsolve(*ii,v[i.index()],d[i.index()],w); 00303 } 00304 } 00305 }; 00306 00307 // recursion end ... 00308 template<> 00309 struct algmeta_bdsolve<0,withrelax> { 00310 template<class M, class X, class Y, class K> 00311 static void bdsolve (const M& A, X& v, const Y& d, const K& w) 00312 { 00313 A.solve(v,d); 00314 v *= w; 00315 } 00316 }; 00317 template<> 00318 struct algmeta_bdsolve<0,norelax> { 00319 template<class M, class X, class Y, class K> 00320 static void bdsolve (const M& A, X& v, const Y& d, const K& w) 00321 { 00322 A.solve(v,d); 00323 } 00324 }; 00325 00326 // user calls 00327 00328 // default block recursion level = 1 00329 00331 template<class M, class X, class Y> 00332 void bdsolve (const M& A, X& v, const Y& d) 00333 { 00334 typename X::field_type w=1; 00335 algmeta_bdsolve<1,norelax>::bdsolve(A,v,d,w); 00336 } 00338 template<class M, class X, class Y, class K> 00339 void bdsolve (const M& A, X& v, const Y& d, const K& w) 00340 { 00341 algmeta_bdsolve<1,withrelax>::bdsolve(A,v,d,w); 00342 } 00343 00344 // general block recursion level >= 0 00345 00347 template<class M, class X, class Y, int l> 00348 void bdsolve (const M& A, X& v, const Y& d, BL<l> bl) 00349 { 00350 typename X::field_type w=1; 00351 algmeta_bdsolve<l,norelax>::bdsolve(A,v,d,w); 00352 } 00354 template<class M, class X, class Y, class K, int l> 00355 void bdsolve (const M& A, X& v, const Y& d, const K& w, BL<l> bl) 00356 { 00357 algmeta_bdsolve<l,withrelax>::bdsolve(A,v,d,w); 00358 } 00359 00360 00361 00362 00363 00364 //============================================================ 00365 // generic steps of iteration methods 00366 // Jacobi, Gauss-Seidel, SOR, SSOR 00367 // work directly on Ax=b, ie solve M(x^{i+1}-x^i) = w (b-Ax^i) 00368 // we can recurse over a fixed number of levels 00369 //============================================================ 00370 00371 00372 // template meta program for iterative solver steps 00373 template<int I> 00374 struct algmeta_itsteps { 00375 00376 #ifdef HAVE_BOOST_FUSION 00377 00378 template<typename T11, typename T12, typename T13, typename T14, 00379 typename T15, typename T16, typename T17, typename T18, 00380 typename T19, typename T21, typename T22, typename T23, 00381 typename T24, typename T25, typename T26, typename T27, 00382 typename T28, typename T29, class K> 00383 static void dbgs (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A, 00384 MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x, 00385 const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b, 00386 const K& w) 00387 { 00388 const int rowcount = 00389 boost::mpl::size<MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19> >::value; 00390 Dune::MultiTypeBlockMatrix_Solver<I,0,rowcount>::dbgs(A, x, b, w); 00391 } 00392 #endif 00393 00394 template<class M, class X, class Y, class K> 00395 static void dbgs (const M& A, X& x, const Y& b, const K& w) 00396 { 00397 typedef typename M::ConstRowIterator rowiterator; 00398 typedef typename M::ConstColIterator coliterator; 00399 typedef typename Y::block_type bblock; 00400 typedef typename X::block_type xblock; 00401 bblock rhs; 00402 00403 X xold(x); // remember old x 00404 00405 rowiterator endi=A.end(); 00406 for (rowiterator i=A.begin(); i!=endi; ++i) 00407 { 00408 rhs = b[i.index()]; 00409 coliterator endj=(*i).end(); 00410 coliterator j=(*i).begin(); 00411 for (; j.index()<i.index(); ++j) 00412 (*j).mmv(x[j.index()],rhs); 00413 coliterator diag=j; 00414 for (; j != endj; ++j) 00415 (*j).mmv(x[j.index()],rhs); 00416 algmeta_itsteps<I-1>::dbgs(*diag,x[i.index()],rhs,w); 00417 } 00418 x *= w; 00419 x.axpy(1-w,xold); 00420 } 00421 00422 #ifdef HAVE_BOOST_FUSION 00423 template<typename T11, typename T12, typename T13, typename T14, 00424 typename T15, typename T16, typename T17, typename T18, 00425 typename T19, typename T21, typename T22, typename T23, 00426 typename T24, typename T25, typename T26, typename T27, 00427 typename T28, typename T29, class K> 00428 static void bsorf (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A, 00429 MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x, 00430 const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b, 00431 const K& w) 00432 { 00433 const int rowcount = 00434 boost::mpl::size<MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19> >::value; 00435 Dune::MultiTypeBlockMatrix_Solver<I,0,rowcount>::bsorf(A, x, b, w); 00436 } 00437 #endif 00438 00439 template<class M, class X, class Y, class K> 00440 static void bsorf (const M& A, X& x, const Y& b, const K& w) 00441 { 00442 typedef typename M::ConstRowIterator rowiterator; 00443 typedef typename M::ConstColIterator coliterator; 00444 typedef typename Y::block_type bblock; 00445 typedef typename X::block_type xblock; 00446 bblock rhs; 00447 xblock v; 00448 00449 // Initialize nested data structure if there are entries 00450 if(A.begin()!=A.end()) 00451 v=x[0]; 00452 00453 rowiterator endi=A.end(); 00454 for (rowiterator i=A.begin(); i!=endi; ++i) 00455 { 00456 rhs = b[i.index()]; 00457 coliterator endj=(*i).end(); 00458 coliterator j=(*i).begin(); 00459 for (; j.index()<i.index(); ++j) 00460 (*j).mmv(x[j.index()],rhs); 00461 coliterator diag=j; 00462 for (; j!=endj; ++j) 00463 (*j).mmv(x[j.index()],rhs); 00464 algmeta_itsteps<I-1>::bsorf(*diag,v,rhs,w); 00465 x[i.index()].axpy(w,v); 00466 } 00467 } 00468 00469 #ifdef HAVE_BOOST_FUSION 00470 00471 template<typename T11, typename T12, typename T13, typename T14, 00472 typename T15, typename T16, typename T17, typename T18, 00473 typename T19, typename T21, typename T22, typename T23, 00474 typename T24, typename T25, typename T26, typename T27, 00475 typename T28, typename T29, class K> 00476 static void bsorb (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A, 00477 MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x, 00478 const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b, 00479 const K& w) 00480 { 00481 const int rowcount = 00482 mpl::size<MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19> >::value; 00483 Dune::MultiTypeBlockMatrix_Solver<I,rowcount-1,rowcount>::bsorb(A, x, b, w); 00484 } 00485 #endif 00486 00487 template<class M, class X, class Y, class K> 00488 static void bsorb (const M& A, X& x, const Y& b, const K& w) 00489 { 00490 typedef typename M::ConstRowIterator rowiterator; 00491 typedef typename M::ConstColIterator coliterator; 00492 typedef typename Y::block_type bblock; 00493 typedef typename X::block_type xblock; 00494 bblock rhs; 00495 xblock v; 00496 00497 // Initialize nested data structure if there are entries 00498 if(A.begin()!=A.end()) 00499 v=x[0]; 00500 00501 rowiterator endi=A.beforeBegin(); 00502 for (rowiterator i=A.beforeEnd(); i!=endi; --i) 00503 { 00504 rhs = b[i.index()]; 00505 coliterator endj=(*i).end(); 00506 coliterator j=(*i).begin(); 00507 for (; j.index()<i.index(); ++j) 00508 (*j).mmv(x[j.index()],rhs); 00509 coliterator diag=j; 00510 for (; j!=endj; ++j) 00511 (*j).mmv(x[j.index()],rhs); 00512 algmeta_itsteps<I-1>::bsorb(*diag,v,rhs,w); 00513 x[i.index()].axpy(w,v); 00514 } 00515 } 00516 00517 #ifdef HAVE_BOOST_FUSION 00518 template<typename T11, typename T12, typename T13, typename T14, 00519 typename T15, typename T16, typename T17, typename T18, 00520 typename T19, typename T21, typename T22, typename T23, 00521 typename T24, typename T25, typename T26, typename T27, 00522 typename T28, typename T29, class K> 00523 static void dbjac (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A, 00524 MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x, 00525 const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b, 00526 const K& w) 00527 { 00528 const int rowcount = 00529 boost::mpl::size<MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19> >::value; 00530 Dune::MultiTypeBlockMatrix_Solver<I,0,rowcount >::dbjac(A, x, b, w); 00531 } 00532 #endif 00533 00534 template<class M, class X, class Y, class K> 00535 static void dbjac (const M& A, X& x, const Y& b, const K& w) 00536 { 00537 typedef typename M::ConstRowIterator rowiterator; 00538 typedef typename M::ConstColIterator coliterator; 00539 typedef typename Y::block_type bblock; 00540 typedef typename X::block_type xblock; 00541 bblock rhs; 00542 00543 X v(x); // allocate with same size 00544 00545 rowiterator endi=A.end(); 00546 for (rowiterator i=A.begin(); i!=endi; ++i) 00547 { 00548 rhs = b[i.index()]; 00549 coliterator endj=(*i).end(); 00550 coliterator j=(*i).begin(); 00551 for (; j.index()<i.index(); ++j) 00552 (*j).mmv(x[j.index()],rhs); 00553 coliterator diag=j; 00554 for (; j!=endj; ++j) 00555 (*j).mmv(x[j.index()],rhs); 00556 algmeta_itsteps<I-1>::dbjac(*diag,v[i.index()],rhs,w); 00557 } 00558 x.axpy(w,v); 00559 } 00560 }; 00561 // end of recursion 00562 template<> 00563 struct algmeta_itsteps<0> { 00564 template<class M, class X, class Y, class K> 00565 static void dbgs (const M& A, X& x, const Y& b, const K& w) 00566 { 00567 A.solve(x,b); 00568 } 00569 template<class M, class X, class Y, class K> 00570 static void bsorf (const M& A, X& x, const Y& b, const K& w) 00571 { 00572 A.solve(x,b); 00573 } 00574 template<class M, class X, class Y, class K> 00575 static void bsorb (const M& A, X& x, const Y& b, const K& w) 00576 { 00577 A.solve(x,b); 00578 } 00579 template<class M, class X, class Y, class K> 00580 static void dbjac (const M& A, X& x, const Y& b, const K& w) 00581 { 00582 A.solve(x,b); 00583 } 00584 }; 00585 00586 00587 // user calls 00588 00590 template<class M, class X, class Y, class K> 00591 void dbgs (const M& A, X& x, const Y& b, const K& w) 00592 { 00593 algmeta_itsteps<1>::dbgs(A,x,b,w); 00594 } 00596 template<class M, class X, class Y, class K, int l> 00597 void dbgs (const M& A, X& x, const Y& b, const K& w, BL<l> bl) 00598 { 00599 algmeta_itsteps<l>::dbgs(A,x,b,w); 00600 } 00602 template<class M, class X, class Y, class K> 00603 void bsorf (const M& A, X& x, const Y& b, const K& w) 00604 { 00605 algmeta_itsteps<1>::bsorf(A,x,b,w); 00606 } 00608 template<class M, class X, class Y, class K, int l> 00609 void bsorf (const M& A, X& x, const Y& b, const K& w, BL<l> bl) 00610 { 00611 algmeta_itsteps<l>::bsorf(A,x,b,w); 00612 } 00614 template<class M, class X, class Y, class K> 00615 void bsorb (const M& A, X& x, const Y& b, const K& w) 00616 { 00617 algmeta_itsteps<1>::bsorb(A,x,b,w); 00618 } 00620 template<class M, class X, class Y, class K, int l> 00621 void bsorb (const M& A, X& x, const Y& b, const K& w, BL<l> bl) 00622 { 00623 algmeta_itsteps<l>::bsorb(A,x,b,w); 00624 } 00626 template<class M, class X, class Y, class K> 00627 void dbjac (const M& A, X& x, const Y& b, const K& w) 00628 { 00629 algmeta_itsteps<1>::dbjac(A,x,b,w); 00630 } 00632 template<class M, class X, class Y, class K, int l> 00633 void dbjac (const M& A, X& x, const Y& b, const K& w, BL<l> bl) 00634 { 00635 algmeta_itsteps<l>::dbjac(A,x,b,w); 00636 } 00637 00638 00641 } // end namespace 00642 00643 #endif
Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].