00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_ODE_IV_SOLVE_H 00024 #define O2SCL_ODE_IV_SOLVE_H 00025 00026 #include <string> 00027 #include <o2scl/collection.h> 00028 #include <o2scl/adapt_step.h> 00029 #include <o2scl/gsl_astep.h> 00030 00031 #ifndef DOXYGENP 00032 namespace o2scl { 00033 #endif 00034 00035 /** 00036 \brief Solve an initial-value ODE problems given an 00037 adaptive ODE stepper 00038 00039 \future Add error information 00040 */ 00041 template<class param_t, class func_t=ode_funct<param_t>, 00042 class vec_t=ovector_base, class alloc_vec_t=ovector, 00043 class alloc_t=ovector_alloc, class mat_row_t=omatrix_row> 00044 class ode_iv_solve { 00045 00046 #ifndef DOXYGEN_INTERNAL 00047 00048 protected: 00049 00050 /// Derivative 00051 alloc_vec_t dydx; 00052 00053 /// Memory allocator 00054 alloc_t ao; 00055 00056 /// The adaptive stepper 00057 adapt_step<param_t,func_t,vec_t,alloc_vec_t,alloc_t> *astp; 00058 00059 /// Print out iteration information 00060 virtual int print_iter(double x, size_t nv, vec_t &y) { 00061 std::cout << type() << " x: " << x << " y: "; 00062 for(size_t i=0;i<nv;i++) std::cout << y[i] << " "; 00063 std::cout << std::endl; 00064 if (verbose>1) { 00065 char ch; 00066 std::cin >> ch; 00067 } 00068 return 0; 00069 } 00070 00071 #endif 00072 00073 public: 00074 00075 ode_iv_solve() { 00076 verbose=0; 00077 ntrial=1000; 00078 nsteps_out=10; 00079 astp=&gsl_astp; 00080 exit_on_fail=true; 00081 } 00082 00083 virtual ~ode_iv_solve() {} 00084 00085 /// Set the adaptive stepper to use 00086 int set_adapt_step(adapt_step<param_t,func_t,vec_t, 00087 alloc_vec_t,alloc_t> &as) { 00088 astp=&as; 00089 return 0; 00090 } 00091 00092 /** 00093 \brief Solve the initial-value problem and output a table 00094 00095 Initially, \c xsol should be a vector of size \c nsol, and \c 00096 ysol should be a two-dimensional array (e.g. omatrix) of 00097 size \c [nsol][n]. On exit, \c nsol will will be the size of 00098 the solution table, less than or equal to the original value 00099 of \c nsol. 00100 00101 If \ref verbose is greater than zero, The solution 00102 at each internal point will be written to \c std::cout. 00103 If \ref verbose is greater than one, a character 00104 will be required after each point. 00105 00106 If the given value of \c h is small enough, the solution 00107 may generate more points than the space initially allocated 00108 and the full solution will not be generated. 00109 00110 \future Should we give the option not to call the error 00111 handler in case running out of space in the table? 00112 \future Consider modifying so that this can handle tables 00113 which are too small by removing half the rows and doubling 00114 the stepsize. Alternatively, make a function which 00115 accepts allocation objects and can re-allocate the vector 00116 and matrix space if required? (This latter option might 00117 be hard for arrays since they are really of fixed size.) 00118 A final alternative is to offer an ovector/omatrix 00119 interface which does the reallocation automatically 00120 if necessary. 00121 \future Some copying is done here and it would be 00122 nice if that could be avoided. 00123 */ 00124 template<class mat_t> 00125 int solve_table(double x0, double x1, double h, size_t n, 00126 vec_t &ystart, size_t &nsol, vec_t &xsol, 00127 mat_t &ysol, param_t &pa, func_t &derivs) { 00128 00129 int ret=0, nmax=nsol, i, first_ret=0; 00130 size_t j; 00131 00132 xsol[0]=x0; 00133 for(j=0;j<n;j++) ysol[0][j]=ystart[j]; 00134 if (verbose>0) { 00135 print_iter(xsol[0],n,ystart); 00136 if (verbose>1) { 00137 char ch; 00138 std::cin >> ch; 00139 } 00140 } 00141 00142 alloc_vec_t yerrl; 00143 ao.allocate(yerrl,n); 00144 ao.allocate(dydx,n); 00145 derivs(x0,n,ystart,dydx,pa); 00146 00147 for(i=1;i<nmax && xsol[i-1]<x1;i++) { 00148 00149 xsol[i]=xsol[i-1]; 00150 for(j=0;j<n;j++) ysol[i][j]=ysol[i-1][j]; 00151 00152 mat_row_t ar(ysol,i); 00153 00154 ret=astp->astep_derivs(xsol[i],h,x1,n,ar,dydx,yerrl,pa,derivs); 00155 if (ret!=0) { 00156 if (exit_on_fail) { 00157 nsol=i+1; 00158 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00159 "ode_iv_solve::solve_table()",ret); 00160 } else if (first_ret!=0) { 00161 first_ret=ret; 00162 } 00163 } 00164 00165 if (verbose>0) { 00166 print_iter(xsol[i],n,ar); 00167 if (verbose>1) { 00168 char ch; 00169 std::cin >> ch; 00170 } 00171 } 00172 } 00173 00174 ao.free(dydx); 00175 ao.free(yerrl); 00176 00177 nsol=i; 00178 00179 if (i==nmax && xsol[i-1]<x1) { 00180 O2SCL_ERR2_RET("Ran out of space in ", 00181 "ode_iv_solve::solve_table().",gsl_etable); 00182 } 00183 return first_ret; 00184 } 00185 00186 /** \brief Solve the initial-value problem from \c x0 00187 to \c x1 over a grid 00188 00189 Initially, \c xsol should be an array of size \c nsol, and 00190 \c ysol should be a \c matrix of size \c [nsol][n]. This 00191 function never takes a step larger than the grid size, 00192 but will take a step smaller than the grid size in order 00193 to ensure accuracy. 00194 00195 If \ref verbose is greater than zero, The at each grid point 00196 will be written to \c std::cout. If \ref verbose is greater 00197 than one, a character will be required after each point. 00198 */ 00199 template<class mat_t> 00200 int solve_grid(double x0, double x1, double h, size_t n, 00201 vec_t &ystart, size_t nsol, vec_t &xsol, mat_t &ysol, 00202 param_t &pa, func_t &derivs) { 00203 00204 double x=x0, xnext; 00205 int ret=0, first_ret=0; 00206 size_t j; 00207 00208 xsol[0]=x0; 00209 for(j=0;j<n;j++) ysol[0][j]=ystart[j]; 00210 if (verbose>0) print_iter(xsol[0],n,ystart); 00211 00212 alloc_vec_t yerrl; 00213 ao.allocate(dydx,n); 00214 ao.allocate(yerrl,n); 00215 00216 if (x0<x1) { 00217 00218 derivs(x0,n,ystart,dydx,pa); 00219 00220 for(size_t i=1;i<nsol && ret==0;i++) { 00221 00222 xsol[i]=xsol[i-1]; 00223 for(j=0;j<n;j++) ysol[i][j]=ysol[i-1][j]; 00224 00225 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00226 00227 mat_row_t ar(ysol,i); 00228 while(xsol[i]<xnext && ret==0) { 00229 00230 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,dydx,yerrl,pa,derivs); 00231 if (ret!=0) { 00232 if (exit_on_fail) { 00233 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00234 "ode_iv_solve::solve_grid()",ret); 00235 } else if (first_ret!=0) { 00236 first_ret=ret; 00237 } 00238 } 00239 } 00240 if (verbose>0) print_iter(xsol[i],n,ar); 00241 00242 } 00243 00244 } else { 00245 00246 derivs(x0,n,ystart,dydx,pa); 00247 00248 for(size_t i=1;i<nsol && ret==0;i++) { 00249 00250 xsol[i]=xsol[i-1]; 00251 for(j=0;j<n;j++) ysol[i][j]=ysol[i-1][j]; 00252 00253 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00254 00255 mat_row_t ar(ysol,i); 00256 while(xsol[i]>xnext && ret==0) { 00257 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,dydx,yerrl,pa,derivs); 00258 if (ret!=0) { 00259 if (exit_on_fail) { 00260 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00261 "ode_iv_solve::solve_grid()",ret); 00262 } else if (first_ret!=0) { 00263 first_ret=ret; 00264 } 00265 } 00266 } 00267 if (verbose>0) print_iter(xsol[i],n,ar); 00268 00269 } 00270 00271 } 00272 00273 ao.free(dydx); 00274 ao.free(yerrl); 00275 00276 return first_ret; 00277 } 00278 00279 /** \brief Solve the initial-value problem from \c x0 00280 to \c x1 over a grid storing derivatives 00281 00282 Initially, \c xsol should be an array of size \c nsol, and 00283 \c ysol should be a \c omatrix of size \c [nsol][n]. This 00284 function never takes a step larger than the grid size. 00285 00286 If \ref verbose is greater than zero, The solution 00287 at each grid point will be written to \c std::cout. 00288 If \ref verbose is greater than one, a character 00289 will be required after each point. 00290 00291 \todo Add error information 00292 */ 00293 template<class mat_t> 00294 int solve_grid_derivs(double x0, double x1, double h, size_t n, 00295 vec_t &ystart, size_t nsol, vec_t &xsol, 00296 mat_t &ysol, mat_t &dydx_sol, 00297 param_t &pa, func_t &derivs) { 00298 00299 double x=x0, xnext; 00300 int ret=0, first_ret=0; 00301 size_t j; 00302 00303 xsol[0]=x0; 00304 for(j=0;j<n;j++) ysol[0][j]=ystart[j]; 00305 if (verbose>0) print_iter(xsol[0],n,ystart); 00306 00307 alloc_vec_t yerrl; 00308 ao.allocate(yerrl,n); 00309 00310 if (x0<x1) { 00311 00312 mat_row_t mr0(dydx_sol,0); 00313 derivs(x0,n,ystart,mr0,pa); 00314 00315 for(size_t i=1;i<nsol && ret==0;i++) { 00316 00317 xsol[i]=xsol[i-1]; 00318 for(j=0;j<n;j++) { 00319 ysol[i][j]=ysol[i-1][j]; 00320 dydx_sol[i][j]=dydx_sol[i-1][j]; 00321 } 00322 00323 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00324 00325 mat_row_t ar(ysol,i); 00326 mat_row_t ar2(dydx_sol,i); 00327 00328 while(xsol[i]<xnext && ret==0) { 00329 00330 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,ar2,yerrl,pa,derivs); 00331 if (ret!=0) { 00332 if (exit_on_fail) { 00333 ao.free(yerrl); 00334 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00335 "ode_iv_solve::solve_grid()",ret); 00336 } else if (first_ret!=0) { 00337 first_ret=ret; 00338 } 00339 } 00340 } 00341 if (verbose>0) print_iter(xsol[i],n,ar); 00342 00343 } 00344 00345 } else { 00346 00347 mat_row_t mr0(dydx_sol,0); 00348 derivs(x0,n,ystart,dydx,pa); 00349 00350 for(size_t i=1;i<nsol && ret==0;i++) { 00351 00352 xsol[i]=xsol[i-1]; 00353 for(j=0;j<n;j++) { 00354 ysol[i][j]=ysol[i-1][j]; 00355 dydx_sol[i][j]=dydx_sol[i-1][j]; 00356 } 00357 00358 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00359 00360 mat_row_t ar(ysol,i); 00361 mat_row_t ar2(dydx_sol,i); 00362 00363 while(xsol[i]>xnext && ret==0) { 00364 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,ar2,yerrl,pa,derivs); 00365 if (ret!=0) { 00366 if (exit_on_fail) { 00367 ao.free(yerrl); 00368 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00369 "ode_iv_solve::solve_grid()",ret); 00370 } else if (first_ret!=0) { 00371 first_ret=ret; 00372 } 00373 } 00374 } 00375 if (verbose>0) print_iter(xsol[i],n,ar); 00376 00377 } 00378 00379 } 00380 00381 ao.free(yerrl); 00382 00383 return first_ret; 00384 } 00385 00386 #ifdef O2SCL_NEVER_DEFINED 00387 00388 /** \brief Solve the initial-value problem from \c x0 00389 to \c x1 over a grid storing derivatives 00390 00391 Initially, \c xsol should be an array of size \c nsol, and 00392 \c ysol should be a \c omatrix of size \c [nsol][n]. This 00393 function never takes a step larger than the grid size. 00394 00395 If \ref verbose is greater than zero, The solution 00396 at each grid point will be written to \c std::cout. 00397 If \ref verbose is greater than one, a character 00398 will be required after each point. 00399 00400 \todo Add error information 00401 */ 00402 template<class mat_t> 00403 int solve_grid_derivs_err(double x0, double x1, double h, size_t n, 00404 vec_t &ystart, size_t nsol, vec_t &xsol, 00405 mat_t &ysol, mat_t &dydx_sol, mat_t &errsol, 00406 param_t &pa, func_t &derivs) { 00407 00408 double x=x0, xnext; 00409 int ret=0, first_ret=0; 00410 size_t j; 00411 00412 xsol[0]=x0; 00413 for(j=0;j<n;j++) ysol[0][j]=ystart[j]; 00414 if (verbose>0) print_iter(xsol[0],n,ystart); 00415 00416 if (x0<x1) { 00417 00418 mat_row_t mr0(dydx_sol,0); 00419 derivs(x0,n,ystart,mr0,pa); 00420 00421 for(size_t i=1;i<nsol && ret==0;i++) { 00422 00423 xsol[i]=xsol[i-1]; 00424 for(j=0;j<n;j++) { 00425 ysol[i][j]=ysol[i-1][j]; 00426 dydx_sol[i][j]=dydx_sol[i-1][j]; 00427 } 00428 00429 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00430 00431 mat_row_t ar(ysol,i); 00432 mat_row_t ar2(dydx_sol,i); 00433 mat_row_t aerr(errsol,i); 00434 00435 while(xsol[i]<xnext && ret==0) { 00436 00437 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,ar2,aerr,pa,derivs); 00438 if (ret!=0) { 00439 if (exit_on_fail) { 00440 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00441 "ode_iv_solve::solve_grid()",ret); 00442 } else if (first_ret!=0) { 00443 first_ret=ret; 00444 } 00445 } 00446 } 00447 if (verbose>0) print_iter(xsol[i],n,ar); 00448 00449 } 00450 00451 } else { 00452 00453 mat_row_t mr0(dydx_sol,0); 00454 derivs(x0,n,ystart,dydx,pa); 00455 00456 for(size_t i=1;i<nsol && ret==0;i++) { 00457 00458 xsol[i]=xsol[i-1]; 00459 for(j=0;j<n;j++) { 00460 ysol[i][j]=ysol[i-1][j]; 00461 dydx_sol[i][j]=dydx_sol[i-1][j]; 00462 } 00463 00464 xnext=x0+(x1-x0)*((double)i)/((double)(nsol-1)); 00465 00466 mat_row_t ar(ysol,i); 00467 mat_row_t ar2(dydx_sol,i); 00468 mat_row_t aerr(errsol,i); 00469 00470 while(xsol[i]>xnext && ret==0) { 00471 ret=astp->astep_derivs(xsol[i],h,xnext,n,ar,ar2,aerr,pa,derivs); 00472 if (ret!=0) { 00473 if (exit_on_fail) { 00474 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00475 "ode_iv_solve::solve_grid()",ret); 00476 } else if (first_ret!=0) { 00477 first_ret=ret; 00478 } 00479 } 00480 } 00481 if (verbose>0) print_iter(xsol[i],n,ar); 00482 00483 } 00484 00485 } 00486 00487 return first_ret; 00488 } 00489 00490 #endif 00491 00492 /** \brief Solve the initial-value problem 00493 to get the final value 00494 00495 If \ref verbose is greater than zero, The solution at less 00496 than or approximately equal to \ref nsteps_out points will be 00497 written to \c std::cout. If \ref verbose is greater than one, 00498 a character will be required after each selected point. 00499 00500 The solution fails if more than \ref ntrial steps are 00501 required. 00502 */ 00503 int solve_final_value(double x0, double x1, double h, size_t n, 00504 vec_t &ystart, vec_t ¥d, 00505 param_t &pa, func_t &derivs) { 00506 00507 ao.allocate(dydx,n); 00508 int ret=solve_final_value_derivs(x0,x1,h,n,ystart,yend,dydx, 00509 dydx,pa,derivs); 00510 ao.free(dydx); 00511 return ret; 00512 } 00513 00514 /** \brief Solve the initial-value problem 00515 to get the final value and derivative 00516 00517 If \ref verbose is greater than zero, The solution at less 00518 than or approximately equal to \ref nsteps_out points will be 00519 written to \c std::cout. If \ref verbose is greater than one, 00520 a character will be required after each selected point. 00521 00522 The solution fails if more than \ref ntrial steps are 00523 required. 00524 00525 This function will also fail if <tt>x1>x0</tt> and <tt>h<0</tt> 00526 or if <tt>x1<x0</tt> but <tt>h>0</tt>. 00527 00528 \todo Add error information 00529 \todo At present, dydx_start is computed, but it should 00530 probably be assumed that the user will provide this. 00531 \future It looks like the x0<x1 code and the x1<x0 code 00532 is the same? 00533 */ 00534 int solve_final_value_derivs(double x0, double x1, double h, size_t n, 00535 vec_t &ystart, vec_t ¥d, 00536 vec_t &dydx_start, vec_t &dydx_end, 00537 param_t &pa, func_t &derivs) { 00538 00539 if ((x1>x0 && h<=0.0) || (x0>x1 && h>=0.0)) { 00540 O2SCL_ERR_RET("Interval ordering and sign of h don't match.", 00541 gsl_einval); 00542 } 00543 00544 double x=x0, xverb=0.0, dxverb=0.0; 00545 int ret=0, it=0, first_ret=0; 00546 00547 // Compute stepsize for verbose output 00548 if (verbose>0) { 00549 print_iter(x0,n,ystart); 00550 if (verbose>1) { 00551 char ch; 00552 std::cin >> ch; 00553 } 00554 if (x1>x0) { 00555 dxverb=(x1-x0)/((double)nsteps_out); 00556 xverb=x0+dxverb; 00557 } else { 00558 dxverb=(x0-x1)/((double)nsteps_out); 00559 xverb=x0-dxverb; 00560 } 00561 } 00562 00563 // Use yend as workspace 00564 for(size_t i=0;i<n;i++) yend[i]=ystart[i]; 00565 00566 // Allocate space for the errors 00567 alloc_vec_t yerrl; 00568 ao.allocate(yerrl,n); 00569 00570 if (x1>x0) { 00571 00572 // Compute initial derivative 00573 derivs(x,n,yend,dydx_start,pa); 00574 for(size_t i=0;i<n;i++) dydx_end[i]=dydx_start[i]; 00575 00576 // Main loop 00577 while(x<x1) { 00578 00579 // Take a step 00580 ret=astp->astep_derivs(x,h,x1,n,yend,dydx_end,yerrl,pa,derivs); 00581 if (ret!=0) { 00582 if (exit_on_fail) { 00583 ao.free(yerrl); 00584 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00585 "ode_iv_solve::solve_grid()",ret); 00586 } else if (first_ret!=0) { 00587 first_ret=ret; 00588 } 00589 } 00590 00591 // Print out verbose info 00592 if (verbose>0 && x>xverb) { 00593 print_iter(x,n,yend); 00594 if (verbose>1) { 00595 char ch; 00596 std::cin >> ch; 00597 } 00598 xverb+=dxverb; 00599 } 00600 00601 // Check number of iterations 00602 it++; 00603 if (it>ntrial) { 00604 ao.free(yerrl); 00605 O2SCL_ERR2_RET("Exceeded max number of iterations in ", 00606 "ode_iv_solve::solve_final_value().", 00607 gsl_emaxiter); 00608 } 00609 } 00610 00611 } else { 00612 00613 // Compute initial derivative 00614 derivs(x,n,yend,dydx_start,pa); 00615 for(size_t i=0;i<n;i++) dydx_end[i]=dydx_start[i]; 00616 00617 // Main loop 00618 while(x>x1 && ret==0) { 00619 00620 // Take a step 00621 ret=astp->astep_derivs(x,h,x1,n,yend,dydx_end,yerrl,pa,derivs); 00622 if (ret!=0) { 00623 if (exit_on_fail) { 00624 ao.free(yerrl); 00625 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00626 "ode_iv_solve::solve_grid()",ret); 00627 } else if (first_ret!=0) { 00628 first_ret=ret; 00629 } 00630 } 00631 00632 // Print out verbose info 00633 if (verbose>0 && x<xverb) { 00634 print_iter(x,n,yend); 00635 if (verbose>1) { 00636 char ch; 00637 std::cin >> ch; 00638 } 00639 xverb-=dxverb; 00640 } 00641 00642 // Check number of iterations 00643 it++; 00644 if (it>ntrial) { 00645 ao.free(yerrl); 00646 O2SCL_ERR2_RET("Exceeded max number of iterations in ", 00647 "ode_iv_solve::solve_final_value().", 00648 gsl_emaxiter); 00649 } 00650 } 00651 } 00652 00653 // Print out final verbose info 00654 if (verbose>0) { 00655 print_iter(x,n,yend); 00656 if (verbose>1) { 00657 char ch; 00658 std::cin >> ch; 00659 } 00660 } 00661 00662 ao.free(yerrl); 00663 return first_ret; 00664 } 00665 00666 /// Set output level 00667 int verbose; 00668 00669 /** 00670 \brief Number of output points if \ref verbose is greater 00671 than zero (default 10) 00672 */ 00673 size_t nsteps_out; 00674 00675 /// Maximum number of steps for solve_final_value() (default 1000) 00676 int ntrial; 00677 00678 /** 00679 \brief If true, stop the solution if the adaptive stepper fails 00680 00681 If this is false, then failures in the adaptive stepper are 00682 ignored. 00683 */ 00684 bool exit_on_fail; 00685 00686 /// The default adaptive stepper 00687 gsl_astep<param_t,func_t,vec_t,alloc_vec_t,alloc_t> gsl_astp; 00688 00689 /// Return the type, \c "ode_iv_solve". 00690 virtual const char *type() { return "ode_iv_solve"; } 00691 00692 }; 00693 00694 #ifdef O2SCL_NEVER_DEFINED 00695 00696 /** 00697 \brief Solve an initial-value ODE problems given an 00698 adaptive ODE stepper 00699 */ 00700 template<class param_t,func_t=ode_funct<param_t> > 00701 class ode_iv_solve_x : public ode_iv_solve<param_t,func_t, 00702 ovector_base,ovector,ovector_alloc,omatrix_row> { 00703 00704 public: 00705 00706 int solve_x_derivs(double x0, double x1, double h, size_t n, 00707 ovector_base &ystart, size_t &nsol, 00708 ovector &xsol, std::vector<ovector> &ysol, 00709 std::vector<ovector> &dydx_sol, 00710 param_t &pa, func_t &derivs) { 00711 00712 // Temporary ovectors 00713 ovector y, dy; 00714 00715 // Proceed with solution 00716 double x=x0, xnext; 00717 int ret=0, first_ret=0; 00718 size_t j; 00719 00720 xsol.push_back(x); 00721 for(j=0;j<n;j++) y[j]=ystart[j]; 00722 ysol.push_back(y); 00723 00724 alloc_vec_t yerrl; 00725 ao.allocate(yerrl,n); 00726 00727 if (x0<x1) { 00728 00729 derivs(x,n,y,dy,pa); 00730 dydx_sol.push_back(dy); 00731 00732 while (x<x1 && ret==0) { 00733 00734 ret=astp->astep_derivs(x,h,x1,n,y,dy,yerrl,pa,derivs); 00735 if (ret!=0) { 00736 if (exit_on_fail) { 00737 ao.free(yerrl); 00738 O2SCL_ERR2_RET("Adaptive stepper failed in ", 00739 "ode_iv_solve::solve_grid()",ret); 00740 } else if (first_ret!=0) { 00741 first_ret=ret; 00742 } 00743 } 00744 if (verbose>0) print_iter(xsol[i],n,ar); 00745 00746 xsol.push_back(x); 00747 ysol.push_back(y); 00748 ysol.push_back(dy); 00749 00750 00751 } 00752 00753 } 00754 00755 ao.free(yerrl); 00756 00757 return first_ret; 00758 } 00759 00760 #endif 00761 00762 00763 #ifndef DOXYGENP 00764 } 00765 #endif 00766 00767 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page