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
|
/* ode-initval/evolve.c
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
/* Author: G. Jungman
*/
#include <config.h>
#include <string.h>
#include <stdlib.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_odeiv2.h>
#include "odeiv_util.h"
gsl_odeiv2_evolve *
gsl_odeiv2_evolve_alloc (size_t dim)
{
gsl_odeiv2_evolve *e =
(gsl_odeiv2_evolve *) malloc (sizeof (gsl_odeiv2_evolve));
if (e == 0)
{
GSL_ERROR_NULL ("failed to allocate space for evolve struct",
GSL_ENOMEM);
}
e->y0 = (double *) malloc (dim * sizeof (double));
if (e->y0 == 0)
{
free (e);
GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
}
e->yerr = (double *) malloc (dim * sizeof (double));
if (e->yerr == 0)
{
free (e->y0);
free (e);
GSL_ERROR_NULL ("failed to allocate space for yerr", GSL_ENOMEM);
}
e->dydt_in = (double *) malloc (dim * sizeof (double));
if (e->dydt_in == 0)
{
free (e->yerr);
free (e->y0);
free (e);
GSL_ERROR_NULL ("failed to allocate space for dydt_in", GSL_ENOMEM);
}
e->dydt_out = (double *) malloc (dim * sizeof (double));
if (e->dydt_out == 0)
{
free (e->dydt_in);
free (e->yerr);
free (e->y0);
free (e);
GSL_ERROR_NULL ("failed to allocate space for dydt_out", GSL_ENOMEM);
}
e->dimension = dim;
e->count = 0;
e->failed_steps = 0;
e->last_step = 0.0;
e->driver = NULL;
return e;
}
int
gsl_odeiv2_evolve_reset (gsl_odeiv2_evolve * e)
{
e->count = 0;
e->failed_steps = 0;
e->last_step = 0.0;
return GSL_SUCCESS;
}
void
gsl_odeiv2_evolve_free (gsl_odeiv2_evolve * e)
{
RETURN_IF_NULL (e);
free (e->dydt_out);
free (e->dydt_in);
free (e->yerr);
free (e->y0);
free (e);
}
/* Evolution framework method.
*
* Uses an adaptive step control object
*/
int
gsl_odeiv2_evolve_apply (gsl_odeiv2_evolve * e,
gsl_odeiv2_control * con,
gsl_odeiv2_step * step,
const gsl_odeiv2_system * dydt,
double *t, double t1, double *h, double y[])
{
const double t0 = *t;
double h0 = *h;
int step_status;
int final_step = 0;
double dt = t1 - t0; /* remaining time, possibly less than h */
if (e->dimension != step->dimension)
{
GSL_ERROR ("step dimension must match evolution size", GSL_EINVAL);
}
if ((dt < 0.0 && h0 > 0.0) || (dt > 0.0 && h0 < 0.0))
{
GSL_ERROR ("step direction must match interval direction", GSL_EINVAL);
}
/* Save y in case of failure in a step */
DBL_MEMCPY (e->y0, y, e->dimension);
/* Calculate initial dydt once or reuse previous value if the method
can benefit. */
if (step->type->can_use_dydt_in)
{
if (e->count == 0)
{
int status = GSL_ODEIV_FN_EVAL (dydt, t0, y, e->dydt_in);
if (status)
{
return status;
}
}
else
{
DBL_MEMCPY (e->dydt_in, e->dydt_out, e->dimension);
}
}
try_step:
if ((dt >= 0.0 && h0 > dt) || (dt < 0.0 && h0 < dt))
{
h0 = dt;
final_step = 1;
}
else
{
final_step = 0;
}
if (step->type->can_use_dydt_in)
{
step_status =
gsl_odeiv2_step_apply (step, t0, h0, y, e->yerr, e->dydt_in,
e->dydt_out, dydt);
}
else
{
step_status =
gsl_odeiv2_step_apply (step, t0, h0, y, e->yerr, NULL, e->dydt_out,
dydt);
}
/* Return if stepper indicates a pointer or user function failure */
if (step_status == GSL_EFAULT || step_status == GSL_EBADFUNC)
{
return step_status;
}
/* Check for stepper internal failure */
if (step_status != GSL_SUCCESS)
{
/* Stepper was unable to calculate step. Try decreasing step size. */
double h_old = h0;
h0 *= 0.5;
/* Check that an actual decrease in h0 occured and the
suggested h0 will change the time by at least 1 ulp */
{
double t_curr = GSL_COERCE_DBL (*t);
double t_next = GSL_COERCE_DBL ((*t) + h0);
if (fabs (h0) < fabs (h_old) && t_next != t_curr)
{
/* Step was decreased. Undo step, and try again with new h0. */
DBL_MEMCPY (y, e->y0, dydt->dimension);
e->failed_steps++;
goto try_step;
}
else
{
*h = h0; /* notify user of step-size which caused the failure */
*t = t0; /* restore original t value */
return step_status;
}
}
}
e->count++;
e->last_step = h0;
if (final_step)
{
*t = t1;
}
else
{
*t = t0 + h0;
}
if (con != NULL)
{
/* Check error and attempt to adjust the step. */
double h_old = h0;
const int hadjust_status
=
gsl_odeiv2_control_hadjust (con, step, y, e->yerr, e->dydt_out, &h0);
if (hadjust_status == GSL_ODEIV_HADJ_DEC)
{
/* Check that the reported status is correct (i.e. an actual
decrease in h0 occured) and the suggested h0 will change
the time by at least 1 ulp */
double t_curr = GSL_COERCE_DBL (*t);
double t_next = GSL_COERCE_DBL ((*t) + h0);
if (fabs (h0) < fabs (h_old) && t_next != t_curr)
{
/* Step was decreased. Undo step, and try again with new h0. */
DBL_MEMCPY (y, e->y0, dydt->dimension);
e->failed_steps++;
goto try_step;
}
else
{
/* Can not obtain required error tolerance, and can not
decrease step-size any further, so give up and return
GSL_FAILURE.
*/
*h = h0; /* notify user of step-size which caused the failure */
return GSL_FAILURE;
}
}
}
/* Suggest step size for next time-step. Change of step size is not
suggested in the final step, because that step can be very
small compared to previous step, to reach t1.
*/
if (final_step == 0)
{
*h = h0;
}
return step_status;
}
/* Evolves the system using the user specified constant step size h.
*/
int
gsl_odeiv2_evolve_apply_fixed_step (gsl_odeiv2_evolve * e,
gsl_odeiv2_control * con,
gsl_odeiv2_step * step,
const gsl_odeiv2_system * dydt,
double *t, const double h, double y[])
{
const double t0 = *t;
int step_status;
if (e->dimension != step->dimension)
{
GSL_ERROR ("step dimension must match evolution size", GSL_EINVAL);
}
/* Save y in case of failure in a step */
DBL_MEMCPY (e->y0, y, e->dimension);
/* Calculate initial dydt once if the method can benefit. */
if (step->type->can_use_dydt_in)
{
int status = GSL_ODEIV_FN_EVAL (dydt, t0, y, e->dydt_in);
if (status)
{
return status;
}
}
if (step->type->can_use_dydt_in)
{
step_status =
gsl_odeiv2_step_apply (step, t0, h, y, e->yerr, e->dydt_in,
e->dydt_out, dydt);
}
else
{
step_status =
gsl_odeiv2_step_apply (step, t0, h, y, e->yerr, NULL, e->dydt_out,
dydt);
}
/* Return the stepper return value in case of an error */
if (step_status != GSL_SUCCESS)
{
return step_status;
}
if (con != NULL)
{
/* Calculate error level. Fail if error level exceeds desired
error level. */
double htemp = h;
const int hadjust_status
= gsl_odeiv2_control_hadjust (con, step, y, e->yerr,
e->dydt_out, &htemp);
if (hadjust_status == GSL_ODEIV_HADJ_DEC)
{
DBL_MEMCPY (y, e->y0, dydt->dimension);
e->failed_steps++;
return GSL_FAILURE;
}
}
/* Step is accepted, update status */
e->count++;
e->last_step = h;
*t = t0 + h;
return GSL_SUCCESS;
}
int
gsl_odeiv2_evolve_set_driver (gsl_odeiv2_evolve * e,
const gsl_odeiv2_driver * d)
{
if (d != NULL)
{
e->driver = d;
}
else
{
GSL_ERROR_NULL ("driver pointer is null", GSL_EFAULT);
}
return GSL_SUCCESS;
}
|