Heidelberg Educational Numerics Library Version 0.24 (from 9 September 2011)
|
00001 // -*- tab-width: 4; indent-tabs-mode: nil -*- 00002 #ifndef HDNUM_ODE_HH 00003 #define HDNUM_ODE_HH 00004 00005 #include<vector> 00006 #include "newton.hh" 00007 00012 namespace hdnum { 00013 00022 template<class M> 00023 class EE 00024 { 00025 public: 00027 typedef typename M::size_type size_type; 00028 00030 typedef typename M::time_type time_type; 00031 00033 typedef typename M::number_type number_type; 00034 00036 EE (const M& model_) 00037 : model(model_), u(model.size()), f(model.size()) 00038 { 00039 model.initialize(t,u); 00040 dt = 0.1; 00041 } 00042 00044 void set_dt (time_type dt_) 00045 { 00046 dt = dt_; 00047 } 00048 00050 void step () 00051 { 00052 model.f(t,u,f); // evaluate model 00053 u.update(dt,f); // advance state 00054 t += dt; // advance time 00055 } 00056 00058 void set_state (time_type t_, const Vector<number_type>& u_) 00059 { 00060 t = t_; 00061 u = u_; 00062 } 00063 00065 const Vector<number_type>& get_state () const 00066 { 00067 return u; 00068 } 00069 00071 time_type get_time () const 00072 { 00073 return t; 00074 } 00075 00077 time_type get_dt () const 00078 { 00079 return dt; 00080 } 00081 00083 size_type get_order () const 00084 { 00085 return 1; 00086 } 00087 00088 private: 00089 const M& model; 00090 time_type t, dt; 00091 Vector<number_type> u; 00092 Vector<number_type> f; 00093 }; 00094 00103 template<class M> 00104 class ModifiedEuler 00105 { 00106 public: 00108 typedef typename M::size_type size_type; 00109 00111 typedef typename M::time_type time_type; 00112 00114 typedef typename M::number_type number_type; 00115 00117 ModifiedEuler (const M& model_) 00118 : model(model_), u(model.size()), w(model.size()), k1(model.size()), k2(model.size()) 00119 { 00120 c2 = 0.5; 00121 a21 = 0.5; 00122 b2 = 1.0; 00123 model.initialize(t,u); 00124 dt = 0.1; 00125 } 00126 00128 void set_dt (time_type dt_) 00129 { 00130 dt = dt_; 00131 } 00132 00134 void step () 00135 { 00136 // stage 1 00137 model.f(t,u,k1); 00138 00139 // stage 2 00140 w = u; 00141 w.update(dt*a21,k1); 00142 model.f(t+c2*dt,w,k2); 00143 00144 // final 00145 u.update(dt*b2,k2); 00146 t += dt; 00147 } 00148 00150 void set_state (time_type t_, const Vector<number_type>& u_) 00151 { 00152 t = t_; 00153 u = u_; 00154 } 00155 00157 const Vector<number_type>& get_state () const 00158 { 00159 return u; 00160 } 00161 00163 time_type get_time () const 00164 { 00165 return t; 00166 } 00167 00169 time_type get_dt () const 00170 { 00171 return dt; 00172 } 00173 00175 size_type get_order () const 00176 { 00177 return 2; 00178 } 00179 00180 private: 00181 const M& model; 00182 time_type t, dt; 00183 time_type c2,a21,b2; 00184 Vector<number_type> u,w; 00185 Vector<number_type> k1,k2; 00186 }; 00187 00188 00197 template<class M> 00198 class Heun2 00199 { 00200 public: 00202 typedef typename M::size_type size_type; 00203 00205 typedef typename M::time_type time_type; 00206 00208 typedef typename M::number_type number_type; 00209 00211 Heun2 (const M& model_) 00212 : model(model_), u(model.size()), w(model.size()), k1(model.size()), k2(model.size()) 00213 { 00214 c2 = 1.0; 00215 a21 = 1.0; 00216 b1 = 0.5; 00217 b2 = 0.5; 00218 model.initialize(t,u); 00219 dt = 0.1; 00220 } 00221 00223 void set_dt (time_type dt_) 00224 { 00225 dt = dt_; 00226 } 00227 00229 void step () 00230 { 00231 // stage 1 00232 model.f(t,u,k1); 00233 00234 // stage 2 00235 w = u; 00236 w.update(dt*a21,k1); 00237 model.f(t+c2*dt,w,k2); 00238 00239 // final 00240 u.update(dt*b1,k1); 00241 u.update(dt*b2,k2); 00242 t += dt; 00243 } 00244 00246 void set_state (time_type t_, const Vector<number_type>& u_) 00247 { 00248 t = t_; 00249 u = u_; 00250 } 00251 00253 const Vector<number_type>& get_state () const 00254 { 00255 return u; 00256 } 00257 00259 time_type get_time () const 00260 { 00261 return t; 00262 } 00263 00265 time_type get_dt () const 00266 { 00267 return dt; 00268 } 00269 00271 size_type get_order () const 00272 { 00273 return 2; 00274 } 00275 00276 private: 00277 const M& model; 00278 time_type t, dt; 00279 time_type c2,a21,b1,b2; 00280 Vector<number_type> u,w; 00281 Vector<number_type> k1,k2; 00282 }; 00283 00284 00293 template<class M> 00294 class Heun3 00295 { 00296 public: 00298 typedef typename M::size_type size_type; 00299 00301 typedef typename M::time_type time_type; 00302 00304 typedef typename M::number_type number_type; 00305 00307 Heun3 (const M& model_) 00308 : model(model_), u(model.size()), w(model.size()), k1(model.size()), 00309 k2(model.size()), k3(model.size()) 00310 { 00311 c2 = time_type(1.0)/time_type(3.0); 00312 c3 = time_type(2.0)/time_type(3.0); 00313 a21 = time_type(1.0)/time_type(3.0); 00314 a32 = time_type(2.0)/time_type(3.0); 00315 b1 = 0.25; 00316 b2 = 0.0; 00317 b3 = 0.75; 00318 model.initialize(t,u); 00319 dt = 0.1; 00320 } 00321 00323 void set_dt (time_type dt_) 00324 { 00325 dt = dt_; 00326 } 00327 00329 void step () 00330 { 00331 // stage 1 00332 model.f(t,u,k1); 00333 00334 // stage 2 00335 w = u; 00336 w.update(dt*a21,k1); 00337 model.f(t+c2*dt,w,k2); 00338 00339 // stage 3 00340 w = u; 00341 w.update(dt*a32,k2); 00342 model.f(t+c3*dt,w,k3); 00343 00344 // final 00345 u.update(dt*b1,k1); 00346 u.update(dt*b3,k3); 00347 t += dt; 00348 } 00349 00351 void set_state (time_type t_, const Vector<number_type>& u_) 00352 { 00353 t = t_; 00354 u = u_; 00355 } 00356 00358 const Vector<number_type>& get_state () const 00359 { 00360 return u; 00361 } 00362 00364 time_type get_time () const 00365 { 00366 return t; 00367 } 00368 00370 time_type get_dt () const 00371 { 00372 return dt; 00373 } 00374 00376 size_type get_order () const 00377 { 00378 return 3; 00379 } 00380 00381 private: 00382 const M& model; 00383 time_type t, dt; 00384 time_type c2,c3,a21,a31,a32,b1,b2,b3; 00385 Vector<number_type> u,w; 00386 Vector<number_type> k1,k2,k3; 00387 }; 00388 00397 template<class M> 00398 class Kutta3 00399 { 00400 public: 00402 typedef typename M::size_type size_type; 00403 00405 typedef typename M::time_type time_type; 00406 00408 typedef typename M::number_type number_type; 00409 00411 Kutta3 (const M& model_) 00412 : model(model_), u(model.size()), w(model.size()), k1(model.size()), 00413 k2(model.size()), k3(model.size()) 00414 { 00415 c2 = 0.5; 00416 c3 = 1.0; 00417 a21 = 0.5; 00418 a31 = -1.0; 00419 a32 = 2.0; 00420 b1 = time_type(1.0)/time_type(6.0); 00421 b2 = time_type(4.0)/time_type(6.0); 00422 b3 = time_type(1.0)/time_type(6.0); 00423 model.initialize(t,u); 00424 dt = 0.1; 00425 } 00426 00428 void set_dt (time_type dt_) 00429 { 00430 dt = dt_; 00431 } 00432 00434 void step () 00435 { 00436 // stage 1 00437 model.f(t,u,k1); 00438 00439 // stage 2 00440 w = u; 00441 w.update(dt*a21,k1); 00442 model.f(t+c2*dt,w,k2); 00443 00444 // stage 3 00445 w = u; 00446 w.update(dt*a31,k1); 00447 w.update(dt*a32,k2); 00448 model.f(t+c3*dt,w,k3); 00449 00450 // final 00451 u.update(dt*b1,k1); 00452 u.update(dt*b2,k2); 00453 u.update(dt*b3,k3); 00454 t += dt; 00455 } 00456 00458 void set_state (time_type t_, const Vector<number_type>& u_) 00459 { 00460 t = t_; 00461 u = u_; 00462 } 00463 00465 const Vector<number_type>& get_state () const 00466 { 00467 return u; 00468 } 00469 00471 time_type get_time () const 00472 { 00473 return t; 00474 } 00475 00477 time_type get_dt () const 00478 { 00479 return dt; 00480 } 00481 00483 size_type get_order () const 00484 { 00485 return 3; 00486 } 00487 00488 private: 00489 const M& model; 00490 time_type t, dt; 00491 time_type c2,c3,a21,a31,a32,b1,b2,b3; 00492 Vector<number_type> u,w; 00493 Vector<number_type> k1,k2,k3; 00494 }; 00495 00504 template<class M> 00505 class RungeKutta4 00506 { 00507 public: 00509 typedef typename M::size_type size_type; 00510 00512 typedef typename M::time_type time_type; 00513 00515 typedef typename M::number_type number_type; 00516 00518 RungeKutta4 (const M& model_) 00519 : model(model_), u(model.size()), w(model.size()), k1(model.size()), 00520 k2(model.size()), k3(model.size()), k4(model.size()) 00521 { 00522 c2 = 0.5; 00523 c3 = 0.5; 00524 c4 = 1.0; 00525 a21 = 0.5; 00526 a32 = 0.5; 00527 a43 = 1.0; 00528 b1 = time_type(1.0)/time_type(6.0); 00529 b2 = time_type(2.0)/time_type(6.0); 00530 b3 = time_type(2.0)/time_type(6.0); 00531 b4 = time_type(1.0)/time_type(6.0); 00532 model.initialize(t,u); 00533 dt = 0.1; 00534 } 00535 00537 void set_dt (time_type dt_) 00538 { 00539 dt = dt_; 00540 } 00541 00543 void step () 00544 { 00545 // stage 1 00546 model.f(t,u,k1); 00547 00548 // stage 2 00549 w = u; 00550 w.update(dt*a21,k1); 00551 model.f(t+c2*dt,w,k2); 00552 00553 // stage 3 00554 w = u; 00555 w.update(dt*a32,k2); 00556 model.f(t+c3*dt,w,k3); 00557 00558 // stage 4 00559 w = u; 00560 w.update(dt*a43,k3); 00561 model.f(t+c4*dt,w,k4); 00562 00563 // final 00564 u.update(dt*b1,k1); 00565 u.update(dt*b2,k2); 00566 u.update(dt*b3,k3); 00567 u.update(dt*b4,k4); 00568 t += dt; 00569 } 00570 00572 void set_state (time_type t_, const Vector<number_type>& u_) 00573 { 00574 t = t_; 00575 u = u_; 00576 } 00577 00579 const Vector<number_type>& get_state () const 00580 { 00581 return u; 00582 } 00583 00585 time_type get_time () const 00586 { 00587 return t; 00588 } 00589 00591 time_type get_dt () const 00592 { 00593 return dt; 00594 } 00595 00597 size_type get_order () const 00598 { 00599 return 4; 00600 } 00601 00602 private: 00603 const M& model; 00604 time_type t, dt; 00605 time_type c2,c3,c4,a21,a32,a43,b1,b2,b3,b4; 00606 Vector<number_type> u,w; 00607 Vector<number_type> k1,k2,k3,k4; 00608 }; 00609 00614 template<class M> 00615 class RKF45 00616 { 00617 public: 00619 typedef typename M::size_type size_type; 00620 00622 typedef typename M::time_type time_type; 00623 00625 typedef typename M::number_type number_type; 00626 00628 RKF45 (const M& model_) 00629 : model(model_), u(model.size()), w(model.size()), ww(model.size()), k1(model.size()), 00630 k2(model.size()), k3(model.size()), k4(model.size()), k5(model.size()), k6(model.size()), 00631 steps(0), rejected(0) 00632 { 00633 TOL = time_type(0.0001); 00634 rho = time_type(0.8); 00635 alpha = time_type(0.25); 00636 beta = time_type(4.0); 00637 dt_min = 1E-12; 00638 00639 c2 = time_type(1.0)/time_type(4.0); 00640 c3 = time_type(3.0)/time_type(8.0); 00641 c4 = time_type(12.0)/time_type(13.0); 00642 c5 = time_type(1.0); 00643 c6 = time_type(1.0)/time_type(2.0); 00644 00645 a21 = time_type(1.0)/time_type(4.0); 00646 00647 a31 = time_type(3.0)/time_type(32.0); 00648 a32 = time_type(9.0)/time_type(32.0); 00649 00650 a41 = time_type(1932.0)/time_type(2197.0); 00651 a42 = time_type(-7200.0)/time_type(2197.0); 00652 a43 = time_type(7296.0)/time_type(2197.0); 00653 00654 a51 = time_type(439.0)/time_type(216.0); 00655 a52 = time_type(-8.0); 00656 a53 = time_type(3680.0)/time_type(513.0); 00657 a54 = time_type(-845.0)/time_type(4104.0); 00658 00659 a61 = time_type(-8.0)/time_type(27.0); 00660 a62 = time_type(2.0); 00661 a63 = time_type(-3544.0)/time_type(2565.0); 00662 a64 = time_type(1859.0)/time_type(4104.0); 00663 a65 = time_type(-11.0)/time_type(40.0); 00664 00665 b1 = time_type(25.0)/time_type(216.0); 00666 b2 = time_type(0.0); 00667 b3 = time_type(1408.0)/time_type(2565.0); 00668 b4 = time_type(2197.0)/time_type(4104.0); 00669 b5 = time_type(-1.0)/time_type(5.0); 00670 00671 bb1 = time_type(16.0)/time_type(135.0); 00672 bb2 = time_type(0.0); 00673 bb3 = time_type(6656.0)/time_type(12825.0); 00674 bb4 = time_type(28561.0)/time_type(56430.0); 00675 bb5 = time_type(-9.0)/time_type(50.0); 00676 bb6 = time_type(2.0)/time_type(55.0); 00677 00678 model.initialize(t,u); 00679 dt = 0.1; 00680 } 00681 00683 void set_dt (time_type dt_) 00684 { 00685 dt = dt_; 00686 } 00687 00689 void set_TOL (time_type TOL_) 00690 { 00691 TOL = TOL_; 00692 } 00693 00695 void step () 00696 { 00697 steps++; 00698 00699 // stage 1 00700 model.f(t,u,k1); 00701 00702 // stage 2 00703 w = u; 00704 w.update(dt*a21,k1); 00705 model.f(t+c2*dt,w,k2); 00706 00707 // stage 3 00708 w = u; 00709 w.update(dt*a31,k1); 00710 w.update(dt*a32,k2); 00711 model.f(t+c3*dt,w,k3); 00712 00713 // stage 4 00714 w = u; 00715 w.update(dt*a41,k1); 00716 w.update(dt*a42,k2); 00717 w.update(dt*a43,k3); 00718 model.f(t+c4*dt,w,k4); 00719 00720 // stage 5 00721 w = u; 00722 w.update(dt*a51,k1); 00723 w.update(dt*a52,k2); 00724 w.update(dt*a53,k3); 00725 w.update(dt*a54,k4); 00726 model.f(t+c5*dt,w,k5); 00727 00728 // stage 6 00729 w = u; 00730 w.update(dt*a61,k1); 00731 w.update(dt*a62,k2); 00732 w.update(dt*a63,k3); 00733 w.update(dt*a64,k4); 00734 w.update(dt*a65,k5); 00735 model.f(t+c6*dt,w,k6); 00736 00737 // compute order 4 approximation 00738 w = u; 00739 w.update(dt*b1,k1); 00740 w.update(dt*b2,k2); 00741 w.update(dt*b3,k3); 00742 w.update(dt*b4,k4); 00743 w.update(dt*b5,k5); 00744 00745 // compute order 5 approximation 00746 ww = u; 00747 ww.update(dt*bb1,k1); 00748 ww.update(dt*bb2,k2); 00749 ww.update(dt*bb3,k3); 00750 ww.update(dt*bb4,k4); 00751 ww.update(dt*bb5,k5); 00752 ww.update(dt*bb6,k6); 00753 00754 // estimate local error 00755 w -= ww; 00756 time_type error(norm(w)); 00757 time_type dt_opt(dt*pow(rho*TOL/error,0.2)); 00758 dt_opt = std::min(beta*dt,std::max(alpha*dt,dt_opt)); 00759 //std::cout << "est. error=" << error << " dt_opt=" << dt_opt << std::endl; 00760 00761 if (error<=TOL) 00762 { 00763 t += dt; 00764 u = ww; 00765 dt = dt_opt; 00766 } 00767 else 00768 { 00769 rejected++; 00770 dt = dt_opt; 00771 if (dt>dt_min) step(); 00772 } 00773 } 00774 00776 const Vector<number_type>& get_state () const 00777 { 00778 return u; 00779 } 00780 00782 time_type get_time () const 00783 { 00784 return t; 00785 } 00786 00788 time_type get_dt () const 00789 { 00790 return dt; 00791 } 00792 00794 size_type get_order () const 00795 { 00796 return 5; 00797 } 00798 00800 void get_info () const 00801 { 00802 std::cout << "RE: steps=" << steps << " rejected=" << rejected << std::endl; 00803 } 00804 00805 private: 00806 const M& model; 00807 time_type t, dt; 00808 time_type TOL,rho,alpha,beta,dt_min; 00809 time_type c2,c3,c4,c5,c6; 00810 time_type a21,a31,a32,a41,a42,a43,a51,a52,a53,a54,a61,a62,a63,a64,a65; 00811 time_type b1,b2,b3,b4,b5; // 4th order 00812 time_type bb1,bb2,bb3,bb4,bb5,bb6; // 5th order 00813 Vector<number_type> u,w,ww; 00814 Vector<number_type> k1,k2,k3,k4,k5,k6; 00815 mutable size_type steps, rejected; 00816 }; 00817 00818 00824 template<class M, class S> 00825 class RE 00826 { 00827 public: 00829 typedef typename M::size_type size_type; 00830 00832 typedef typename M::time_type time_type; 00833 00835 typedef typename M::number_type number_type; 00836 00838 RE (const M& model_, S& solver_) 00839 : model(model_), solver(solver_), u(model.size()), 00840 wlow(model.size()), whigh(model.size()), ww(model.size()), 00841 steps(0), rejected(0) 00842 { 00843 model.initialize(t,u); // initialize state 00844 dt = 0.1; // set initial time step 00845 two_power_m = 1.0; 00846 for (size_type i=0; i<solver.get_order(); i++) 00847 two_power_m *= 2.0; 00848 TOL = time_type(0.0001); 00849 rho = time_type(0.8); 00850 alpha = time_type(0.25); 00851 beta = time_type(4.0); 00852 dt_min = 1E-12; 00853 } 00854 00856 void set_dt (time_type dt_) 00857 { 00858 dt = dt_; 00859 } 00860 00862 void set_TOL (time_type TOL_) 00863 { 00864 TOL = TOL_; 00865 } 00866 00868 void step () 00869 { 00870 // count steps done 00871 steps++; 00872 00873 // do 1 step with 2*dt 00874 time_type H(2.0*dt); 00875 solver.set_state(t,u); 00876 solver.set_dt(H); 00877 solver.step(); 00878 wlow = solver.get_state(); 00879 00880 // do 2 steps with dt 00881 solver.set_state(t,u); 00882 solver.set_dt(dt); 00883 solver.step(); 00884 solver.step(); 00885 whigh = solver.get_state(); 00886 00887 // estimate local error 00888 ww = wlow; 00889 ww -= whigh; 00890 time_type error(norm(ww)/(pow(H,1.0+solver.get_order())*(1.0-1.0/two_power_m))); 00891 time_type dt_opt(pow(rho*TOL/error,1.0/((time_type)solver.get_order()))); 00892 dt_opt = std::min(beta*dt,std::max(alpha*dt,dt_opt)); 00893 //std::cout << "est. error=" << error << " dt_opt=" << dt_opt << std::endl; 00894 00895 if (dt<=dt_opt) 00896 { 00897 t += H; 00898 u = whigh; 00899 u *= two_power_m; 00900 u -= wlow; 00901 u /= two_power_m-1.0; 00902 dt = dt_opt; 00903 } 00904 else 00905 { 00906 rejected++; 00907 dt = dt_opt; 00908 if (dt>dt_min) step(); 00909 } 00910 } 00911 00913 const Vector<number_type>& get_state () const 00914 { 00915 return u; 00916 } 00917 00919 time_type get_time () const 00920 { 00921 return t; 00922 } 00923 00925 time_type get_dt () const 00926 { 00927 return dt; 00928 } 00929 00931 size_type get_order () const 00932 { 00933 return solver.get_order()+1; 00934 } 00935 00937 void get_info () const 00938 { 00939 std::cout << "RE: steps=" << steps << " rejected=" << rejected << std::endl; 00940 } 00941 00942 private: 00943 const M& model; 00944 S& solver; 00945 time_type t, dt; 00946 time_type two_power_m; 00947 Vector<number_type> u,wlow,whigh,ww; 00948 time_type TOL,rho,alpha,beta,dt_min; 00949 mutable size_type steps, rejected; 00950 }; 00951 00952 00962 template<class M, class S> 00963 class IE 00964 { 00966 // h_n f(t_n, y_n) - y_n + y_{n-1} = 0 00967 class NonlinearProblem 00968 { 00969 public: 00971 typedef typename M::size_type size_type; 00972 00974 typedef typename M::number_type number_type; 00975 00977 NonlinearProblem (const M& model_, const Vector<number_type>& yold_, 00978 typename M::time_type tnew_, typename M::time_type dt_) 00979 : model(model_), yold(yold_), tnew(tnew_), dt(dt_) 00980 {} 00981 00983 std::size_t size () const 00984 { 00985 return model.size(); 00986 } 00987 00989 void F (const Vector<number_type>& x, Vector<number_type>& result) const 00990 { 00991 model.f(tnew,x,result); 00992 result *= dt; 00993 result -= x; 00994 result += yold; 00995 } 00996 00998 void F_x (const Vector<number_type>& x, DenseMatrix<number_type>& result) const 00999 { 01000 model.f_x(tnew,x,result); 01001 result *= dt; 01002 for (size_type i=0; i<model.size(); i++) result[i][i] -= number_type(1.0); 01003 } 01004 01005 void set_tnew_dt (typename M::time_type tnew_, typename M::time_type dt_) 01006 { 01007 tnew = tnew_; 01008 dt = dt_; 01009 } 01010 01011 private: 01012 const M& model; 01013 const Vector<number_type>& yold; 01014 typename M::time_type tnew; 01015 typename M::time_type dt; 01016 }; 01017 01018 public: 01020 typedef typename M::size_type size_type; 01021 01023 typedef typename M::time_type time_type; 01024 01026 typedef typename M::number_type number_type; 01027 01029 IE (const M& model_, const S& newton_) 01030 : verbosity(0), model(model_), newton(newton_), u(model.size()), unew(model.size()) 01031 { 01032 model.initialize(t,u); 01033 dt = dtmax = 0.1; 01034 } 01035 01037 void set_dt (time_type dt_) 01038 { 01039 dt = dtmax = dt_; 01040 } 01041 01043 void set_verbosity (size_type verbosity_) 01044 { 01045 verbosity = verbosity_; 01046 } 01047 01049 void step () 01050 { 01051 if (verbosity>=2) 01052 std::cout << "IE: step" << " t=" << t << " dt=" << dt << std::endl; 01053 NonlinearProblem nlp(model,u,t+dt,dt); 01054 bool reduced = false; 01055 error = false; 01056 while (1) 01057 { 01058 unew = u; 01059 newton.solve(nlp,unew); 01060 if (newton.has_converged()) 01061 { 01062 u = unew; 01063 t += dt; 01064 if (!reduced && dt<dtmax-1e-13) 01065 { 01066 dt = std::min(2.0*dt,dtmax); 01067 if (verbosity>0) 01068 std::cout << "IE: increasing time step to " << dt << std::endl; 01069 } 01070 return; 01071 } 01072 else 01073 { 01074 if (dt<1e-12) 01075 { 01076 HDNUM_ERROR("time step too small in implicit Euler"); 01077 error = true; 01078 break; 01079 } 01080 dt *= 0.5; 01081 reduced = true; 01082 nlp.set_tnew_dt(t+dt,dt); 01083 if (verbosity>0) std::cout << "IE: reducing time step to " << dt << std::endl; 01084 } 01085 } 01086 } 01087 01089 bool get_error () const 01090 { 01091 return error; 01092 } 01093 01095 void set_state (time_type t_, const Vector<number_type>& u_) 01096 { 01097 t = t_; 01098 u = u_; 01099 } 01100 01102 const Vector<number_type>& get_state () const 01103 { 01104 return u; 01105 } 01106 01108 time_type get_time () const 01109 { 01110 return t; 01111 } 01112 01114 time_type get_dt () const 01115 { 01116 return dt; 01117 } 01118 01120 size_type get_order () const 01121 { 01122 return 1; 01123 } 01124 01126 void get_info () const 01127 { 01128 } 01129 01130 private: 01131 size_type verbosity; 01132 const M& model; 01133 const S& newton; 01134 time_type t, dt, dtmax; 01135 number_type reduction; 01136 size_type linesearchsteps; 01137 Vector<number_type> u; 01138 Vector<number_type> unew; 01139 mutable bool error; 01140 }; 01141 01152 template<class M, class S> 01153 class DIRK 01154 { 01155 public: 01156 01158 typedef typename M::size_type size_type; 01159 01161 typedef typename M::time_type time_type; 01162 01164 typedef typename M::number_type number_type; 01165 01167 typedef DenseMatrix<number_type> ButcherTableau; 01168 01169 private: 01170 01173 static ButcherTableau initTableau(const std::string method) 01174 { 01175 if(method.find("Implicit Euler") != std::string::npos){ 01176 ButcherTableau butcher(2,2,0.0); 01177 butcher[1][1] = 1; 01178 butcher[0][1] = 1; 01179 butcher[0][0] = 1; 01180 01181 return butcher; 01182 } 01183 else if(method.find("Alexander") != std::string::npos){ 01184 ButcherTableau butcher(3,3,0.0); 01185 const number_type alpha = 1. - sqrt(2.)/2.; 01186 butcher[0][0] = alpha; 01187 butcher[0][1] = alpha; 01188 01189 butcher[1][0] = 1.; 01190 butcher[1][1] = 1. - alpha; 01191 butcher[1][2] = alpha; 01192 01193 butcher[2][1] = 1. - alpha; 01194 butcher[2][2] = alpha; 01195 01196 return butcher; 01197 } 01198 else if(method.find("Crouzieux") != std::string::npos){ 01199 ButcherTableau butcher(3,3,0.0); 01200 const number_type beta = 1./2./sqrt(3); 01201 butcher[0][0] = 0.5 + beta; 01202 butcher[0][1] = 0.5 + beta; 01203 01204 butcher[1][0] = 0.5 - beta; 01205 butcher[1][1] = -1. / sqrt(3); 01206 butcher[1][2] = 0.5 + beta; 01207 01208 butcher[2][1] = 0.5; 01209 butcher[2][2] = 0.5; 01210 01211 return butcher; 01212 } 01213 else if(method.find("Midpoint Rule") != std::string::npos){ 01214 ButcherTableau butcher(2,2,0.0); 01215 butcher[0][0] = 0.5; 01216 butcher[0][1] = 0.5; 01217 butcher[1][1] = 1; 01218 01219 return butcher; 01220 } 01221 else if(method.find("Fractional Step Theta") != std::string::npos){ 01222 ButcherTableau butcher(5,5,0.0); 01223 const number_type theta = 1 - sqrt(2.)/2.; 01224 const number_type alpha = 2. - sqrt(2.); 01225 const number_type beta = 1. - alpha; 01226 butcher[1][0] = theta; 01227 butcher[1][1] = beta * theta; 01228 butcher[1][2] = alpha * theta; 01229 01230 butcher[2][0] = 1.-theta; 01231 butcher[2][1] = beta * theta; 01232 butcher[2][2] = alpha * (1.-theta); 01233 butcher[2][3] = alpha * theta; 01234 01235 butcher[3][0] = 1.; 01236 butcher[3][1] = beta * theta; 01237 butcher[3][2] = alpha * (1.-theta); 01238 butcher[3][3] = (alpha + beta) * theta; 01239 butcher[3][4] = alpha * theta; 01240 01241 butcher[4][1] = beta * theta; 01242 butcher[4][2] = alpha * (1.-theta); 01243 butcher[4][3] = (alpha + beta) * theta; 01244 butcher[4][4] = alpha * theta; 01245 01246 return butcher; 01247 } 01248 else{ 01249 HDNUM_ERROR("Order not available for Runge Kutta solver."); 01250 } 01251 } 01252 01253 static int initOrder(const std::string method) 01254 { 01255 if(method.find("Implicit Euler") != std::string::npos){ 01256 return 1; 01257 } 01258 else if(method.find("Alexander") != std::string::npos){ 01259 return 2; 01260 } 01261 else if(method.find("Crouzieux") != std::string::npos){ 01262 return 3; 01263 } 01264 else if(method.find("Midpoint Rule") != std::string::npos){ 01265 return 2; 01266 } 01267 else if(method.find("Fractional Step Theta") != std::string::npos){ 01268 return 2; 01269 } 01270 else{ 01271 HDNUM_ERROR("Order not available for Runge Kutta solver."); 01272 } 01273 } 01274 01275 01277 // h_n f(t_n, y_n) - y_n + y_{n-1} = 0 01278 class NonlinearProblem 01279 { 01280 public: 01282 typedef typename M::size_type size_type; 01283 01285 typedef typename M::number_type number_type; 01286 01288 NonlinearProblem (const M& model_, const Vector<number_type>& yold_, 01289 typename M::time_type told_, typename M::time_type dt_, 01290 const ButcherTableau & butcher_, const int rk_step_, 01291 const std::vector< Vector<number_type> > & k_) 01292 : model(model_), yold(yold_), told(told_), 01293 dt(dt_), butcher(butcher_), rk_step(rk_step_), k_old(model.size(),0) 01294 { 01295 for(int i=0; i<rk_step; ++i) 01296 k_old.update(butcher[rk_step][1+i] * dt, k_[i]); 01297 } 01298 01300 std::size_t size () const 01301 { 01302 return model.size(); 01303 } 01304 01306 void F (const Vector<number_type>& x, Vector<number_type>& result) const 01307 { 01308 result = k_old; 01309 01310 Vector<number_type> current_z(x); 01311 current_z.update(1.,yold); 01312 01313 const number_type tnew = told + butcher[rk_step][0] * dt; 01314 01315 Vector<number_type> current_k(model.size(),0.); 01316 model.f(tnew,current_z,current_k); 01317 result.update(butcher[rk_step][rk_step+1] * dt, current_k); 01318 01319 result.update(-1.,x); 01320 } 01321 01323 void F_x (const Vector<number_type>& x, DenseMatrix<number_type>& result) const 01324 { 01325 const number_type tnew = told + butcher[rk_step][0] * dt; 01326 01327 Vector<number_type> current_z(x); 01328 current_z.update(1.,yold); 01329 01330 model.f_x(tnew,current_z,result); 01331 01332 result *= dt * butcher[rk_step][rk_step+1]; 01333 01334 for (size_type i=0; i<model.size(); i++) result[i][i] -= number_type(1.0); 01335 } 01336 01337 void set_told_dt (typename M::time_type told_, typename M::time_type dt_) 01338 { 01339 told = told_; 01340 dt = dt_; 01341 } 01342 01343 private: 01344 const M& model; 01345 const Vector<number_type>& yold; 01346 typename M::time_type told; 01347 typename M::time_type dt; 01348 const ButcherTableau & butcher; 01349 const int rk_step; 01350 Vector<number_type> k_old; 01351 }; 01352 01353 public: 01354 01357 DIRK (const M& model_, const S& newton_, const ButcherTableau & butcher_, const int order_) 01358 : verbosity(0), butcher(butcher_), model(model_), newton(newton_), 01359 u(model.size()), order(order_) 01360 { 01361 model.initialize(t,u); 01362 dt = dtmax = 0.1; 01363 } 01364 01367 DIRK (const M& model_, const S& newton_, const std::string method) 01368 : verbosity(0), butcher(initTableau(method)), model(model_), newton(newton_), u(model.size()), 01369 order(initOrder(method)) 01370 { 01371 model.initialize(t,u); 01372 dt = dtmax = 0.1; 01373 } 01374 01375 01377 void set_dt (time_type dt_) 01378 { 01379 dt = dtmax = dt_; 01380 } 01381 01383 void set_verbosity (size_type verbosity_) 01384 { 01385 verbosity = verbosity_; 01386 } 01387 01389 void step () 01390 { 01391 01392 const size_type R = butcher.colsize()-1; 01393 01394 bool reduced = false; 01395 error = false; 01396 if(verbosity>=2) 01397 std::cout << "DIRK: step to" << " t+dt=" << t+dt << " dt=" << dt << std::endl; 01398 01399 while (1) 01400 { 01401 bool converged = true; 01402 01403 // Perform R Runge-Kutta steps 01404 std::vector< Vector<number_type> > k; 01405 for(size_type i=0; i<R; ++i) { 01406 if (verbosity>=2) 01407 std::cout << "DIRK: step nr "<< i << std::endl; 01408 01409 Vector<number_type> current_z(model.size(),0.0); 01410 01411 // Set starting value of k_i 01412 // model.f(t,u,current_k); 01413 01414 // Solve nonlinear problem 01415 NonlinearProblem nlp(model,u,t,dt,butcher,i,k); 01416 01417 newton.solve(nlp,current_z); 01418 01419 converged = converged && newton.has_converged(); 01420 if(!converged) 01421 break; 01422 01423 current_z.update(1., u); 01424 const number_type t_i = t + butcher[i][0] * dt; 01425 Vector<number_type>current_k(model.size(),0.); 01426 model.f(t_i,current_z,current_k); 01427 01428 k.push_back( current_k ); 01429 } 01430 01431 if (converged) 01432 { 01433 if(verbosity >= 2) 01434 std::cout << "DIRK finished"<< std::endl; 01435 01436 // Update to new solution 01437 for(size_type i=0; i<R; ++i) 01438 u.update(dt*butcher[R][1+i],k[i]); 01439 01440 t += dt; 01441 if (!reduced && dt<dtmax-1e-13) 01442 { 01443 dt = std::min(2.0*dt,dtmax); 01444 if (verbosity>0) 01445 std::cout << "DIRK: increasing time step to " << dt << std::endl; 01446 } 01447 return; 01448 } 01449 else 01450 { 01451 if (dt<1e-12) 01452 { 01453 HDNUM_ERROR("time step too small in implicit Euler"); 01454 error = true; 01455 break; 01456 } 01457 dt *= 0.5; 01458 reduced = true; 01459 if (verbosity>0) std::cout << "DIRK: reducing time step to " << dt << std::endl; 01460 } 01461 } 01462 } 01463 01465 bool get_error () const 01466 { 01467 return error; 01468 } 01469 01471 void set_state (time_type t_, const Vector<number_type>& u_) 01472 { 01473 t = t_; 01474 u = u_; 01475 } 01476 01478 const Vector<number_type>& get_state () const 01479 { 01480 return u; 01481 } 01482 01484 time_type get_time () const 01485 { 01486 return t; 01487 } 01488 01490 time_type get_dt () const 01491 { 01492 return dt; 01493 } 01494 01496 size_type get_order () const 01497 { 01498 return order; 01499 } 01500 01502 void get_info () const 01503 { 01504 } 01505 01506 private: 01507 size_type verbosity; 01508 const DenseMatrix<number_type> butcher; 01509 const M& model; 01510 const S& newton; 01511 time_type t, dt, dtmax; 01512 number_type reduction; 01513 size_type linesearchsteps; 01514 Vector<number_type> u; 01515 int order; 01516 mutable bool error; 01517 }; 01518 01520 template<class T, class N> 01521 inline void gnuplot (const std::string& fname, const std::vector<T> t, const std::vector<Vector<N> > u) 01522 { 01523 std::fstream f(fname.c_str(),std::ios::out); 01524 for (typename std::vector<T>::size_type n=0; n<t.size(); n++) 01525 { 01526 f << std::scientific << std::showpoint 01527 << std::setprecision(16) << t[n]; 01528 for (typename Vector<N>::size_type i=0; i<u[n].size(); i++) 01529 f << " " << std::scientific << std::showpoint 01530 << std::setprecision(u[n].precision()) << u[n][i]; 01531 f << std::endl; 01532 } 01533 f.close(); 01534 } 01535 01537 template<class T, class N> 01538 inline void gnuplot (const std::string& fname, const std::vector<T> t, const std::vector<Vector<N> > u, const std::vector<T> dt) 01539 { 01540 std::fstream f(fname.c_str(),std::ios::out); 01541 for (typename std::vector<T>::size_type n=0; n<t.size(); n++) 01542 { 01543 f << std::scientific << std::showpoint 01544 << std::setprecision(16) << t[n]; 01545 for (typename Vector<N>::size_type i=0; i<u[n].size(); i++) 01546 f << " " << std::scientific << std::showpoint 01547 << std::setprecision(u[n].precision()) << u[n][i]; 01548 f << " " << std::scientific << std::showpoint 01549 << std::setprecision(16) << dt[n]; 01550 f << std::endl; 01551 } 01552 f.close(); 01553 } 01554 01555 } // namespace hdnum 01556 01557 #endif