module Common:sig
..end
Common definitions that are included in each of the time-stepping modules.
type
solver_result =
| |
Success |
(* | The solution was advanced. (ARK_SUCCESS) | *) |
| |
RootsFound |
(* | A root was found. See | *) |
| |
StopTimeReached |
(* | The stop time was reached. See
| *) |
Values returned by the step functions. Failures are indicated by exceptions.
type
step_stats = {
|
num_steps : |
(* | Cumulative number of internal solver steps. | *) |
|
actual_init_step : |
(* | Integration step sized used on the first step. | *) |
|
last_step : |
(* | Integration step size of the last internal step. | *) |
|
current_step : |
(* | Integration step size to attempt on the next internal step. | *) |
|
current_time : |
(* | Current internal time reached by the solver. | *) |
}
Summaries of integrator statistics.
type
interpolant_type =
| |
Hermite |
(* | Polynomial interpolants of Hermite form. (ARK_INTERP_HERMITE) | *) |
| |
Lagrange |
(* | Polynomial interpolants of Lagrange form, for stiff problems. (ARK_INTERP_LAGRANGE) | *) |
Used to specify the nterpolation method used for output values and
implicit method predictors. See, for example,
Arkode.MRIStep.set_interpolant_type
.
type
linearity =
| |
Linear of |
(* | Implicit portion is linear. Specifies whether $f_I(t, y)$ or the preconditioner is time-dependent. | *) |
| |
Nonlinear |
(* | Implicit portion is nonlinear. | *) |
The linearity of the implicit portion of the problem.
type'data
error_weight_fun ='data -> 'data -> unit
Functions that set the multiplicative error weights for use in the weighted
RMS norm. The call efun y ewt
takes the dependent variable vector y
and
fills the error-weight vector ewt
with positive values or raises
Sundials.NonPositiveEwt
. Other exceptions are eventually propagated, but
should be avoided (efun
is not allowed to abort the solver).
type ('data, 'kind)
tolerance =
| |
SStolerances of |
(* |
| *) |
| |
SVtolerances of |
(* |
| *) |
| |
WFtolerances of |
(* | Set the multiplicative error weights for the weighted RMS norm. | *) |
Tolerance specifications.
val default_tolerances : ('data, 'kind) tolerance
A default relative tolerance of 1.0e-4 and absolute tolerance of 1.0e-9.
type
predictor_method =
| |
TrivialPredictor |
(* | Piece-wise constant interpolant $p_0(\tau) = y_{n-1}$ % | *) |
| |
MaximumOrderPredictor |
(* | An interpolant $p_q(t)$ of polynomial order up to $q=3$. | *) |
| |
VariableOrderPredictor |
(* | Decrease the polynomial degree for later RK stages. | *) |
| |
CutoffOrderPredictor |
(* | Maximum order for early RK stages and first-order for later ones. | *) |
Method choices for predicting implicit solutions.
type 'd
nonlin_system_data = {
|
tcur : |
(* | Independent variable value for slow stage $t^S_{n,i}$ . | *) |
|
zpred : |
(* | Predicted nonlinear solution $z_{\mathit{pred}}$ . This data must not be changed. | *) |
|
zi : |
(* | Stage vector $z_i$ . This data may not be current and may need to be filled. | *) |
|
fi : |
(* | Memory available for evaluating the slow right-hand side $f^S(t^S_{n,i}, z_i)$ . This data may not be current and may need to be filled. | *) |
|
gamma : |
(* | Current $\gamma$ for slow-stage calculation. | *) |
|
sdata : |
(* | Accumulated data from previous solution and stages $\tilde{a}_i$ . This data must not be changed. | *) |
}
Internal data required to construct the current nonlinear implicit system within a nonlinear solver.
type'd
rootsfn =float -> 'd -> Sundials.RealArray.t -> unit
Called by the solver to calculate the values of root functions. These ‘zero-crossings’ are used to detect significant events. The function is passed three arguments:
t
, the value of the independent variable, i.e., the simulation time,y
, the vector of dependent-variable values, i.e., $y(t)$, and,gout
, a vector for storing the value of $g(t, y)$.y
and gout
should not be accessed after the function has
returned.val no_roots : int * 'd rootsfn
A convenience value for signalling that there are no roots to monitor.
type
adaptivity_args = {
|
h1 : |
(* | the current step size, $t_m - t_{m-1}$. | *) |
|
h2 : |
(* | the previous step size, $t_{m-1} - t_{m-2}$. | *) |
|
h3 : |
(* | the step size $t_{m-2} - t_{m-3}$. | *) |
|
e1 : |
(* | the error estimate from the current step, $m$. | *) |
|
e2 : |
(* | the error estimate from the previous step, $m-1$. | *) |
|
e3 : |
(* | the error estimate from the step $m-2$. | *) |
|
q : |
(* | the global order of accuracy for the integration method. | *) |
|
p : |
(* | the global order of accuracy for the embedding. | *) |
}
Arguments for Arkode.Common.adaptivity_fn
s.
type'd
adaptivity_fn =float -> 'd -> adaptivity_args -> float
A function implementing a time step adaptivity algorithm that chooses an
$h$ that satisfies the error tolerances. The call hnew = adapt_fn t y args
has as arguments
t
, the value of the independent variable,y
, the value of the dependent variable vector $y(t)$, andargs
, information on step sizes, error estimates, and accuracies.
and returns the next step size hnew
. The function should raise an
exception if it cannot set the next step size. The step size should be the
maximum value where the error estimates remain below 1.This function should focus on accuracy-based time step estimation; for
stability based time steps, Arkode.ARKStep.set_stability_fn
and
Arkode.ERKStep.set_stability_fn
should be used.
type
adaptivity_params = {
|
ks : |
|
method_order : |
}
Parameters for the standard adaptivity algorithms. There are two:
adaptivity_ks
, the k1
, k2
, and k3
parameters, or None
to use the defaults, andadaptivity_method_order
, true
specifies the method order of
accuracy $q$ and false
specifies the embedding order of
accuracy $p$.type 'd
adaptivity_method =
| |
PIDcontroller of |
(* | The default time adaptivity controller. | *) |
| |
PIcontroller of |
(* | Uses the two most recent step sizes. | *) |
| |
Icontroller of |
(* | Standard time adaptivity control algorithm. | *) |
| |
ExplicitGustafsson of |
(* | Primarily used with explicit RK methods. | *) |
| |
ImplicitGustafsson of |
(* | Primarily used with implicit RK methods. | *) |
| |
ImExGustafsson of |
(* | An ImEx version of the two preceding controllers. | *) |
| |
AdaptivityFn of |
(* | A custom time-step adaptivity function. | *) |
Asymptotic error control algorithms.
type'd
rhsfn =float -> 'd -> 'd -> unit
Right-hand side functions for calculating ODE derivatives. They are passed three arguments:
t
, the value of the independent variable, i.e., the simulation time,y
, the vector of dependent-variable values, i.e., $y(t)$, and,y'
, a vector for storing the value of $f(t, y)$.Within the function, raising a Sundials.RecoverableFailure
exception
indicates a recoverable error. Any other exception is treated as an
unrecoverable error.
y
and y'
should not be accessed after the function returns.type'd
stability_fn =float -> 'd -> float
A function that predicts the maximum stable step size for the explicit
portions of an ImEx ODE system. The call hstab = stab_fn t y
has as arguments
t
, the value of the independent variable, andy
, the value of the dependent variable vector $y(t)$.
and returns the absolute value of the maximum stable step size hstab
.
Returning $\mathtt{hstab}\le 0$ indicates that there is no explicit
stability restriction on the time step size. The function should raise
an exception if it cannot set the next step size.type'd
resize_fn ='d -> 'd -> unit
Called to resize a vector to match the dimensions of another. The call
resizefn y ytemplate
must resize y
to match the size of ytemplate
.
y
and ytemplate
should not be accessed after the function
has returned.type'd
postprocess_step_fn =float -> 'd -> unit
A function to process the results of each timestep solution. The arguments are
t
, the value of the independent variable, andy
, the value of the dependent variable vector $y(t)$.Definitions common to linear solvers within the time-stepping modules.
type'd
triple ='d * 'd * 'd
Workspaces with three temporary vectors.
type('t, 'd)
jacobian_arg =('t, 'd) jacobian_arg
= {
|
jac_t : |
(* | The independent variable. | *) |
|
jac_y : |
(* | The dependent variable vector. | *) |
|
jac_fy : |
(* | The derivative vector $f_I(t, y)$. | *) |
|
jac_tmp : |
(* | Workspace data. | *) |
}
Arguments common to Jacobian callback functions.