23#define MAX_EXIT_ITERATIONS 400
24#define ITER_FACTOR 0.995
28void new_scale_factors_dp(
39 for (i = 0; i < n_eq; i++) {
46 (y0[i] + dydx0[i] * h_start + tiny[i]) * accuracy[i];
50 (dydx0[i] * h_start + tiny[i]) * accuracy[i];
53 yscale[i] = accuracy[i];
56 yscale[i] = accuracy[i] * h_start;
59 printf(
"error: accmode[%d] = %ld (new_scale_factors_dp)\n", i, accmode[i]);
64 printf(
"error: yscale[%d] = %e\n", i, yscale[i]);
70void initial_scale_factors_dp(
83 for (i = 0; i < n_eq; i++) {
84 if ((accur[i] = accuracy[i]) <= 0) {
85 printf(
"error: accuracy[%d] = %e (initial_scale_factors_dp)\n", i, accuracy[i]);
93 yscale[i] = (y0[i] + dydx0[i] * h_start + tiny[i]) * accur[i];
96 yscale[i] = (dydx0[i] * h_start + tiny[i]) * accur[i];
102 yscale[i] = (accur[i] /= (xf - x0)) * h_start;
105 printf(
"error: accmode[%d] = %ld (initial_scale_factors_dp)\n", i, accmode[i]);
109 if (yscale[i] <= 0) {
110 printf(
"error: yscale[%d] = %e (initial_scale_factors_dp)\n", i, yscale[i]);
116void report_state_dp(FILE *fp,
double *y,
double *dydx,
double *yscale,
long *misses,
117 double x,
double h,
long n_eq) {
119 fputs(
"integration state:\n", fp);
120 fprintf(fp,
"%ld equations, indep.var.=%e, step size=%e",
122 fprintf(fp,
"\ny : ");
123 for (i = 0; i < n_eq; i++)
124 fprintf(fp,
"%e ", y[i]);
125 fprintf(fp,
"\ndydx : ");
126 for (i = 0; i < n_eq; i++)
127 fprintf(fp,
"%e ", dydx[i]);
128 fprintf(fp,
"\ntol.scale: ");
129 for (i = 0; i < n_eq; i++)
130 fprintf(fp,
"%e ", yscale[i]);
131 fprintf(fp,
"\nmisses : ");
132 for (i = 0; i < n_eq; i++)
133 fprintf(fp,
"%ld ", misses[i]);
148 void (*derivs)(
double *dydx,
double *y,
double x)
151 static long last_n_eq = 0;
152 static double *k1, *k2, *k3, *yTemp, *dydxTemp;
157 if (last_n_eq < n_eq) {
158 if (last_n_eq != 0) {
166 k1 =
tmalloc(
sizeof(*k1) * n_eq);
167 k2 =
tmalloc(
sizeof(*k2) * n_eq);
168 k3 =
tmalloc(
sizeof(*k3) * n_eq);
169 yTemp =
tmalloc(
sizeof(*yTemp) * n_eq);
170 dydxTemp =
tmalloc(
sizeof(*dydxTemp) * n_eq);
181 for (i = 0; i < n_eq; i++) {
183 yTemp[i] = yi[i] + k1[i] / 2;
188 (*derivs)(dydxTemp, yTemp, x1);
189 for (i = 0; i < n_eq; i++) {
190 k2[i] = h * dydxTemp[i];
191 yTemp[i] = yi[i] + k2[i] / 2;
195 (*derivs)(dydxTemp, yTemp, x1);
196 for (i = 0; i < n_eq; i++) {
197 k3[i] = h * dydxTemp[i];
198 yTemp[i] = yi[i] + k3[i];
203 (*derivs)(dydxTemp, yTemp, x1);
204 for (i = 0; i < n_eq; i++)
205 yf[i] = yi[i] + (k1[i] / 2 + k2[i] + k3[i] + h * dydxTemp[i] / 2) / 3;
215static double safetyMargin = 0.9;
216static double increasePower = 0.2;
217static double decreasePower = 0.25;
218static double maxIncreaseFactor = 4.0;
220void rk4_qctune(
double newSafetyMargin,
double newIncreasePower,
221 double newDecreasePower,
double newMaxIncreaseFactor) {
222 if (newSafetyMargin > 0 && newSafetyMargin < 1)
223 safetyMargin = newSafetyMargin;
224 if (newIncreasePower > 0)
225 increasePower = newIncreasePower;
226 if (newDecreasePower > 0)
227 decreasePower = newDecreasePower;
228 if (newMaxIncreaseFactor > 1)
229 maxIncreaseFactor = newMaxIncreaseFactor;
239 double *hRecommended,
242 void (*derivs)(
double *dydx,
double *y,
double x),
248 static long last_equations = 0;
249 static double *dydxTemp, *yTemp;
250 double hOver2, h, xTemp;
251 double error, maxError, hFactor;
252 long i, iWorst = 0, minStepped = 0, noAdaptation = 0;
255 if (last_equations < equations) {
256 if (last_equations != 0) {
260 last_equations = equations;
261 dydxTemp =
tmalloc(
sizeof(*dydxTemp) * equations);
262 yTemp =
tmalloc(
sizeof(*yTemp) * equations);
265 for (i = 0; i < equations; i++)
266 if (yScale[i] != DBL_MAX)
274 printf(
"x = %e, h = %e\n", *x, h);
280 puts(
"warning: step-size underflow in rk_qcstep()");
281 report_state_dp(stdout, yInitial, dydxInitial, yScale, misses, *x, h, equations);
282 if ((h = 2 * fabs(*x) * DBL_EPSILON) == 0)
294 rk4_step(yTemp, xTemp, yInitial, dydxInitial, hOver2, equations, derivs);
297 (*derivs)(dydxTemp, yTemp, xTemp);
299 rk4_step(yFinal, xTemp, yTemp, dydxTemp, hOver2, equations, derivs);
302 rk4_step(yTemp, xTemp, yInitial, dydxInitial, h, equations, derivs);
305 for (i = 0; i < equations; i++)
306 yTemp[i] = yFinal[i] - yTemp[i];
310 for (i = 0; i < equations; i++) {
312 printf(
"%ld: yTemp=%e yFinal=%e yScale=%e, ", i, yTemp[i], yFinal[i], yScale[i]);
314 if (maxError < (error = fabs(yTemp[i]) / yScale[i])) {
319 printf(
" error = %e\n", error);
334 hFactor = safetyMargin * pow(maxError, -increasePower);
337 hFactor = maxIncreaseFactor;
339 printf(
"maxError = %e, hFactor = %e\n", maxError, hFactor);
341 if (hFactor > maxIncreaseFactor)
342 hFactor = maxIncreaseFactor;
343 else if (hFactor < 1)
348 *hRecommended = hFactor * h;
353 for (i = 0; i < equations; i++)
354 yFinal[i] += yTemp[i] / 15;
367 hOver2 = (h = safetyMargin * h * pow(maxError, -decreasePower)) / 2;
400 void (*derivs)(
double *dydx,
double *y,
double x),
416 double (*exit_func)(
double *dydx,
double *y,
double x),
417 double exit_accuracy,
421 void (*store_data)(
double *dydx,
double *y,
double x,
double exval)) {
423 double *dydx0, *y1, *dydx1, *dydx2, *y2;
424 double ex0, ex1, ex2, *yscale, *accur;
425 double h_used, h_next, x1, x2, xdiff;
426 long i, n_exit_iterations, n_step_ups = 0, is_zero;
427#define MAX_N_STEP_UPS 10
430 return (DIFFEQ_XI_GT_XF);
431 if (FABS(*x0 - xf) < x_accuracy)
432 return (DIFFEQ_SOLVED_ALREADY);
446 for (i = 0; i < n_eq; i++) {
447 if (accmode[i] < -1 || accmode[i] > 3)
448 bomb(
"accmode must be on [-1, 3] (rk_odeint)", NULL);
449 if (accmode[i] < 2 && tiny[i] < TINY)
455 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
456 y1 =
tmalloc(
sizeof(
double) * n_eq);
457 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
458 y2 =
tmalloc(
sizeof(
double) * n_eq);
459 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
460 yscale =
tmalloc(
sizeof(
double) * n_eq);
463 (*derivs)(dydx0, y0, *x0);
468 accur =
tmalloc(
sizeof(
double) * n_eq);
469 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
470 accuracy, accur, *x0, xf, n_eq);
472 ex0 = exit_func ? (*exit_func)(dydx0, y0, *x0) : 0;
474 (*store_data)(dydx0, y0, *x0, ex0);
479 if (exit_func && FABS(ex0) < exit_accuracy) {
481 if (n_to_skip == 0) {
483 (*store_data)(dydx0, y0, *x0, ex0);
484 for (i = 0; i < n_eq; i++)
498 return (DIFFEQ_ZERO_FOUND);
507 if ((xdiff = xf - *x0) < h_start)
511 if (!rk_qcstep(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
512 yscale, n_eq, derivs, misses)) {
513 if (n_step_ups++ > MAX_N_STEP_UPS)
514 bomb(
"cannot take initial step (rk_odeint--1)", NULL);
515 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
519 (*derivs)(dydx1, y1, x1);
520 ex1 = exit_func ? (*exit_func)(dydx1, y1, x1) : 0;
522 (*store_data)(dydx1, y1, x1, ex1);
524 if (exit_func && SIGN(ex0) != SIGN(ex1) && !is_zero) {
534 if (FABS(xdiff = xf - x1) < x_accuracy) {
537 (*derivs)(dydx1, y1, x1);
538 ex1 = exit_func ? (*exit_func)(dydx1, y1, x1) : 0;
539 (*store_data)(dydx1, y1, x1, ex1);
541 for (i = 0; i < n_eq; i++)
556 return (DIFFEQ_END_OF_INTERVAL);
560 SWAP_PTR(dydx0, dydx1);
565 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
567 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
574 printf(
"failure in rk_odeint(): solution stepped outside interval\n");
586 return (DIFFEQ_OUTSIDE_INTERVAL);
590 n_exit_iterations = MAX_EXIT_ITERATIONS;
593 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
596 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
598 if (!rk_qcstep(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
599 yscale, n_eq, derivs, misses))
600 bomb(
"step size too small (rk_odeint--2)", NULL);
602 (*derivs)(dydx2, y2, x2);
603 ex2 = (*exit_func)(dydx2, y2, x2);
604 if (FABS(ex2) < exit_accuracy) {
605 for (i = 0; i < n_eq; i++)
619 return (DIFFEQ_ZERO_FOUND);
622 if (SIGN(ex1) == SIGN(ex2)) {
624 SWAP_PTR(dydx1, dydx2);
629 SWAP_PTR(dydx0, dydx2);
633 }
while (n_exit_iterations--);
634 return (DIFFEQ_EXIT_COND_FAILED);
662 void (*derivs)(
double *dydx,
double *y,
double x),
679 double *dydx0, *y1, *dydx1, *dydx2, *y2;
680 double *yscale, *accur;
682 double h_used, h_next, xdiff;
683 long i, n_step_ups = 0;
684#define MAX_N_STEP_UPS 10
687 return (DIFFEQ_XI_GT_XF);
688 if (FABS(*x0 - xf) < x_accuracy)
689 return (DIFFEQ_SOLVED_ALREADY);
703 for (i = 0; i < n_eq; i++) {
704 if (accmode[i] < -1 || accmode[i] > 3)
705 bomb(
"accmode must be on [-1, 3] (rk_odeint)", NULL);
706 if (accmode[i] < 2 && tiny[i] < TINY)
712 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
713 y1 =
tmalloc(
sizeof(
double) * n_eq);
714 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
715 y2 =
tmalloc(
sizeof(
double) * n_eq);
716 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
717 yscale =
tmalloc(
sizeof(
double) * n_eq);
720 (*derivs)(dydx0, y0, *x0);
725 accur =
tmalloc(
sizeof(
double) * n_eq);
726 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
727 accuracy, accur, *x0, xf, n_eq);
731 if ((xdiff = xf - *x0) < h_start)
735 if (!rk_qcstep(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
736 yscale, n_eq, derivs, misses)) {
737 if (n_step_ups++ > MAX_N_STEP_UPS) {
738 puts(
"error: cannot take step (rk_odeint1--1)");
739 printf(
"xf = %.16e x0 = %.16e\n", xf, *x0);
740 printf(
"h_start = %.16e h_used = %.16e\n", h_start, h_used);
741 puts(
"dump of integration state:");
742 puts(
" variable value derivative scale misses");
743 puts(
"---------------------------------------------------------------");
744 for (i = 0; i < n_eq; i++)
745 printf(
" %5ld %13.6e %13.6e %13.6e %5ld \n",
746 i, y0[i], dydx0[i], yscale[i], misses[i]);
749 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
753 if (FABS(xdiff = xf - x1) < x_accuracy) {
755 for (i = 0; i < n_eq; i++)
770 return (DIFFEQ_END_OF_INTERVAL);
773 (*derivs)(dydx1, y1, x1);
775 SWAP_PTR(dydx0, dydx1);
779 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
781 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
815 void (*derivs)(
double *dydx,
double *y,
double x),
833 double exit_accuracy,
836 double *y_return, *accur;
837 double *dydx0, *y1, *dydx1, *dydx2, *y2;
838 double ex0, ex1, ex2, x1, x2, *yscale;
839 double h_used, h_next, xdiff;
840 long i, n_exit_iterations, n_step_ups = 0, is_zero;
841#define MAX_N_STEP_UPS 10
844 return (DIFFEQ_XI_GT_XF);
845 if (FABS(*x0 - xf) < x_accuracy)
846 return (DIFFEQ_SOLVED_ALREADY);
847 if (i_exit_value < 0 || i_exit_value >= n_eq)
848 bomb(
"index of variable for exit testing is out of range (rk_odeint2)", NULL);
862 for (i = 0; i < n_eq; i++) {
863 if (accmode[i] < -1 || accmode[i] > 3)
864 bomb(
"accmode must be on [-1, 3] (rk_odeint2)", NULL);
865 if (accmode[i] < 2 && tiny[i] < TINY)
871 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
872 y1 =
tmalloc(
sizeof(
double) * n_eq);
873 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
874 y2 =
tmalloc(
sizeof(
double) * n_eq);
875 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
876 yscale =
tmalloc(
sizeof(
double) * n_eq);
879 (*derivs)(dydx0, y0, *x0);
884 accur =
tmalloc(
sizeof(
double) * n_eq);
885 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
886 accuracy, accur, *x0, xf, n_eq);
888 ex0 = exit_value - y0[i_exit_value];
892 if (FABS(ex0) < exit_accuracy) {
894 if (n_to_skip == 0) {
895 for (i = 0; i < n_eq; i++)
909 return (DIFFEQ_ZERO_FOUND);
918 if ((xdiff = xf - *x0) < h_start)
922 if (!rk_qcstep(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
923 yscale, n_eq, derivs, misses)) {
924 if (n_step_ups++ > MAX_N_STEP_UPS) {
925 bomb(
"error: cannot take initial step (rk_odeint2--1)", NULL);
927 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
931 (*derivs)(dydx1, y1, x1);
932 ex1 = exit_value - y1[i_exit_value];
934 if (SIGN(ex0) != SIGN(ex1) && !is_zero) {
943 if (FABS(xdiff = xf - x1) < x_accuracy) {
945 for (i = 0; i < n_eq; i++)
960 return (DIFFEQ_END_OF_INTERVAL);
963 SWAP_PTR(dydx0, dydx1);
968 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
970 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
977 n_exit_iterations = MAX_EXIT_ITERATIONS;
980 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
983 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
985 if (!rk_qcstep(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
986 yscale, n_eq, derivs, misses))
987 bomb(
"step size too small (rk_odeint2--2)", NULL);
989 (*derivs)(dydx2, y2, x2);
990 ex2 = exit_value - y2[i_exit_value];
991 if (FABS(ex2) < exit_accuracy) {
992 for (i = 0; i < n_eq; i++)
1006 return (DIFFEQ_ZERO_FOUND);
1009 if (SIGN(ex1) == SIGN(ex2)) {
1011 SWAP_PTR(dydx1, dydx2);
1016 SWAP_PTR(dydx0, dydx2);
1020 }
while (n_exit_iterations--);
1021 return (DIFFEQ_EXIT_COND_FAILED);
1049 void (*derivs)(
double *dydx,
double *y,
double x),
1064 double *dydx2, *dydx1, *ytemp, *dydx0;
1065 double xh, hh, h6, x;
1069 return (DIFFEQ_XI_GT_XF);
1071 return (DIFFEQ_ZERO_STEPSIZE);
1072 if (!(n = (xf - x) / h + 0.5))
1076 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
1077 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
1078 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
1079 ytemp =
tmalloc(
sizeof(
double) * n_eq);
1085 for (j = 0; j < n; j++) {
1094 (*derivs)(dydx0, y0, x);
1095 for (i = n_eq; i >= 0; i--)
1096 ytemp[i] = y0[i] + hh * dydx0[i];
1099 (*derivs)(dydx1, ytemp, xh);
1100 for (i = n_eq; i >= 0; i--)
1101 ytemp[i] = y0[i] + hh * dydx1[i];
1104 (*derivs)(dydx2, ytemp, xh);
1105 for (i = n_eq; i >= 0; i--) {
1106 ytemp[i] = y0[i] + h * dydx2[i];
1107 dydx2[i] += dydx1[i];
1111 (*derivs)(dydx1, ytemp, x = xh + hh);
1112 for (i = n_eq; i >= 0; i--)
1113 y0[i] += h6 * (dydx0[i] + dydx1[i] + 2 * dydx2[i]);
1121 return (DIFFEQ_END_OF_INTERVAL);
1149 void (*derivs)(
double *dydx,
double *y,
double x),
1165 double (*exit_func)(
double *dydx,
double *y,
double x),
1167 double exit_accuracy
1169 static double *y0, *yscale;
1170 static double *dydx0, *y1, *dydx1, *dydx2, *y2, *accur;
1171 static long last_neq = 0;
1172 double ex0, ex1, ex2, x1, x2;
1173 double h_used, h_next, xdiff;
1174 long i, n_exit_iterations, n_step_ups = 0;
1175#define MAX_N_STEP_UPS 10
1178 return (DIFFEQ_XI_GT_XF);
1179 if (FABS(*x0 - xf) < x_accuracy)
1180 return (DIFFEQ_SOLVED_ALREADY);
1194 for (i = 0; i < n_eq; i++) {
1195 if (accmode[i] < -1 || accmode[i] > 3)
1196 bomb(
"accmode must be on [-1, 3] (rk_odeint)", NULL);
1197 if (accmode[i] < 2 && tiny[i] < TINY)
1202 if (last_neq < n_eq) {
1203 if (last_neq != 0) {
1213 y0 =
tmalloc(
sizeof(
double) * n_eq);
1214 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
1215 y1 =
tmalloc(
sizeof(
double) * n_eq);
1216 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
1217 y2 =
tmalloc(
sizeof(
double) * n_eq);
1218 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
1219 yscale =
tmalloc(
sizeof(
double) * n_eq);
1220 accur =
tmalloc(
sizeof(
double) * n_eq);
1224 for (i = 0; i < n_eq; i++)
1228 (*derivs)(dydx0, y0, *x0);
1233 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1234 accuracy, accur, *x0, xf, n_eq);
1236 ex0 = (*exit_func)(dydx0, y0, *x0);
1240 if (FABS(ex0) < exit_accuracy) {
1241 for (i = 0; i < n_eq; i++)
1244 return (DIFFEQ_ZERO_FOUND);
1248 if ((xdiff = xf - *x0) < h_start)
1252 if (!rk_qcstep(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
1253 yscale, n_eq, derivs, misses)) {
1254 if (n_step_ups++ > MAX_N_STEP_UPS)
1255 bomb(
"error: cannot take initial step (rk_odeint3--1)", NULL);
1256 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
1260 (*derivs)(dydx1, y1, x1);
1261 ex1 = (*exit_func)(dydx1, y1, x1);
1262 if (SIGN(ex0) != SIGN(ex1))
1265 if (FABS(xdiff = xf - x1) < x_accuracy) {
1267 for (i = 0; i < n_eq; i++)
1271 return (DIFFEQ_END_OF_INTERVAL);
1274 SWAP_PTR(dydx0, dydx1);
1279 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
1281 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1287 printf(
"failure in rk_odeint3(): solution stepped outside interval\n");
1288 return (DIFFEQ_OUTSIDE_INTERVAL);
1291 if (FABS(ex1) < exit_accuracy) {
1292 for (i = 0; i < n_eq; i++)
1295 return (DIFFEQ_ZERO_FOUND);
1299 n_exit_iterations = MAX_EXIT_ITERATIONS;
1302 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
1305 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1307 if (!rk_qcstep(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
1308 yscale, n_eq, derivs, misses))
1309 bomb(
"step size too small (rk_odeint3--2)", NULL);
1311 (*derivs)(dydx2, y2, x2);
1312 ex2 = (*exit_func)(dydx2, y2, x2);
1313 if (FABS(ex2) < exit_accuracy) {
1314 for (i = 0; i < n_eq; i++)
1317 return (DIFFEQ_ZERO_FOUND);
1320 if (SIGN(ex1) == SIGN(ex2)) {
1322 SWAP_PTR(dydx1, dydx2);
1327 SWAP_PTR(dydx0, dydx2);
1331 }
while (n_exit_iterations--);
1332 return (DIFFEQ_EXIT_COND_FAILED);
1365 void (*derivs)(
double *dydx,
double *y,
double x),
1383 double exit_accuracy,
1385 void (*store_data)(
double *dydx,
double *y,
double x,
double exf)
1387 double *y_return, *accur;
1388 double *dydx0, *y1, *dydx1, *dydx2, *y2;
1389 double ex0, ex1, ex2, x1, x2, *yscale;
1390 double h_used, h_next, xdiff;
1391 long i, n_exit_iterations, n_step_ups = 0, is_zero;
1392#define MAX_N_STEP_UPS 10
1395 return (DIFFEQ_XI_GT_XF);
1396 if (fabs(*x0 - xf) < x_accuracy)
1397 return (DIFFEQ_SOLVED_ALREADY);
1398 if (i_exit_value < 0 || i_exit_value >= n_eq)
1399 bomb(
"index of variable for exit testing is out of range (rk_odeint4)", NULL);
1413 for (i = 0; i < n_eq; i++) {
1414 if (accmode[i] < -1 || accmode[i] > 3)
1415 bomb(
"accmode must be on [-1, 3] (rk_odeint4)", NULL);
1416 if (accmode[i] < 2 && tiny[i] < TINY)
1422 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
1423 y1 =
tmalloc(
sizeof(
double) * n_eq);
1424 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
1425 y2 =
tmalloc(
sizeof(
double) * n_eq);
1426 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
1427 yscale =
tmalloc(
sizeof(
double) * n_eq);
1430 (*derivs)(dydx0, y0, *x0);
1435 accur =
tmalloc(
sizeof(
double) * n_eq);
1436 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1437 accuracy, accur, *x0, xf, n_eq);
1439 ex0 = exit_value - y0[i_exit_value];
1441 (*store_data)(dydx0, y0, *x0, ex0);
1445 if (fabs(ex0) < exit_accuracy) {
1447 if (n_to_skip == 0) {
1449 (*store_data)(dydx0, y0, *x0, ex0);
1450 for (i = 0; i < n_eq; i++)
1451 y_return[i] = y0[i];
1464 return (DIFFEQ_ZERO_FOUND);
1473 if ((xdiff = xf - *x0) < h_start)
1477 if (!rk_qcstep(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
1478 yscale, n_eq, derivs, misses)) {
1479 if (n_step_ups++ > MAX_N_STEP_UPS) {
1480 bomb(
"error: cannot take initial step (rk_odeint4--1)", NULL);
1482 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
1486 (*derivs)(dydx1, y1, x1);
1487 ex1 = exit_value - y1[i_exit_value];
1489 (*store_data)(dydx1, y1, x1, ex1);
1491 if (SIGN(ex0) != SIGN(ex1) && !is_zero) {
1500 if (fabs(xdiff = xf - x1) < x_accuracy) {
1503 (*derivs)(dydx1, y1, x1);
1504 ex1 = exit_value - y0[i_exit_value];
1505 (*store_data)(dydx1, y1, x1, ex1);
1507 for (i = 0; i < n_eq; i++)
1508 y_return[i] = y1[i];
1522 return (DIFFEQ_END_OF_INTERVAL);
1525 SWAP_PTR(dydx0, dydx1);
1530 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
1532 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1538 n_exit_iterations = MAX_EXIT_ITERATIONS;
1541 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
1544 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1546 if (!rk_qcstep(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
1547 yscale, n_eq, derivs, misses))
1548 bomb(
"step size too small (rk_odeint4--2)", NULL);
1550 (*derivs)(dydx2, y2, x2);
1551 ex2 = exit_value - y2[i_exit_value];
1552 if (fabs(ex2) < exit_accuracy) {
1553 for (i = 0; i < n_eq; i++)
1554 y_return[i] = y2[i];
1567 return (DIFFEQ_ZERO_FOUND);
1570 if (SIGN(ex1) == SIGN(ex2)) {
1572 SWAP_PTR(dydx1, dydx2);
1577 SWAP_PTR(dydx0, dydx2);
1581 }
while (n_exit_iterations--);
1582 return (DIFFEQ_EXIT_COND_FAILED);
1611 void (*derivs)(
double *dydx,
double *y,
double x),
1626 double (*exit_func)(
double *dydx,
double *y,
double x),
1628 double exit_accuracy,
1630 void (*stochastic)(
double *y,
double x,
double h)) {
1631 static double *y0, *yscale;
1632 static double *dydx0, *y1, *dydx1, *dydx2, *y2, *accur;
1633 static long last_neq = 0;
1634 double ex0, ex1, ex2, x1, x2;
1636 long i, n_exit_iterations;
1637#define MAX_N_STEP_UPS 10
1640 return (DIFFEQ_XI_GT_XF);
1641 if (FABS(*x0 - xf) < x_accuracy)
1642 return (DIFFEQ_SOLVED_ALREADY);
1644 if (last_neq < n_eq) {
1645 if (last_neq != 0) {
1655 y0 =
tmalloc(
sizeof(
double) * n_eq);
1656 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
1657 y1 =
tmalloc(
sizeof(
double) * n_eq);
1658 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
1659 y2 =
tmalloc(
sizeof(
double) * n_eq);
1660 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
1664 for (i = 0; i < n_eq; i++)
1668 (*derivs)(dydx0, y0, *x0);
1670 ex0 = (*exit_func)(dydx0, y0, *x0);
1674 if (FABS(ex0) < exit_accuracy) {
1675 for (i = 0; i < n_eq; i++)
1677 return (DIFFEQ_ZERO_FOUND);
1681 if ((xdiff = xf - *x0) < h_step)
1685 rk4_step(y1, x1, y0, dydx0, h_step, n_eq, derivs);
1688 (*stochastic)(y1, x1, h_step);
1691 (*derivs)(dydx1, y1, x1);
1692 ex1 = (*exit_func)(dydx1, y1, x1);
1693 if (SIGN(ex0) != SIGN(ex1))
1696 if (FABS(xdiff = xf - x1) < x_accuracy) {
1698 for (i = 0; i < n_eq; i++)
1701 return (DIFFEQ_END_OF_INTERVAL);
1704 SWAP_PTR(dydx0, dydx1);
1711 printf(
"failure in rk_odeint3_na(): solution stepped outside interval\n");
1712 return (DIFFEQ_OUTSIDE_INTERVAL);
1715 if (FABS(ex1) < exit_accuracy) {
1716 for (i = 0; i < n_eq; i++)
1719 return (DIFFEQ_ZERO_FOUND);
1723 n_exit_iterations = MAX_EXIT_ITERATIONS;
1727 h_step = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
1729 rk4_step(y2, x2, y0, dydx0, h_step, n_eq, derivs);
1732 (*derivs)(dydx2, y2, x2);
1733 ex2 = (*exit_func)(dydx2, y2, x2);
1734 if (FABS(ex2) < exit_accuracy) {
1735 for (i = 0; i < n_eq; i++)
1738 return (DIFFEQ_ZERO_FOUND);
1741 if (SIGN(ex1) == SIGN(ex2)) {
1743 SWAP_PTR(dydx1, dydx2);
1748 SWAP_PTR(dydx0, dydx2);
1752 }
while (n_exit_iterations--);
1753 return (DIFFEQ_EXIT_COND_FAILED);
int tfree(void *ptr)
Frees a memory block and records the deallocation if tracking is enabled.
void * tmalloc(uint64_t size_of_block)
Allocates a memory block of the specified size with zero initialization.
void bomb(char *error, char *usage)
Reports error messages to the terminal and aborts the program.
long rk_odeint1(double *y0, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_start, double h_max, double *h_rec)
Integrate ODEs without exit conditions or intermediate output.
long rk_odeint(double *y0, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_start, double h_max, double *h_rec, double(*exit_func)(double *dydx, double *y, double x), double exit_accuracy, long n_to_skip, void(*store_data)(double *dydx, double *y, double x, double exval))
Integrate a set of ODEs until the upper limit or an exit condition is met.
long rk_odeint_na(double *y0, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h, double h_max, double *h_rec)
Integrate ODEs without adaptive step-size control.
long rk_odeint2(double *y0, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_start, double h_max, double *h_rec, double exit_value, long i_exit_value, double exit_accuracy, long n_to_skip)
Integrate ODEs until a specific component reaches a target value.
long rk_odeint4(double *y0, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_start, double h_max, double *h_rec, double exit_value, long i_exit_value, double exit_accuracy, long n_to_skip, void(*store_data)(double *dydx, double *y, double x, double exf))
Integrate ODEs until a specific component reaches a target value with intermediate storage.
long rk_odeint3(double *yif, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_start, double h_max, double *h_rec, double(*exit_func)(double *dydx, double *y, double x), double exit_accuracy)
Integrate ODEs with an exit condition.
long rk_odeint3_na(double *yif, void(*derivs)(double *dydx, double *y, double x), long n_eq, double *accuracy, long *accmode, double *tiny, long *misses, double *x0, double xf, double x_accuracy, double h_step, double h_max, double *h_rec, double(*exit_func)(double *dydx, double *y, double x), double exit_accuracy, void(*stochastic)(double *y, double x, double h))
Integrate ODEs without adaptive step-size and with stochastic processes.