module ERKStep:sig..end
ERKStep Time-Stepping Module for nonstiff initial value problems.
This module solves problems of the form $\dot{y} = f(t, y)$, $y(t_0) = y_0$.
Its interface is structured as follows.
include Arkode.Common
type('d, 'k)session =('d, 'k, Arkode.erkstep) session
A session with the ERKStep time-stepping solver.
An example session with ERKStep (arkode_erk_skel.ml):
open Sundials
module ERKStep = Arkode.ERKStep
(* 1. Define right-hand-side functions. *)
let f _t y yd = yd.{0} <- y.{1};
yd.{1} <- -9.81
(* 2. Optionally define a root function. *)
let g _t y gout = gout.{0} <- 1.0 -. y.{0}
(* 3. Set vector of initial values.
The length of this vector determines the problem size. *)
let yd = RealArray.of_list [ 10.0; 0.0 ]
let y = Nvector_serial.wrap yd
(* 4. Create and initialize a solver session. *)
let s = ERKStep.(init (SStolerances (1e-4, 1e-9)) f ~roots:(1, g) 0.0 y);;
(* 5. Set optional inputs, e.g.,
call [set_*] functions to change solver parameters. *)
ERKStep.set_stop_time s 10.0;;
ERKStep.set_all_root_directions s RootDirs.Increasing;;
(* 6. Advance the solution in time,
by repeatedly calling [evolve_normal] or [evolve_one_step]. *)
let rec go (t, r) =
Printf.printf "% .10e\t% .10e\t% .10e\n" t yd.{0} yd.{1};
match r with
| ERKStep.Success -> go (ERKStep.evolve_normal s (t +. 0.5) y)
| ERKStep.RootsFound -> begin
yd.{1} <- -0.8 *. yd.{1};
ERKStep.reinit s t y;
go (t, ERKStep.Success)
end
| ERKStep.StopTimeReached -> ();;
Printf.printf "time\ty\ty'\n";;
go (0.0, ERKStep.Success);;
(* 7. Get optional outputs,
call the [get_*] functions to examine solver statistics. *)
let ns = ERKStep.get_num_steps s
val init : ?context:Sundials.Context.t ->
('data, 'kind) tolerance ->
?order:int ->
'data rhsfn ->
?roots:int * 'data rootsfn ->
float -> ('data, 'kind) Nvector.t -> ('data, 'kind) sessionCreates and initializes a session with the solver. The call
init tol ~order f ~roots:(nroots, g) t0 y0has as arguments:
tol, the integration tolerances,order, the order of accuracy for the integration method,f, the ODE right-hand side function,nroots, the number of root functions,g, the root function ((nroots, g) defaults to
Arkode.Common.no_roots),t0, the initial value of the independent variable, andy0, a vector of initial values that also determines the number
of equations.This function does everything necessary to initialize a session, i.e.,
it makes the calls referenced below. The Arkode.ERKStep.evolve_normal and
Arkode.ERKStep.evolve_one_step functions may be called directly.
By default, the session is created using the context returned by
Sundials.Context.default, but this can be overridden by passing
an optional context argument.
val evolve_normal : ('d, 'k) session ->
float -> ('d, 'k) Nvector.t -> float * solver_resultIntegrates an ODE system over an interval. The call
tret, r = evolve_normal s tout yout has as arguments
s, a solver session,tout, the next time at which a solution is desired, and,yout, a vector to store the computed solution.It returns tret, the time reached by the solver, which will be equal to
tout if no errors occur, and, r, a Arkode.Common.solver_result.
IllInput Missing or illegal solver inputs.TooMuchWork The requested time could not be reached in mxstep internal steps.TooMuchAccuracy The requested accuracy could not be satisfied.ErrFailure Too many error test failures within a step or at the minimum step size.RhsFuncFailure Unrecoverable failure in one of the RHS functions.FirstRhsFuncFailure Initial unrecoverable failure in one of the RHS functions.RepeatedRhsFuncFailure Too many convergence test failures, or unable to estimate the initial step size, due to repeated recoverable errors in one of the right-hand side functions.UnrecoverableRhsFuncFailure One of the right-hand side functions had a recoverable error, but no recovery was possible. This error can only occur after an error test failure at order one.RootFuncFailure Failure in the rootfinding function g.val evolve_one_step : ('d, 'k) session ->
float -> ('d, 'k) Nvector.t -> float * solver_resultLike Arkode.ERKStep.evolve_normal but returns after one internal solver step.
val get_dky : ('d, 'k) session -> ('d, 'k) Nvector.t -> float -> int -> unitReturns the interpolated solution or derivatives.
get_dky s dky t k computes the kth derivative of the function
at time t, i.e.,
$\frac{d^\mathtt{k}}{\mathit{dt}^\mathtt{k}}y(\mathtt{t})$,
and stores it in dky. The arguments must satisfy
$t_n - h_n \leq \mathtt{t} \leq t_n$—where $t_n$
denotes Arkode.ERKStep.get_current_time and $h_n$ denotes Arkode.ERKStep.get_last_step,—
and $0 \leq \mathtt{k} \leq 3$.
This function may only be called after a successful return from either
Arkode.ERKStep.evolve_normal or Arkode.ERKStep.evolve_one_step.
BadT t is not in the interval $[t_n - h_n, t_n]$.BadK k is not in the range {0, 1, ..., dord}.val reinit : ('d, 'k) session ->
?order:int -> ?roots:int * 'd rootsfn -> float -> ('d, 'k) Nvector.t -> unitReinitializes the solver with new parameters and state values. The
values of the independent variable, i.e., the simulation time, and the
state variables must be given. If given, order changes the order of
accuracy, and roots specifies a new root finding function. The new
problem must have the same size as the previous one.
The allowed values for order are as for Arkode.ERKStep.init.
val reset : ('d, 'k) session -> float -> ('d, 'k) Nvector.tResets the state to the given independent variable value and dependent variable vector. All previously set options, internal counter values, and step-size/error histories are retained.
val resize : ('d, 'k) session ->
?resize_nvec:'d resize_fn ->
('d, 'k) tolerance -> float -> ('d, 'k) Nvector.t -> float -> unitChange the number of equations and unknowns between integrator steps.
The call
resize s ~resize_nvec:rfn tol hscale ynew t0
has as arguments:
s, the solver session to resize,rfn, a resize function that transforms nvectors in place-otherwise
they are simply destroyed and recloned,tol, tolerance values (ensures that any tolerance vectors are resized),hscale, the next step will be of size $h \mathtt{hscale}$,ynew, the newly-sized solution vector with the value $y(t_0)$, andt0, the current value of the independent variable $t_0$.A call to this function disables inequality constraint checking. It can
be reenabled by calling Arkode.ERKStep.set_constraints.
val set_tolerances : ('d, 'k) session -> ('d, 'k) tolerance -> unitSets the integration tolerances.
val set_defaults : ('d, 'k) session -> unitResets all optional input parameters to their default values. Neither the problem-defining functions nor the root-finding functions are changed.
val set_diagnostics : ?logfile:Sundials.Logfile.t -> ('d, 'k) session -> unitWrite step adaptivity and solver diagnostics on the standard output (or given file).
val set_interpolant_type : ('d, 'k) session -> interpolant_type -> unitSpecifies the interpolation module used for output value interpolation and implicit method predictors.
val set_interpolant_degree : ('d, 'k) session -> int -> unitSpecifies the degree of the polynomial interpolant used for output values and implicit method predictors.
val clear_diagnostics : ('d, 'k) session -> unitDo not write step adaptivity or solver diagnostics of a file.
val set_error_file : ('d, 'k) session -> Sundials.Logfile.t -> unitConfigure the default error handler to write messages to a file. By default it writes to Logfile.stderr.
val set_err_handler_fn : ('d, 'k) session ->
(Sundials.Util.error_details -> unit) -> unitSpecifies a custom function for handling error messages. The handler must not fail: any exceptions are trapped and discarded.
val clear_err_handler_fn : ('d, 'k) session -> unitRestores the default error handling function.
val set_init_step : ('d, 'k) session -> float -> unitSpecifies the initial step size.
val set_fixed_step : ('d, 'k) session -> float option -> unitDisables time step adaptivity and fix the step size for all internal
steps. Pass None to restore time step adaptivity. This function is not
recommended since there is no assurance of the validity of the computed
solutions. It is provided primarily for code-to-code verification
testing. Use in conjunction with Arkode.ERKStep.set_min_step or Arkode.ERKStep.set_max_step.
val set_max_num_constr_fails : ('d, 'k) session -> int -> unitSpecifies the maximum number of constraint failures in a step before an error is signalled.
val set_max_hnil_warns : ('d, 'k) session -> int -> unitSpecifies the maximum number of messages warning that t + h = t on
the next internal step.
val set_max_num_steps : ('d, 'k) session -> int -> unitSpecifies the maximum number of steps taken in attempting to reach a given output time.
val set_max_err_test_fails : ('d, 'k) session -> int -> unitSpecifies the maximum number of error test failures permitted in attempting one step.
val set_min_step : ('d, 'k) session -> float -> unitSpecifies a lower bound on the magnitude of the step size.
val set_max_step : ('d, 'k) session -> float -> unitSpecifies an upper bound on the magnitude of the step size.
val set_stop_time : ('d, 'k) session -> float -> unitLimits the value of the independent variable t when solving.
By default no stop time is imposed.
val set_table : ('d, 'k) session -> Arkode.ButcherTable.t -> unitSpecifies a customized Butcher table.
val set_table_num : ('d, 'k) session -> Arkode.ButcherTable.erk_table -> unitUse a specific built-in Butcher table for integration.
val set_adaptivity_method : ('d, 'k) session -> 'd adaptivity_method -> unitSpecifies the method and associated parameters used for time step adaptivity.
val set_cfl_fraction : ('d, 'k) session -> float -> unitSpecifies the fraction of the estimated explicitly stable step to use. Any non-positive argument resets to the default value (0.5).
val set_error_bias : ('d, 'k) session -> float -> unitSpecifies the bias to apply to the error estimates within accuracy-based adaptivity strategies.
val set_fixed_step_bounds : ('d, 'k) session -> float -> float -> unitSpecifies the step growth interval in which the step size will remain
unchanged. In the call set_fixed_step_bounds s lb ub, lb specifies a
lower bound on the window to leave the step size fixed and ub
specifies an upper bound. Any interval not containing 1.0 resets to
default values.
val set_max_efail_growth : ('d, 'k) session -> float -> unitSpecifies the maximum step size growth factor upon multiple successive accuracy-based error failures in the solver. Any value outside the interval $(0, 1]$ resets to the default value.
val set_max_first_growth : ('d, 'k) session -> float -> unitSpecifies the maximum allowed step size change following the very first integration step. Any value $\le 1$ resets to the default value.
val set_max_growth : ('d, 'k) session -> float -> unitSpecifies the maximum growth of the step size between consecutive time steps. Any value $\le 1$ resets to the default value.
val set_min_reduction : ('d, 'k) session -> float -> unitSpecifies the minimum allowed reduction factor in step size between step attempts that result from a temporal error failure in the integration process.
val set_safety_factor : ('d, 'k) session -> float -> unitSpecifies the safety factor to be applied to the accuracy-based estimated step. Any non-positive value resets to the default value.
val set_small_num_efails : ('d, 'k) session -> float -> unitSpecifies the threshold for “multiple” successive error failures
before the factor from Arkode.ERKStep.set_max_efail_growth is applied. Any
non-positive value resets to the default value.
val set_stability_fn : ('d, 'k) session -> 'd stability_fn -> unitSets a problem-dependent function to estimate a stable time step size for the explicit portion of the ODE system.
val clear_stability_fn : ('d, 'k) session -> unitClears the problem-dependent function that estimates a stable time step size for the explicit portion of the ODE system.
val set_postprocess_step_fn : ('d, 'k) session -> 'd postprocess_step_fn -> unitSet a post processing step function.
val clear_postprocess_step_fn : ('d, 'k) session -> unitClear the post processing step function.
val set_constraints : ('d, 'k) session -> ('d, 'k) Nvector.t -> unitSpecifies a vector defining inequality constraints for each
component of the solution vector. See Sundials.Constraint.
val get_work_space : ('d, 'k) session -> int * intReturns the real and integer workspace sizes.
real_size, integer_size)val get_num_steps : ('d, 'k) session -> intReturns the cumulative number of internal steps taken by the solver.
val get_num_exp_steps : ('d, 'k) session -> intReturns the cumulative number of stability-limited steps taken by the solver.
val get_num_acc_steps : ('d, 'k) session -> intReturns the cumulative number of accuracy-limited steps taken by the solver.
val get_num_step_attempts : ('d, 'k) session -> intReturns the cumulative number of steps attempted by the solver.
val get_num_rhs_evals : ('d, 'k) session -> intReturns the number of calls to the right-hand side function.
val get_num_err_test_fails : ('d, 'k) session -> intReturns the number of local error test failures that have occurred.
val get_last_step : ('d, 'k) session -> floatReturns the integration step size taken on the last successful internal step.
val get_current_step : ('d, 'k) session -> floatReturns the integration step size to be attempted on the next internal step.
val get_actual_init_step : ('d, 'k) session -> floatReturns the the value of the integration step size used on the first step.
val get_current_time : ('d, 'k) session -> floatReturns the the current internal time reached by the solver.
val get_current_butcher_table : ('d, 'k) session -> Arkode.ButcherTable.tReturns the Butcher table in use by the solver.
val get_tol_scale_factor : ('d, 'k) session -> floatReturns a suggested factor by which the user's tolerances should be scaled when too much accuracy has been requested for some internal step.
val get_err_weights : ('d, 'k) session -> ('d, 'k) Nvector.t -> unitReturns the solution error weights at the current time.
val get_est_local_errors : ('d, 'k) session -> ('d, 'k) Nvector.t -> unitReturns the vector of estimated local errors.
type timestepper_stats = {
|
exp_steps : |
(* | Cumulative number of stability-limited solver steps. | *) |
|
acc_steps : |
(* | Cumulative number of accuracy-limited solver steps. | *) |
|
step_attempts : |
(* | Cumulative number of steps attempted by the solver. | *) |
|
num_nf_evals : |
(* | Number of calls to the right-hand side function. | *) |
|
num_err_test_fails : |
(* | Number of error test failures. | *) |
}
Summaries of time-stepper statistics.
val get_timestepper_stats : ('d, 'k) session -> timestepper_statsReturns a grouped set of time-stepper statistics.
val get_step_stats : ('d, 'k) session -> step_statsReturns a grouped set of integrator statistics.
val print_timestepper_stats : ('d, 'k) session -> Stdlib.out_channel -> unitPrints time-stepper statistics on the given channel.
val print_step_stats : ('d, 'k) session -> Stdlib.out_channel -> unitPrints integrator statistics on the given channel.
val write_session : ?logfile:Sundials.Logfile.t -> ('d, 'k) session -> unitSummarize the session on the standard output (or given file).
val set_root_direction : ('d, 'k) session -> Sundials.RootDirs.d array -> unitset_root_direction s dir specifies the direction of zero-crossings to
be located and returned. dir may contain one entry for each root
function.
val set_all_root_directions : ('d, 'k) session -> Sundials.RootDirs.d -> unitLike Arkode.ERKStep.set_root_direction but specifies a single direction for all root
functions.
val set_no_inactive_root_warn : ('d, 'k) session -> unitDisables issuing a warning if some root function appears to be identically zero at the beginning of the integration.
val get_num_roots : ('d, 'k) session -> intReturns the number of root functions.
val get_root_info : ('d, 'k) session -> Sundials.Roots.t -> unitFills an array showing which functions were found to have a root.
val get_num_g_evals : ('d, 'k) session -> intReturns the cumulative number of calls made to the user-supplied root function g.
val get_num_constr_fails : ('d, 'k) session -> intReturns the cumulative number of test failures.
val write_parameters : ?logfile:Sundials.Logfile.t -> ('d, 'k) session -> unitOutputs all the solver parameters on the standard output (or given file).
val write_butcher : ?logfile:Sundials.Logfile.t -> ('d, 'k) session -> unitOutputs the current butcher table on the standard output (or given file).