51 void (*derivs)(
double *_dydx,
double *_y,
double _x)) {
53 double x = 0, ynSave, h, hTimes2;
54 static double *ym, *yn;
55 static long last_equations = 0;
58 if (equations > last_equations) {
64 ym =
tmalloc(
sizeof(*ym) * equations);
65 yn =
tmalloc(
sizeof(*yn) * equations);
66 last_equations = equations;
69 hTimes2 = (h = interval / steps) * 2;
72 for (i = 0; i < equations; i++) {
74 yn[i] = yInitial[i] + h * dydxInitial[i];
78 for (j = 1; j < steps; j++) {
80 (*derivs)(dydxTemp, yn, x);
81 for (i = 0; i < equations; i++) {
83 yn[i] = ym[i] + hTimes2 * dydxTemp[i];
89 (*derivs)(dydxTemp, yn, x + interval);
90 for (i = 0; i < equations; i++)
91 yFinal[i] = (ym[i] + yn[i] + h * dydxTemp[i]) / 2;
120 void (*derivs)(
double *_dydx,
double *_y,
double _x)) {
121 static double *yFinal2;
122 static long i, last_equations = 0;
129 if (equations > last_equations) {
130 if (last_equations) {
134 yFinal2 =
tmalloc(
sizeof(*yFinal2) * equations);
135 last_equations = equations;
138 mmid(y, dydx, equations, x0, interval, steps, yFinal, derivs);
139 mmid(y, dydx, equations, x0, interval, steps / 2, yFinal2, derivs);
140 for (i = 0; i < equations; i++)
141 yFinal[i] = (4 * yFinal[i] - yFinal2[i]) / 3;
175 void (*derivs)(
double *dydx,
double *y,
double x),
190 double (*exit_func)(
double *dydx,
double *y,
double x),
194 static double *y0, *yscale;
195 static double *dydx0, *y1, *dydx1, *dydx2, *y2, *accur;
196 static long last_neq = 0;
197 double ex0, ex1, ex2, x1, x2;
199 long i, n_exit_iterations;
200#define MAX_N_STEP_UPS 10
203 return (DIFFEQ_XI_GT_XF);
204 if (FABS(*x0 - xf) < x_accuracy)
205 return (DIFFEQ_SOLVED_ALREADY);
207 if (last_neq < n_eq) {
218 y0 =
tmalloc(
sizeof(
double) * n_eq);
219 dydx0 =
tmalloc(
sizeof(
double) * n_eq);
220 y1 =
tmalloc(
sizeof(
double) * n_eq);
221 dydx1 =
tmalloc(
sizeof(
double) * n_eq);
222 y2 =
tmalloc(
sizeof(
double) * n_eq);
223 dydx2 =
tmalloc(
sizeof(
double) * n_eq);
227 for (i = 0; i < n_eq; i++)
231 (*derivs)(dydx0, y0, *x0);
232 ex0 = (*exit_func)(dydx0, y0, *x0);
236 if (FABS(ex0) < exit_accuracy) {
237 for (i = 0; i < n_eq; i++)
239 return (DIFFEQ_ZERO_FOUND);
243 if ((xdiff = xf - *x0) < h_step)
247 mmid2(y0, dydx0, n_eq, x1, h_step, 8, y1, derivs);
250 (*derivs)(dydx1, y1, x1);
251 ex1 = (*exit_func)(dydx1, y1, x1);
252 if (SIGN(ex0) != SIGN(ex1))
255 if (FABS(xdiff = xf - x1) < x_accuracy) {
257 for (i = 0; i < n_eq; i++)
260 return (DIFFEQ_END_OF_INTERVAL);
263 SWAP_PTR(dydx0, dydx1);
270 printf(
"failure in mmid_odeint3_na(): solution stepped outside interval\n");
271 return (DIFFEQ_OUTSIDE_INTERVAL);
274 if (FABS(ex1) < exit_accuracy) {
275 for (i = 0; i < n_eq; i++)
278 return (DIFFEQ_ZERO_FOUND);
282 n_exit_iterations = MAX_EXIT_ITERATIONS;
285 h_step = -ex0 * (x1 - *x0) / (ex1 - ex0) * ITER_FACTOR;
287 mmid2(y0, dydx0, n_eq, x2, h_step, 8, y2, derivs);
290 (*derivs)(dydx2, y2, x2);
291 ex2 = (*exit_func)(dydx2, y2, x2);
292 if (FABS(ex2) < exit_accuracy) {
293 for (i = 0; i < n_eq; i++)
296 return (DIFFEQ_ZERO_FOUND);
299 if (SIGN(ex1) == SIGN(ex2)) {
301 SWAP_PTR(dydx1, dydx2);
306 SWAP_PTR(dydx0, dydx2);
310 }
while (n_exit_iterations--);
311 return (DIFFEQ_EXIT_COND_FAILED);
long mmid_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)
Integrates ODEs until a condition is met or the interval is reached.
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.
void mmid2(double *y, double *dydx, long equations, double x0, double interval, long steps, double *yFinal, void(*derivs)(double *_dydx, double *_y, double _x))
Enhances the modified midpoint method with error correction.