1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467
|
/* DOP853
------
This code computes the numerical solution of a system of first order ordinary
differential equations y'=f(x,y). It uses an explicit Runge-Kutta method of
order 8(5,3) due to Dormand & Prince with step size control and dense output.
Authors : E. Hairer & G. Wanner
Universite de Geneve, dept. de Mathematiques
CH-1211 GENEVE 4, SWITZERLAND
E-mail : HAIRER@DIVSUN.UNIGE.CH, WANNER@DIVSUN.UNIGE.CH
The code is described in : E. Hairer, S.P. Norsett and G. Wanner, Solving
ordinary differential equations I, nonstiff problems, 2nd edition,
Springer Series in Computational Mathematics, Springer-Verlag (1993).
Version of Mai 2, 1994.
Remarks about the C version : this version allocates memory by itself, the
iwork array (among the initial FORTRAN parameters) has been splitted into
independant initial parameters, the statistical variables and last step size
and x have been encapsulated in the module and are now accessible through
dedicated functions; the variable names have been kept to maintain a kind
of reading compatibility between the C and FORTRAN codes; adaptation made by
J.Colinge (COLINGE@DIVSUN.UNIGE.CH).
INPUT PARAMETERS
----------------
n Dimension of the system (n < UINT_MAX).
fcn A pointer the the function definig the differential equation, this
function must have the following prototype
void fcn (unsigned n, double x, double *y, double *f)
where the array f will be filled with the function result.
x Initial x value.
*y Initial y values (double y[n]).
xend Final x value (xend-x may be positive or negative).
*rtoler Relative and absolute error tolerances. They can be both scalars or
*atoler vectors of length n (in the scalar case pass the addresses of
variables where you have placed the tolerance values).
itoler Switch for atoler and rtoler :
itoler=0 : both atoler and rtoler are scalars, the code keeps
roughly the local error of y[i] below
rtoler*abs(y[i])+atoler.
itoler=1 : both rtoler and atoler are vectors, the code keeps
the local error of y[i] below
rtoler[i]*abs(y[i])+atoler[i].
solout A pointer to the output function called during integration.
If iout >= 1, it is called after every successful step. If iout = 0,
pass a pointer equal to NULL. solout must must have the following
prototype
solout (long nr, double xold, double x, double* y, unsigned n, int* irtrn)
where y is the solution the at nr-th grid point x, xold is the
previous grid point and irtrn serves to interrupt the integration
(if set to a negative value).
Continuous output : during the calls to solout, a continuous solution
for the interval (xold,x) is available through the function
contd8(i,s)
which provides an approximation to the i-th component of the solution
at the point s (s must lie in the interval (xold,x)).
iout Switch for calling solout :
iout=0 : no call,
iout=1 : solout only used for output,
iout=2 : dense output is performed in solout (in this case nrdens
must be greater than 0).
fileout A pointer to the stream used for messages, if you do not want any
message, just pass NULL.
icont An array containing the indexes of components for which dense
output is required. If no dense output is required, pass NULL.
licont The number of cells in icont.
Sophisticated setting of parameters
-----------------------------------
Several parameters have a default value (if set to 0) but, to better
adapt the code to your problem, you can specify particular initial
values.
uround The rounding unit, default 2.3E-16 (this default value can be
replaced in the code by DBL_EPSILON providing float.h defines it
in your system).
safe Safety factor in the step size prediction, default 0.9.
fac1 Parameters for step size selection; the new step size is chosen
fac2 subject to the restriction fac1 <= hnew/hold <= fac2.
Default values are fac1=0.333 and fac2=6.0.
beta The "beta" for stabilized step size control (see section IV.2 of our
book). Larger values for beta ( <= 0.1 ) make the step size control
more stable. Negative initial value provoke beta=0; default beta=0.
hmax Maximal step size, default xend-x.
h Initial step size, default is a guess computed by the function hinit.
nmax Maximal number of allowed steps, default 100000.
meth Switch for the choice of the method coefficients; at the moment the
only possibility and default value are 1.
nstiff Test for stiffness is activated when the current step number is a
multiple of nstiff. A negative value means no test and the default
is 1000.
nrdens Number of components for which dense outpout is required, default 0.
For 0 < nrdens < n, the components have to be specified in icont[0],
icont[1], ... icont[nrdens-1]. Note that if nrdens=0 or nrdens=n, no
icont is needed, pass NULL.
Memory requirements
-------------------
The function dop853 allocates dynamically 11*n doubles for the method
stages, 8*nrdens doubles for the interpolation if dense output is
performed and n unsigned if 0 < nrdens < n.
OUTPUT PARAMETERS
-----------------
y numerical solution at x=xRead() (see below).
dopri5 returns the following values
1 : computation successful,
2 : computation successful interrupted by solout,
-1 : input is not consistent,
-2 : larger nmax is needed,
-3 : step size becomes too small,
-4 : the problem is probably stff (interrupted).
Several functions provide access to different values :
xRead x value for which the solution has been computed (x=xend after
successful return).
hRead Predicted step size of the last accepted step (useful for a subsequent
call to dop853).
nstepRead Number of used steps.
naccptRead Number of accepted steps.
nrejctRead Number of rejected steps.
nfcnRead Number of function calls.
*/
/* DOPRI5
------
This code computes the numerical solution of a system of first order ordinary
differential equations y'=f(x,y). It uses an explicit Runge-Kutta method of
order (4)5 due to Dormand & Prince with step size control and dense output.
Authors : E. Hairer & G. Wanner
Universite de Geneve, dept. de Mathematiques
CH-1211 GENEVE 4, SWITZERLAND
E-mail : HAIRER@DIVSUN.UNIGE.CH, WANNER@DIVSUN.UNIGE.CH
The code is described in : E. Hairer, S.P. Norsett and G. Wanner, Solving
ordinary differential equations I, nonstiff problems, 2nd edition,
Springer Series in Computational Mathematics, Springer-Verlag (1993).
Version of April 28, 1994.
Remarks about the C version : this version allocates memory by itself, the
iwork array (among the initial FORTRAN parameters) has been splitted into
independant initial parameters, the statistical variables and last step size
and x have been encapsulated in the module and are now accessible through
dedicated functions, the variable names have been kept to maintain a kind
of reading compatibility between the C and FORTRAN codes; adaptation made by
J.Colinge (COLINGE@DIVSUN.UNIGE.CH).
INPUT PARAMETERS
----------------
n Dimension of the system (n < UINT_MAX).
fcn A pointer the the function definig the differential equation, this
function must have the following prototype
void fcn (unsigned n, double x, double *y, double *f)
where the array f will be filled with the function result.
x Initial x value.
*y Initial y values (double y[n]).
xend Final x value (xend-x may be positive or negative).
*rtoler Relative and absolute error tolerances. They can be both scalars or
*atoler vectors of length n (in the scalar case pass the addresses of
variables where you have placed the tolerance values).
itoler Switch for atoler and rtoler :
itoler=0 : both atoler and rtoler are scalars, the code keeps
roughly the local error of y[i] below
rtoler*abs(y[i])+atoler.
itoler=1 : both rtoler and atoler are vectors, the code keeps
the local error of y[i] below
rtoler[i]*abs(y[i])+atoler[i].
solout A pointer to the output function called during integration.
If iout >= 1, it is called after every successful step. If iout = 0,
pass a pointer equal to NULL. solout must must have the following
prototype
solout (long nr, double xold, double x, double* y, unsigned n, int* irtrn)
where y is the solution the at nr-th grid point x, xold is the
previous grid point and irtrn serves to interrupt the integration
(if set to a negative value).
Continuous output : during the calls to solout, a continuous solution
for the interval (xold,x) is available through the function
contd5(i,s)
which provides an approximation to the i-th component of the solution
at the point s (s must lie in the interval (xold,x)).
iout Switch for calling solout :
iout=0 : no call,
iout=1 : solout only used for output,
iout=2 : dense output is performed in solout (in this case nrdens
must be greater than 0).
fileout A pointer to the stream used for messages, if you do not want any
message, just pass NULL.
icont An array containing the indexes of components for which dense
output is required. If no dense output is required, pass NULL.
licont The number of cells in icont.
Sophisticated setting of parameters
-----------------------------------
Several parameters have a default value (if set to 0) but, to better
adapt the code to your problem, you can specify particular initial
values.
uround The rounding unit, default 2.3E-16 (this default value can be
replaced in the code by DBL_EPSILON providing float.h defines it
in your system).
safe Safety factor in the step size prediction, default 0.9.
fac1 Parameters for step size selection; the new step size is chosen
fac2 subject to the restriction fac1 <= hnew/hold <= fac2.
Default values are fac1=0.2 and fac2=10.0.
beta The "beta" for stabilized step size control (see section IV.2 of our
book). Larger values for beta ( <= 0.1 ) make the step size control
more stable. dopri5 needs a larger beta than Higham & Hall. Negative
initial value provoke beta=0; default beta=0.04.
hmax Maximal step size, default xend-x.
h Initial step size, default is a guess computed by the function hinit.
nmax Maximal number of allowed steps, default 100000.
meth Switch for the choice of the method coefficients; at the moment the
only possibility and default value are 1.
nstiff Test for stiffness is activated when the current step number is a
multiple of nstiff. A negative value means no test and the default
is 1000.
nrdens Number of components for which dense outpout is required, default 0.
For 0 < nrdens < n, the components have to be specified in icont[0],
icont[1], ... icont[nrdens-1]. Note that if nrdens=0 or nrdens=n, no
icont is needed, pass NULL.
Memory requirements
-------------------
The function dopri5 allocates dynamically 8*n doubles for the method
stages, 5*nrdens doubles for the interpolation if dense output is
performed and n unsigned if 0 < nrdens < n.
OUTPUT PARAMETERS
-----------------
y numerical solution at x=xRead() (see below).
dopri5 returns the following values
1 : computation successful,
2 : computation successful interrupted by solout,
-1 : input is not consistent,
-2 : larger nmax is needed,
-3 : step size becomes too small,
-4 : the problem is probably stff (interrupted).
Several functions provide access to different values :
xRead x value for which the solution has been computed (x=xend after
successful return).
hRead Predicted step size of the last accepted step (useful for a
subsequent call to dopri5).
nstepRead Number of used steps.
naccptRead Number of accepted steps.
nrejctRead Number of rejected steps.
nfcnRead Number of function calls.
*/
#include "my_rhs.h"
#include <stdio.h>
#include <limits.h>
typedef void (*FcnEqDiff)(unsigned n, double x, double *y, double *f);
typedef void (*SolTrait)(long nr, double xold, double x, double* y, unsigned n, int* irtrn);
extern int dop853
(unsigned n, /* dimension of the system <= UINT_MAX-1*/
FcnEqDiff fcn, /* function computing the value of f(x,y) */
double x, /* initial x-value */
double* y, /* initial values for y */
double xend, /* final x-value (xend-x may be positive or negative) */
double* rtoler, /* relative error tolerance */
double* atoler, /* absolute error tolerance */
int itoler, /* switch for rtoler and atoler */
SolTrait solout, /* function providing the numerical solution during integration */
int iout, /* switch for calling solout */
FILE* fileout, /* messages stream */
double uround, /* rounding unit */
double safe, /* safety factor */
double fac1, /* parameters for step size selection */
double fac2,
double beta, /* for stabilized step size control */
double hmax, /* maximal step size */
double h, /* initial step size */
long nmax, /* maximal number of allowed steps */
int meth, /* switch for the choice of the coefficients */
long nstiff, /* test for stiffness */
unsigned nrdens, /* number of components for which dense outpout is required */
unsigned* icont, /* indexes of components for which dense output is required, >= nrdens */
unsigned licont, /* declared length of icon */
double *work
);
extern double contd8
(unsigned ii, /* index of desired component */
double x /* approximation at x */
);
extern int dopri5
(unsigned n, /* dimension of the system <= UINT_MAX-1*/
FcnEqDiff fcn, /* function computing the value of f(x,y) */
double x, /* initial x-value */
double* y, /* initial values for y */
double xend, /* final x-value (xend-x may be positive or negative) */
double* rtoler, /* relative error tolerance */
double* atoler, /* absolute error tolerance */
int itoler, /* switch for rtoler and atoler */
SolTrait solout, /* function providing the numerical solution during integration */
int iout, /* switch for calling solout */
FILE* fileout, /* messages stream */
double uround, /* rounding unit */
double safe, /* safety factor */
double fac1, /* parameters for step size selection */
double fac2,
double beta, /* for stabilized step size control */
double hmax, /* maximal step size */
double h, /* initial step size */
long nmax, /* maximal number of allowed steps */
int meth, /* switch for the choice of the coefficients */
long nstiff, /* test for stiffness */
unsigned nrdens, /* number of components for which dense outpout is required */
unsigned* icont, /* indexes of components for which dense output is required, >= nrdens */
unsigned licont , /* declared length of icon */
double *work
);
extern double contd5
(unsigned ii, /* index of desired component */
double x /* approximation at x */
);
void dprhs(unsigned n, double t, double *y, double *f);
void dp_err(int k);
int dp(int *istart, double *y, double *t, int n, double tout, double *tol, double *atol, int flag, int *kflag);
int dormprin(int *istart, double *y, double *t, int n, double tout, double *tol, double *atol, int flag, int *kflag);
long nfcnRead(void);
long nstepRead(void);
long naccptRead(void);
long nrejctRead(void);
double hRead(void);
double xRead(void);
int dop853(unsigned n, FcnEqDiff fcn, double x, double *y, double xend, double *rtoler, double *atoler, int itoler, SolTrait solout, int iout, FILE *fileout, double uround, double safe, double fac1, double fac2, double beta, double hmax, double h, long nmax, int meth, long nstiff, unsigned nrdens, unsigned *icont, unsigned licont, double *work);
double contd8(unsigned ii, double x);
int dopri5(unsigned n, FcnEqDiff fcn, double x, double *y, double xend, double *rtoler, double *atoler, int itoler, SolTrait solout, int iout, FILE *fileout, double uround, double safe, double fac1, double fac2, double beta, double hmax, double h, long nmax, int meth, long nstiff, unsigned nrdens, unsigned *icont, unsigned licont, double *work);
double contd5(unsigned ii, double x);
extern long nfcnRead (void); /* encapsulation of statistical data */
extern long nstepRead (void);
extern long naccptRead (void);
extern long nrejctRead (void);
extern double hRead (void);
extern double xRead (void);
|