23 #ifndef O2SCL_ODE_IV_SOLVE_H 24 #define O2SCL_ODE_IV_SOLVE_H 32 #include <boost/numeric/ublas/vector.hpp> 33 #include <boost/numeric/ublas/matrix.hpp> 34 #include <boost/numeric/ublas/matrix_proxy.hpp> 36 #include <o2scl/astep.h> 37 #include <o2scl/astep_gsl.h> 39 #ifndef DOXYGEN_NO_O2NS 97 #ifndef DOXYGEN_INTERNAL 103 vec_t vtemp, vtemp2, vtemp3, vtemp4;
114 std::cout <<
type() <<
" x: " << x <<
" y: ";
115 for(
size_t i=0;i<nv;i++) std::cout << y[i] <<
" ";
116 std::cout << std::endl;
182 vec_t &ystart, vec_t ¥d, func_t &derivs) {
205 vec_t &ystart, vec_t ¥d, vec_t &yerr,
240 vec_t &ystart, vec_t ¥d, vec_t &yerr,
241 vec_t &dydx_end, func_t &derivs) {
243 if ((x1>x0 && h<=0.0) || (x0>x1 && h>=0.0)) {
244 std::string str=
"Interval direction (x1-x0="+
o2scl::dtos(x1-x0)+
245 ") does not match step direction (h="+
o2scl::dtos(h)+
246 " in ode_iv_solve::solve_final_value().";
250 O2SCL_ERR2(
"Starting and final endpoints identical in ",
251 "ode_iv_solve::solve_final_value().",
exc_einval);
261 double x=x0, xverb=0.0, dxverb=0.0, xnext;
262 int ret=0, first_ret=0;
267 vec_t &dydx_start=vtemp;
287 for(
size_t i=0;i<n;i++) yend[i]=ystart[i];
290 derivs(x,n,yend,dydx_start);
293 while (done==
false) {
296 ret=astp->
astep_full(x,x1,xnext,h,n,ystart,dydx_start,
297 yend,yerr,dydx_end,derivs);
302 "ode_iv_solve::solve_final_value()",ret);
303 }
else if (first_ret!=0) {
309 if (
verbose>0 && xnext>xverb) {
321 std::string str=
"Too many steps required (ntrial="+
itos(
ntrial)+
322 ", x="+
o2scl::dtos(x)+
") in ode_iv_solve::solve_final_value().";
330 if (xnext>=x1) done=
true;
332 if (xnext<=x1) done=
true;
339 ret=astp->
astep_full(xnext,x1,x,h,n,yend,dydx_end,ystart,yerr,
345 "ode_iv_solve::solve_final_value()",ret);
346 }
else if (first_ret!=0) {
364 std::string str=
"Too many steps required (ntrial="+
itos(
ntrial)+
365 ", x="+
o2scl::dtos(x)+
") in ode_iv_solve::solve_final_value().";
375 for(
size_t j=0;j<n;j++) {
377 dydx_end[j]=dydx_start[j];
383 for(
size_t j=0;j<n;j++) {
385 dydx_end[j]=dydx_start[j];
437 template<
class mat_t>
439 vec_t &x_sol, mat_t &y_sol, mat_t &yerr_sol,
440 mat_t &dydx_sol, func_t &derivs,
size_t istart=0) {
449 double dx_verb=(x1-x0)/((
double)(
nsteps_out-1));
451 double dx_tab=(x1-x0)/((
double)(n_sol-istart-1));
453 double x_verb=x0+dx_verb;
454 double x_tab=x0+dx_tab;
461 vec_t &ystart=vtemp4;
462 vec_t &dydx_start=vtemp2;
465 for(
size_t j=0;j<n;j++) ystart[j]=y_sol(istart,j);
477 derivs(x0,n,ystart,dydx_start);
481 for(
size_t j=0;j<n;j++) {
482 dydx_sol(istart,j)=dydx_start[j];
483 yerr_sol(istart,j)=0.0;
487 size_t icurr=istart+1;
489 for(
size_t j=0;j<n;j++) y_sol(icurr,j)=ystart[j];
494 while (done==
false) {
499 vec_t &dydx_out=vtemp3;
504 for(
size_t i=0;i<n;i++) yrow[i]=y_sol(icurr,i);
505 ret=astp->
astep_full(x0,x1,xnext,h,n,yrow,dydx_start,ystart,yerr,
512 O2SCL_ERR2(
"Adaptive stepper returned non-zero in ",
514 }
else if (first_ret==0) {
522 std::string str=
"Too many steps required (ntrial="+
itos(
ntrial)+
523 ", x="+
o2scl::dtos(x0)+
") in ode_iv_solve::solve_store().";
528 if (xnext>=x_verb &&
verbose>0) {
545 for(
size_t j=0;j<n;j++) {
546 y_sol(icurr,j)=ystart[j];
547 dydx_sol(icurr,j)=dydx_out[j];
548 yerr_sol(icurr,j)=yerr[j];
556 if (xnext>=x_tab && icurr<nmax) {
560 for(
size_t j=0;j<n;j++) {
561 y_sol(icurr,j)=ystart[j];
562 dydx_sol(icurr,j)=dydx_out[j];
563 yerr_sol(icurr,j)=yerr[j];
566 y_sol(icurr+1,j)=ystart[j];
567 dydx_start[j]=dydx_out[j];
580 for(
size_t j=0;j<n;j++) {
581 dydx_start[j]=dydx_out[j];
582 y_sol(icurr,j)=ystart[j];
611 template<
class mat_t,
class mat_row_t>
613 mat_t &ysol, mat_t &err_sol, mat_t &dydx_sol,
617 double x1=xsol[nsol-1];
620 int ret=0, first_ret=0;
632 vec_t &dydx_start=vtemp2;
635 for(
size_t j=0;j<n;j++) {
643 derivs(x0,n,ystart,dydx_start);
646 for(
size_t j=0;j<n;j++) {
647 dydx_sol(0,j)=dydx_start[j];
651 for(
size_t i=1;i<nsol && ret==0;i++) {
653 mat_row_t y_row(ysol,i);
654 mat_row_t dydx_row(dydx_sol,i);
655 mat_row_t yerr_row(err_sol,i);
659 while(done==
false && ret==0) {
661 ret=astp->
astep_full(x,xsol[i],xnext,h,n,ystart,dydx_start,
662 y_row,yerr_row,dydx_row,derivs);
670 "ode_iv_solve::solve_grid()",ret);
671 }
else if (first_ret!=0) {
677 std::string str=
"Too many steps required (ntrial="+
itos(
ntrial)+
678 ", x="+
o2scl::dtos(x)+
") in ode_iv_solve::solve_grid().";
685 for(
size_t j=0;j<n;j++) {
687 dydx_start[j]=dydx_row[j];
692 if (x>=xsol[i]) done=
true;
694 if (x<=xsol[i]) done=
true;
745 virtual const char *
type() {
return "ode_iv_solve"; }
749 #ifndef DOXYGEN_NO_O2NS Solve an initial-value ODE problems given an adaptive ODE stepper.
astep_base< vec_t, vec_t, vec_t, func_t > * astp
The adaptive stepper.
The main O<span style='position: relative; top: 0.3em; font-size: 0.8em'>2</span>scl O$_2$scl names...
astep_gsl< vec_t, vec_t, vec_t, func_t > gsl_astp
The default adaptive stepper.
invalid argument supplied by user
virtual const char * type()
Return the type, "ode_iv_solve".
size_t nsteps
Number of adaptive ste!ps employed.
exceeded max number of iterations
bool exit_on_fail
If true, stop the solution if the adaptive stepper fails (default true)
int verbose
Set output level.
int solve_final_value(double x0, double x1, double h, size_t n, vec_t &ystart, vec_t ¥d, vec_t &yerr, func_t &derivs)
Solve the initial-value problem to get the final value with errors.
int solve_final_value(double x0, double x1, double h, size_t n, vec_t &ystart, vec_t ¥d, vec_t &yerr, vec_t &dydx_end, func_t &derivs)
Solve the initial-value problem to get the final value, derivatives, and errors.
int solve_final_value(double x0, double x1, double h, size_t n, vec_t &ystart, vec_t ¥d, func_t &derivs)
Solve the initial-value problem to get the final value.
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
std::string dtos(double x, int prec=6, bool auto_prec=false)
Convert a double to a string.
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
virtual int astep_full(double x, double xlimit, double &x_out, double &h, size_t n, vec_y_t &y, vec_dydx_t &dydx, vec_y_t &yout, vec_yerr_t &yerr, vec_dydx_t &dydx_out, func_t &derivs)=0
Make an adaptive integration step of the system derivs with derivatives.
size_t ntrial
Maximum number of applications of the adaptive stepper (default 1000)
int set_astep(astep_base< vec_t, vec_t, vec_t, func_t > &as)
Set the adaptive stepper to use.
int solve_store(double x0, double x1, double h, size_t n, size_t &n_sol, vec_t &x_sol, mat_t &y_sol, mat_t &yerr_sol, mat_t &dydx_sol, func_t &derivs, size_t istart=0)
Solve the initial-value problem and store the associated output.
void allocate(size_t n)
Allocate space for temporary vectors.
virtual int print_iter(double x, size_t nv, vec_t &y)
Print out iteration information.
void free()
Free allocated memory.
int solve_grid(double h, size_t n, size_t nsol, vec_t &xsol, mat_t &ysol, mat_t &err_sol, mat_t &dydx_sol, func_t &derivs)
Solve the initial-value problem from x0 to x1 over a grid storing derivatives and errors...
size_t nsteps_out
Number of output points for verbose output (default 10)
std::string itos(int x)
Convert an integer to a string.
static const double x1[5]
std::function< int(double, size_t, const boost::numeric::ublas::vector< double > &, boost::numeric::ublas::vector< double > &)> ode_funct11
Ordinary differential equation function.