20void new_scale_factors_dp(
double *yscale,
double *y0,
double *dydx0,
21 double h_start,
double *tiny,
long *accmode,
double *accuracy,
23void initial_scale_factors_dp(
double *yscale,
double *y0,
double *dydx0,
24 double h_start,
double *tiny,
long *accmode,
double *accuracy,
25 double *accur,
double x0,
double xf,
long n_eq);
27static double stepIncreaseFactor = 0.50;
28static double stepDecreaseFactor = 0.95;
30void bs_qctune(
double newStepIncreaseFactor,
double newStepDecreaseFactor) {
31 if (newStepIncreaseFactor > 0)
32 stepIncreaseFactor = newStepIncreaseFactor;
33 if (newStepDecreaseFactor > 0)
34 stepDecreaseFactor = newStepDecreaseFactor;
54 double *stepRecommended,
57 void (*derivs)(
double *dydx,
double *y,
double x),
61 static double **solution, *hSqr, *yLast, *yError;
62 long i, j, iWorst = 0, code, nuse;
63 double maxError, error, yInterp;
64 static long mmidSteps[IMAX] = {2, 4, 6, 8, 12, 16, 24, 32, 48, 64, 96};
65 static long lastEquations = 0;
67 if (equations > lastEquations) {
68 if (lastEquations != 0) {
72 free_array_2d((
void **)solution,
sizeof(*solution), 0, lastEquations - 1, 0, NUSE - 1);
74 yLast =
tmalloc(
sizeof(
double) * equations);
75 yError =
tmalloc(
sizeof(
double) * equations);
76 hSqr =
tmalloc(
sizeof(
double) * IMAX);
77 solution = (
double **)
array_2d(
sizeof(
double), 0, equations - 1, 0, NUSE - 1);
78 lastEquations = equations;
83 printf(
"step = %e\n", step);
85 for (i = 0; i < IMAX; i++) {
86 mmid(yInitial, dydxInitial, equations, *x, step, mmidSteps[i], yFinal, derivs);
87 hSqr[i % NUSE] = sqr(step / mmidSteps[i]);
88 nuse = i > NUSE ? NUSE : i;
89 for (j = 0; j < equations; j++) {
91 solution[j][i % NUSE] = yFinal[j];
100 yError[j] = yInterp - yLast[j];
104 printf(
"i = %ld, j = %ld: yFinal = %.10e, yError = %.10e, yScale = %.10e\n",
105 i, j, yFinal[j], yInterp, yError[j], yScale[j]);
110 for (j = 0; j < equations; j++)
111 if (maxError < (error = FABS(yError[j] / yScale[j]))) {
116 printf(
"maxError = %e, iWorst = %ld\n", maxError, iWorst);
118 if (maxError < 1.0) {
120 *stepRecommended = *stepUsed = step;
123 *stepRecommended *= stepDecreaseFactor;
126 *stepRecommended *= stepIncreaseFactor / sqrt(maxError);
129 printf(
"returning with i=%ld, stepUsed=%e, stepRec=%e\n", i, *stepUsed, *stepRecommended);
131 for (j = 0; j < equations; j++)
132 yFinal[j] = yLast[j];
141 for (i = 0; i < (IMAX - NUSE) / 2; i++)
143 }
while ((*x + step) != *x);
144 fprintf(stderr,
"error: step size underflow in bs_step()--step reduced to %e\n", step);
181 void (*derivs)(
double *dydx,
double *y,
double x),
197 double (*exit_func)(
double *dydx,
double *y,
double x),
199 double exit_accuracy,
202 void (*store_data)(
double *dydx,
double *y,
double x,
double exf)
204 double *y_return, *accur;
205 double *dydx0, *y1, *dydx1, *dydx2, *y2;
206 double ex0, ex1, ex2, x1, x2, *yscale;
207 double h_used, h_next, xdiff;
208 long i, n_step_ups = 0, is_zero;
209#define MAX_N_STEP_UPS 10
212 return (DIFFEQ_XI_GT_XF);
213 if (FABS(*x0 - xf) < x_accuracy)
214 return (DIFFEQ_SOLVED_ALREADY);
227 for (i = 0; i < n_eq; i++) {
228 if (accmode[i] < 0 || accmode[i] > 3)
229 bomb(
"accmode must be on [0, 3] (bs_odeint)", NULL);
230 if (accmode[i] < 2 && tiny[i] < TINY)
236 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
237 y1 =
tmalloc(
sizeof(
double) * n_eq);
238 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
239 y2 =
tmalloc(
sizeof(
double) * n_eq);
240 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
241 yscale =
tmalloc(
sizeof(
double) * n_eq);
244 (*derivs)(dydx0, y0, *x0);
249 accur =
tmalloc(
sizeof(
double) * n_eq);
250 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
251 accuracy, accur, *x0, xf, n_eq);
253 ex0 = exit_func ? (*exit_func)(dydx0, y0, *x0) : 0;
255 (*store_data)(dydx0, y0, *x0, ex0);
260 if (exit_func && FABS(ex0) < exit_accuracy) {
262 if (n_to_skip == 0) {
264 (*store_data)(dydx0, y0, *x0, ex0);
265 for (i = 0; i < n_eq; i++)
279 return (DIFFEQ_ZERO_FOUND);
288 if ((xdiff = xf - *x0) < h_start)
292 if (!bs_step(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
293 yscale, n_eq, derivs, misses)) {
294 if (n_step_ups++ > MAX_N_STEP_UPS)
295 bomb(
"error: cannot take initial step (bs_odeint--1)", NULL);
296 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
300 (*derivs)(dydx1, y1, x1);
301 ex1 = exit_func ? (*exit_func)(dydx1, y1, x1) : 0;
303 (*store_data)(dydx1, y1, x1, ex1);
305 if (exit_func && SIGN(ex0) != SIGN(ex1) && !is_zero) {
314 if (FABS(xdiff = xf - x1) < x_accuracy) {
317 (*derivs)(dydx1, y1, x1);
318 ex1 = exit_func ? (*exit_func)(dydx1, y1, x1) : 0;
319 (*store_data)(dydx1, y1, x1, ex1);
321 for (i = 0; i < n_eq; i++)
336 return (DIFFEQ_END_OF_INTERVAL);
339 SWAP_PTR(dydx0, dydx1);
344 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
346 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
352 printf(
"failure in bs_odeint(): solution stepped outside interval\n");
364 return (DIFFEQ_OUTSIDE_INTERVAL);
367 if (FABS(ex1) < exit_accuracy) {
368 for (i = 0; i < n_eq; i++)
382 return (DIFFEQ_ZERO_FOUND);
388 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0 + TINY);
391 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
393 if (!bs_step(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
394 yscale, n_eq, derivs, misses))
395 bomb(
"step size too small (bs_odeint--2)", NULL);
397 (*derivs)(dydx2, y2, x2);
398 ex2 = (*exit_func)(dydx2, y2, x2);
399 if (FABS(ex2) < exit_accuracy) {
400 for (i = 0; i < n_eq; i++)
414 return (DIFFEQ_ZERO_FOUND);
417 if (SIGN(ex1) == SIGN(ex2)) {
419 SWAP_PTR(dydx1, dydx2);
424 SWAP_PTR(dydx0, dydx2);
475 double *dydx0, *y1, *dydx1, *dydx2, *y2;
476 double x1, *yscale, *accur;
477 double h_used, h_next, xdiff;
478 long i, n_step_ups = 0;
479#define MAX_N_STEP_UPS 10
482 return (DIFFEQ_XI_GT_XF);
483 if (fabs(*x0 - xf) < x_accuracy)
484 return (DIFFEQ_SOLVED_ALREADY);
497 for (i = 0; i < n_eq; i++) {
498 if (accmode[i] < 0 || accmode[i] > 3)
499 bomb(
"accmode must be on [0, 3] (bs_odeint)", NULL);
500 if (accmode[i] < 2 && tiny[i] < TINY)
506 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
507 y1 =
tmalloc(
sizeof(
double) * n_eq);
508 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
509 y2 =
tmalloc(
sizeof(
double) * n_eq);
510 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
511 yscale =
tmalloc(
sizeof(
double) * n_eq);
514 (*derivs)(dydx0, y0, *x0);
519 accur =
tmalloc(
sizeof(
double) * n_eq);
520 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
521 accuracy, accur, *x0, xf, n_eq);
525 if ((xdiff = xf - *x0) < h_start)
529 if (!bs_step(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
530 yscale, n_eq, derivs, misses)) {
531 if (n_step_ups++ > MAX_N_STEP_UPS)
532 bomb(
"error: cannot take initial step (bs_odeint1--1)", NULL);
533 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
537 if (fabs(xdiff = xf - x1) < x_accuracy) {
539 for (i = 0; i < n_eq; i++)
554 return (DIFFEQ_END_OF_INTERVAL);
557 (*derivs)(dydx1, y1, x1);
559 SWAP_PTR(dydx0, dydx1);
563 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
565 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
602 void (*derivs)(
double *dydx,
double *y,
double x),
620 double exit_accuracy,
623 double *y_return, *accur;
624 double *dydx0, *y1, *dydx1, *dydx2, *y2;
625 double ex0, ex1, ex2, x1, x2, *yscale;
626 double h_used, h_next, xdiff;
627 long i, n_step_ups = 0, is_zero;
628#define MAX_N_STEP_UPS 10
631 return (DIFFEQ_XI_GT_XF);
632 if (fabs(*x0 - xf) < x_accuracy)
633 return (DIFFEQ_SOLVED_ALREADY);
634 if (i_exit_value < 0 || i_exit_value >= n_eq)
635 bomb(
"index of variable for exit testing is out of range (bs_odeint2)", NULL);
648 for (i = 0; i < n_eq; i++) {
649 if (accmode[i] < 0 || accmode[i] > 3)
650 bomb(
"accmode must be on [0, 3] (bs_odeint2)", NULL);
651 if (accmode[i] < 2 && tiny[i] < TINY)
657 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
658 y1 =
tmalloc(
sizeof(
double) * n_eq);
659 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
660 y2 =
tmalloc(
sizeof(
double) * n_eq);
661 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
662 yscale =
tmalloc(
sizeof(
double) * n_eq);
665 (*derivs)(dydx0, y0, *x0);
670 accur =
tmalloc(
sizeof(
double) * n_eq);
671 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
672 accuracy, accur, *x0, xf, n_eq);
674 ex0 = exit_value - y0[i_exit_value];
678 if (fabs(ex0) < exit_accuracy) {
680 if (n_to_skip == 0) {
681 for (i = 0; i < n_eq; i++)
695 return (DIFFEQ_ZERO_FOUND);
704 if ((xdiff = xf - *x0) < h_start)
708 if (!bs_step(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
709 yscale, n_eq, derivs, misses)) {
710 if (n_step_ups++ > MAX_N_STEP_UPS) {
711 bomb(
"error: cannot take initial step (bs_odeint2--1)", NULL);
713 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
717 (*derivs)(dydx1, y1, x1);
718 ex1 = exit_value - y1[i_exit_value];
720 if (SIGN(ex0) != SIGN(ex1) && !is_zero) {
729 if (fabs(xdiff = xf - x1) < x_accuracy) {
731 for (i = 0; i < n_eq; i++)
746 return (DIFFEQ_END_OF_INTERVAL);
749 SWAP_PTR(dydx0, dydx1);
754 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
756 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
764 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0 + TINY);
767 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
769 if (!bs_step(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
770 yscale, n_eq, derivs, misses))
771 bomb(
"step size too small (bs_odeint2--2)", NULL);
773 (*derivs)(dydx2, y2, x2);
774 ex2 = exit_value - y2[i_exit_value];
775 if (fabs(ex2) < exit_accuracy) {
776 for (i = 0; i < n_eq; i++)
790 return (DIFFEQ_ZERO_FOUND);
793 if (SIGN(ex1) == SIGN(ex2)) {
795 SWAP_PTR(dydx1, dydx2);
800 SWAP_PTR(dydx0, dydx2);
838 void (*derivs)(
double *dydx,
double *y,
double x),
854 double (*exit_func)(
double *dydx,
double *y,
double x),
858 static double *yscale;
859 static double *dydx0, *y1, *dydx1, *dydx2, *y2, *accur, *y0;
860 static long last_neq = 0;
861 double ex0, ex1, ex2, x1, x2;
862 double h_used, h_next, xdiff;
863 long i, n_step_ups = 0;
864#define MAX_N_STEP_UPS 10
867 return (DIFFEQ_XI_GT_XF);
868 if (FABS(*x0 - xf) < x_accuracy)
869 return (DIFFEQ_SOLVED_ALREADY);
882 for (i = 0; i < n_eq; i++) {
883 if (accmode[i] < 0 || accmode[i] > 3)
884 bomb(
"accmode must be on [0, 3] (bs_odeint)", NULL);
885 if (accmode[i] < 2 && tiny[i] < TINY)
890 if (last_neq < n_eq) {
901 y0 =
tmalloc(
sizeof(
double) * n_eq);
902 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
903 y1 =
tmalloc(
sizeof(
double) * n_eq);
904 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
905 y2 =
tmalloc(
sizeof(
double) * n_eq);
906 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
907 yscale =
tmalloc(
sizeof(
double) * n_eq);
908 accur =
tmalloc(
sizeof(
double) * n_eq);
912 for (i = 0; i < n_eq; i++)
916 (*derivs)(dydx0, y0, *x0);
921 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
922 accuracy, accur, *x0, xf, n_eq);
924 ex0 = (*exit_func)(dydx0, y0, *x0);
928 if (FABS(ex0) < exit_accuracy) {
929 for (i = 0; i < n_eq; i++)
932 return (DIFFEQ_ZERO_FOUND);
936 if ((xdiff = xf - *x0) < h_start)
940 if (!bs_step(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
941 yscale, n_eq, derivs, misses)) {
942 if (n_step_ups++ > MAX_N_STEP_UPS)
943 bomb(
"error: cannot take initial step (bs_odeint3--1)", NULL);
944 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
948 (*derivs)(dydx1, y1, x1);
949 ex1 = (*exit_func)(dydx1, y1, x1);
950 if (SIGN(ex0) != SIGN(ex1))
953 if (FABS(xdiff = xf - x1) < x_accuracy) {
955 for (i = 0; i < n_eq; i++)
959 return (DIFFEQ_END_OF_INTERVAL);
962 SWAP_PTR(dydx0, dydx1);
967 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
969 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
975 printf(
"failure in bs_odeint3(): solution stepped outside interval\n");
976 return (DIFFEQ_OUTSIDE_INTERVAL);
979 if (FABS(ex1) < exit_accuracy) {
980 for (i = 0; i < n_eq; i++)
983 return (DIFFEQ_ZERO_FOUND);
989 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0 + TINY);
992 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
994 if (!bs_step(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
995 yscale, n_eq, derivs, misses))
996 bomb(
"step size too small (bs_odeint3--2)", NULL);
998 (*derivs)(dydx2, y2, x2);
999 ex2 = (*exit_func)(dydx2, y2, x2);
1000 if (FABS(ex2) < exit_accuracy) {
1001 for (i = 0; i < n_eq; i++)
1004 return (DIFFEQ_ZERO_FOUND);
1007 if (SIGN(ex1) == SIGN(ex2)) {
1009 SWAP_PTR(dydx1, dydx2);
1014 SWAP_PTR(dydx0, dydx2);
1054 void (*derivs)(
double *dydx,
double *y,
double x),
1072 double exit_accuracy,
1074 void (*store_data)(
double *dydx,
double *y,
double x,
double exf)
1076 double *y_return, *accur;
1077 double *dydx0, *y1, *dydx1, *dydx2, *y2;
1078 double ex0, ex1, ex2, x1, x2, *yscale;
1079 double h_used, h_next, xdiff;
1080 long i, n_step_ups = 0, is_zero;
1081#define MAX_N_STEP_UPS 10
1084 return (DIFFEQ_XI_GT_XF);
1085 if (fabs(*x0 - xf) < x_accuracy)
1086 return (DIFFEQ_SOLVED_ALREADY);
1087 if (i_exit_value < 0 || i_exit_value >= n_eq)
1088 bomb(
"index of variable for exit testing is out of range (bs_odeint4)", NULL);
1101 for (i = 0; i < n_eq; i++) {
1102 if (accmode[i] < 0 || accmode[i] > 3)
1103 bomb(
"accmode must be on [0, 3] (bs_odeint4)", NULL);
1104 if (accmode[i] < 2 && tiny[i] < TINY)
1110 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
1111 y1 =
tmalloc(
sizeof(
double) * n_eq);
1112 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
1113 y2 =
tmalloc(
sizeof(
double) * n_eq);
1114 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
1115 yscale =
tmalloc(
sizeof(
double) * n_eq);
1118 (*derivs)(dydx0, y0, *x0);
1123 accur =
tmalloc(
sizeof(
double) * n_eq);
1124 initial_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1125 accuracy, accur, *x0, xf, n_eq);
1127 ex0 = exit_value - y0[i_exit_value];
1129 (*store_data)(dydx0, y0, *x0, ex0);
1133 if (fabs(ex0) < exit_accuracy) {
1135 if (n_to_skip == 0) {
1137 (*store_data)(dydx0, y0, *x0, ex0);
1138 for (i = 0; i < n_eq; i++)
1139 y_return[i] = y0[i];
1152 return (DIFFEQ_ZERO_FOUND);
1161 if ((xdiff = xf - *x0) < h_start)
1165 if (!bs_step(y1, &x1, y0, dydx0, h_start, &h_used, &h_next,
1166 yscale, n_eq, derivs, misses)) {
1167 if (n_step_ups++ > MAX_N_STEP_UPS) {
1168 bomb(
"error: cannot take initial step (bs_odeint4--1)", NULL);
1170 h_start = (n_step_ups - 1 ? h_start * 10 : h_used * 10);
1174 (*derivs)(dydx1, y1, x1);
1175 ex1 = exit_value - y1[i_exit_value];
1177 (*store_data)(dydx1, y1, x1, ex1);
1179 if (SIGN(ex0) != SIGN(ex1) && !is_zero) {
1188 if (fabs(xdiff = xf - x1) < x_accuracy) {
1191 (*derivs)(dydx1, y1, x1);
1192 ex1 = exit_value - y0[i_exit_value];
1193 (*store_data)(dydx1, y1, x1, ex1);
1195 for (i = 0; i < n_eq; i++)
1196 y_return[i] = y1[i];
1210 return (DIFFEQ_END_OF_INTERVAL);
1213 SWAP_PTR(dydx0, dydx1);
1218 h_start = (h_next > h_max ? (h_max ? h_max : h_next) : h_next);
1220 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1228 h_start = -ex0 * (x1 - *x0) / (ex1 - ex0 + TINY);
1231 new_scale_factors_dp(yscale, y0, dydx0, h_start, tiny, accmode,
1233 if (!bs_step(y2, &x2, y0, dydx0, h_start, &h_used, &h_next,
1234 yscale, n_eq, derivs, misses))
1235 bomb(
"step size too small (bs_odeint4--2)", NULL);
1237 (*derivs)(dydx2, y2, x2);
1238 ex2 = exit_value - y2[i_exit_value];
1239 if (fabs(ex2) < exit_accuracy) {
1240 for (i = 0; i < n_eq; i++)
1241 y_return[i] = y2[i];
1254 return (DIFFEQ_ZERO_FOUND);
1257 if (SIGN(ex1) == SIGN(ex2)) {
1259 SWAP_PTR(dydx1, dydx2);
1264 SWAP_PTR(dydx0, dydx2);
void ** array_2d(uint64_t size, uint64_t lower1, uint64_t upper1, uint64_t lower2, uint64_t upper2)
Allocates a 2D array with specified lower and upper indices for both dimensions.
int tfree(void *ptr)
Frees a memory block and records the deallocation if tracking is enabled.
int free_array_2d(void **array, uint64_t size, uint64_t lower1, uint64_t upper1, uint64_t lower2, uint64_t upper2)
Frees a 2D array and its associated memory.
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 bs_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)
Integrates a system of ordinary differential equations using the Bulirsch-Stoer method with optimized...
long bs_odeint1(double *y0, void(*derivs)(), 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)
Integrates a system of ordinary differential equations to a specified upper limit without exit condit...
long bs_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))
Integrates a system of ordinary differential equations until a specified component reaches a target v...
long bs_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)
Integrates a system of ordinary differential equations until a specified component reaches a target v...
long bs_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 exf))
Integrates a system of ordinary differential equations using the Bulirsch-Stoer method.
double LagrangeInterp(double *x, double *f, long order1, double x0, long *returnCode)
Performs Lagrange interpolation of data.
void mmid(double *yInitial, double *dydxInitial, long equations, double xInitial, double interval, long steps, double *yFinal, void(*derivs)(double *_dydx, double *_y, double _x))
Integrates a system of ODEs using the modified midpoint method.