module Adjoint:sig
..end
(Adjoint) Sensitivity analysis of ODEs with respect to their parameters.
Provides an alternative to forward sensitivity analysis, which can become prohibitively expensive. This technique does not calculate sensitivities, but rather gradients with respect to the parameters of a relatively few derived functionals of the solution, that is the gradient $\frac{\mathrm{d}G}{\mathrm{d}p}$ of $G(p) = \int_{t_0}^T \! g(t, y, p)\,\mathrm{d}t$. The gradients are evaluated by first calculating forward and checkpointing certain intermediate state values, and then integrating backward to $t_0$.
This documented interface is structured as follows.
type('data, 'kind)
bsession =('data, 'kind) AdjointTypes.bsession
A backward session with the CVODES solver. Multiple backward sessions may be associated with a single parent session.
type[> Nvector_serial.kind ]
serial_bsession =(Nvector_serial.data, [> Nvector_serial.kind ] as 'a) bsession
Alias for backward sessions based on serial nvectors.
type
interpolation =
| |
IPolynomial |
(* | (CV_POLYNOMIAL) | *) |
| |
IHermite |
(* | (CV_HERMITE) | *) |
Specifies the type of interpolation to use between checkpoints.
val init : ('d, 'k) Cvode.session -> int -> interpolation -> unit
Activates the forward-backward problem. The arguments specify the number of integration steps between consecutive checkpoints, and the type of variable-degree interpolation.
val forward_normal : ('d, 'k) Cvode.session ->
float -> ('d, 'k) Nvector.t -> float * int * Cvode.solver_result
Integrates the forward problem over an interval and saves
checkpointing data. The arguments are the next time at which a solution
is desired (tout
) and a vector to receive the computed result
(y
). The function returns a triple tret, ncheck, sr
: the time
reached by the solver, the cumulative number of checkpoints stored, and
whether tout
was reached. The solver takes internal steps until it
has reached or just passed the tout
parameter (CV_NORMAL),
it then interpolates to approximate y(tout)
.
Cvode.IllInput
One of the inputs is invalid.Cvode.TooMuchWork
Could not reach tout
in mxstep
stepsCvode.TooMuchAccuracy
Could not satisfy the demanded accuracyCvode.ErrFailure
Too many error test failures.Cvode.ConvergenceFailure
Too many convergence test failures.Cvode.LinearSetupFailure
Unrecoverable failure in linear solver setup function.Cvode.LinearSolveFailure
Unrecoverable failure in linear solver solve function.AdjointNotInitialized
The Cvodes.Adjoint.init
function has not been called.val forward_one_step : ('d, 'k) Cvode.session ->
float -> ('d, 'k) Nvector.t -> float * int * Cvode.solver_result
Integrates the forward problem over an interval and saves
checkpointing data. The arguments are the next time at which a solution
is desired (tout
) and a vector to receive the computed result
(yret
). The function returns a triple tret, ncheck, sr
: the time
reached by the solver, the cumulative number of checkpoints stored, and
whether tout
was reached. The solver takes one step
(CV_ONE_STEP) and returns the solution reached.
Cvode.IllInput
One of the inputs is invalid.Cvode.TooMuchWork
Could not reach tout
in mxstep
stepsCvode.TooMuchAccuracy
Could not satisfy the demanded accuracyCvode.ErrFailure
Too many error test failures.Cvode.ConvergenceFailure
Too many convergence test failures.Cvode.LinearSetupFailure
Unrecoverable failure in linear solver setup function.Cvode.LinearSolveFailure
Unrecoverable failure in linear solver solve function.AdjointNotInitialized
The Cvodes.Adjoint.init
function has not been called.type('data, 'kind)
linear_solver =('data, 'kind) AdjointTypes.linear_solver
Linear solvers used in backward problems.
type[> Nvector_serial.kind ]
serial_linear_solver =(Nvector_serial.data, [> Nvector_serial.kind ] as 'a)
linear_solver
Alias for linear solvers that are restricted to serial nvectors.
type'd
triple ='d * 'd * 'd
Workspaces with three temporary vectors.
type('t, 'd)
jacobian_arg =('t, 'd) AdjointTypes.jacobian_arg
= {
|
jac_t : |
(* | The independent variable. | *) |
|
jac_y : |
(* | The forward solution vector. | *) |
|
jac_yb : |
(* | The backward solution vector. | *) |
|
jac_fyb : |
(* | The backward right-hand side function | *) |
|
jac_tmp : |
(* | Temporary storage vectors. | *) |
}
Arguments common to Jacobian callback functions.
module Diag:sig
..end
Diagonal approximation of Jacobians by difference quotients.
module Dls:sig
..end
Direct Linear Solvers operating on dense, banded, and sparse matrices.
module Spils:sig
..end
Scaled Preconditioned Iterative Linear Solvers.
val matrix_embedded_solver : (unit, 'data, 'kind, [> `MatE ]) Sundials.LinearSolver.t ->
('data, 'kind) linear_solver
Create a CVode-specific linear solver from a generic matrix embedded solver.
type'd
brhsfn_args ='d AdjointTypes.brhsfn_args
= {
|
t : |
(* | The value of the independent variable. | *) |
|
y : |
(* | The vector of dependent-variable values $y(t)$. | *) |
|
yb : |
(* | The vector of backward dependent-variable values $y_B(t)$. | *) |
}
Arguments common to Cvodes.Adjoint.brhsfn_no_sens
and Cvodes.Adjoint.brhsfn_with_sens
.
type'd
brhsfn_no_sens ='d brhsfn_args -> 'd -> unit
Backward functions without forward sensitivities. They are passed the arguments:
args
, the current values of forward and backward variables, and,yb'
, a vector for storing the values
$\dot{y}_B = f_B(t, y, y_B)$.Within the function, raising a Sundials.RecoverableFailure
exception
indicates a recoverable error. Any other exception is treated as an
unrecoverable error.
type'd
brhsfn_with_sens ='d brhsfn_args -> 'd array -> 'd -> unit
Backward functions with forward sensitivities. They are passed the arguments:
args
, the current values of state and backward sensitivity variables,s
, an array holding the values of forward sensitivity vectors, and,yb'
, a vector for storing the values
$\dot{y}_B = f_B(t, y, y_S, y_B)$.Within the function, raising a Sundials.RecoverableFailure
exception
indicates a recoverable error. Any other exception is treated as an
unrecoverable error.
type 'd
brhsfn =
| |
NoSens of |
(* | No dependency on forward sensitivities. | *) |
| |
WithSens of |
(* | Dependency on forward sensitivities. | *) |
Functions that evaluate the right-hand side of a backward ODE system with or without forward sensitivities.
type ('d, 'k)
tolerance =
| |
SStolerances of |
(* |
| *) |
| |
SVtolerances of |
(* |
| *) |
Tolerance specifications.
val init_backward : ('d, 'k) Cvode.session ->
Cvode.lmm ->
('d, 'k) tolerance ->
?nlsolver:('d, 'k, ('d, 'k) Cvode.session, [ `Nvec ])
Sundials.NonlinearSolver.t ->
?lsolver:('d, 'k) linear_solver ->
'd brhsfn ->
float -> ('d, 'k) Nvector.t -> ('d, 'k) bsession
Creates and initializes a backward session attached to an existing (forward) session. The call
init_backward s lmm iter tol fb tb0 yb0
has as arguments:
s
, the parent (forward) session,lmm
, the linear multistep method (see Cvode.lmm
),tol
, the integration tolerances,nlsolver
, the solver to use to calculate integration steps,lsolver
, used by nlsolver
s based on Newton interation,fb
, the backward right-hand side function,tb0
, specifies the endpoint where final conditions are provided
for the backward problem, which is normally the endpoint of
forward integration, and,yb0
, a vector of final values that also determines the number
of equations.This function does everything necessary to initialize a backward
session, i.e., it makes the calls referenced below. The
Cvodes.Adjoint.backward_normal
and Cvodes.Adjoint.backward_one_step
functions may be called
directly.
If an nlsolver
is not specified, then the
Newton module is used by default.
In this case only, lsolver
defaults to Cvodes.Adjoint.Diag.solver
if not otherwise
specified. Specifying an nlsolver
that requires a linear solver without
specifying an lsolver
results in a Cvode.NonlinearInitFailure
(or
Cvode.IllInput
for Sundials < 4.0.0) exception on the first call to
Cvodes.Adjoint.backward_normal
or Cvodes.Adjoint.backward_one_step
.
AdjointNotInitialized
The Cvodes.Adjoint.init
function has not been called.BadFinalTime
The final time is outside the interval over which the forward problem was solved.module Quadrature:sig
..end
Support for backward quadrature equations that may or may not depend on forward sensitivities.
val backward_normal : ('d, 'k) Cvode.session -> float -> unit
Integrates a backward ODE system over an interval. The solver takes internal steps until it has reached or just passed the specified value.
AdjointNotInitialized
The Cvodes.Adjoint.init
function has not been called.NoBackwardProblem
The Cvodes.Adjoint.init_backward
function has not been called.NoForwardCall
Neither Cvodes.Adjoint.forward_normal
nor Cvodes.Adjoint.forward_one_step
has been called.Cvode.IllInput
One of the inputs is invalid.Cvode.TooMuchWork
Could not reach tout
in mxstep
stepsCvode.TooMuchAccuracy
Could not satisfy the demanded accuracyCvode.ErrFailure
Too many error test failures.Cvode.ConvergenceFailure
Too many convergence test failures.Cvode.LinearSetupFailure
Unrecoverable failure in linear solver setup function.Cvode.LinearSolveFailure
Unrecoverable failure in linear solver solve function.BadOutputTime
The requested output time is outside the interval over which the forward problem was solved.ForwardReinitializationFailed
Reinitialization of the forward problem failed at the first checkpoint (corresponding to the initial time of the forward problem).ForwardFail
An error occurred during the integration of the forward problem.val backward_one_step : ('d, 'k) Cvode.session -> float -> unit
Like Cvodes.Adjoint.backward_normal
but returns after one internal solver step.
val get : ('d, 'k) bsession -> ('d, 'k) Nvector.t -> float
Fills the given vector with the solution of the backward ODE problem at the returned time, interpolating if necessary.
val get_dky : ('d, 'k) bsession ->
('d, 'k) Nvector.t -> float -> int -> unit
Returns the interpolated solution or derivatives.
get_dky s dkyb t k
computes the k
th derivative of the backward
function at time t
, i.e.,
$\frac{d^\mathtt{k}y_B(\mathtt{t})}{\mathit{dt}^\mathtt{k}}$,
and stores it in dkyb
. The arguments must satisfy
$t_n - h_u \leq \mathtt{t} \leq t_n$—where $t_n$
denotes Cvodes.Adjoint.get_current_time
and $h_u$ denotes Cvodes.Adjoint.get_last_step
,—
and $0 \leq \mathtt{k} \leq q_u$—where $q_u$ denotes
Cvodes.Adjoint.get_last_order
.
This function may only be called after a successful return from either
Cvodes.Adjoint.backward_normal
or Cvodes.Adjoint.backward_one_step
.
BadT
t
is not in the interval $[t_n - h_u, t_n]$.BadK
k
is not in the range 0, 1, ..., $q_u$.val get_y : ('d, 'k) Cvode.session -> ('d, 'k) Nvector.t -> float -> unit
Fills the vector with the interpolated forward solution at the given time during a backward simulation.
val reinit : ('d, 'k) bsession ->
?nlsolver:('d, 'k, ('d, 'k) Cvode.session, [ `Nvec ])
Sundials.NonlinearSolver.t ->
?lsolver:('d, 'k) linear_solver ->
float -> ('d, 'k) Nvector.t -> unit
Reinitializes the backward problem with new parameters and state values. The values of the independent variable, i.e., the simulation time, and the state variables must be given. It is also possible to change the solution method (and linear solver).
AdjointNotInitialized
The Cvodes.Adjoint.init
function has not been called.BadFinalTime
The final time is not within the forward problem solution interval.val set_no_sensitivity : ('d, 'k) Cvode.session -> unit
Cancels the storage of sensitivity checkpointing data during forward
solution (with Cvodes.Adjoint.forward_normal
or Cvodes.Adjoint.forward_one_step
).
val set_tolerances : ('d, 'k) bsession -> ('d, 'k) tolerance -> unit
Sets the integration tolerances for the backward problem.
val set_max_ord : ('d, 'k) bsession -> int -> unit
Specifies the maximum order of the linear multistep method.
val set_max_num_steps : ('d, 'k) bsession -> int -> unit
Specifies the maximum number of steps taken in attempting to reach a given output time.
val set_init_step : ('d, 'k) bsession -> float -> unit
Specifies the initial step size.
val set_min_step : ('d, 'k) bsession -> float -> unit
Specifies a lower bound on the magnitude of the step size.
val set_max_step : ('d, 'k) bsession -> float -> unit
Specifies an upper bound on the magnitude of the step size.
val set_stab_lim_det : ('d, 'k) bsession -> bool -> unit
Indicates whether the BDF stability limit detection algorithm should be used.
val set_constraints : ('d, 'k) bsession -> ('d, 'k) Nvector.t -> unit
Specifies a vector defining inequality constraints for each
component of the solution vector y
. See Sundials.Constraint
.
val clear_constraints : ('d, 'k) bsession -> unit
Disables constraint checking.
val get_work_space : ('d, 'k) bsession -> int * int
Returns the real and integer workspace sizes.
real_size
, integer_size
)val get_num_steps : ('d, 'k) bsession -> int
Returns the cumulative number of internal steps taken by the solver.
val get_num_rhs_evals : ('d, 'k) bsession -> int
Returns the number of calls to the backward right-hand side function.
val get_num_lin_solv_setups : ('d, 'k) bsession -> int
Returns the number of calls made to the linear solver's setup function.
val get_num_err_test_fails : ('d, 'k) bsession -> int
Returns the number of local error test failures that have occurred.
val get_last_order : ('d, 'k) bsession -> int
Returns the integration method order used during the last internal step.
val get_current_order : ('d, 'k) bsession -> int
Returns the integration method order to be used on the next internal step.
val get_last_step : ('d, 'k) bsession -> float
Returns the integration step size taken on the last internal step.
val get_current_step : ('d, 'k) bsession -> float
Returns the integration step size to be attempted on the next internal step.
val get_actual_init_step : ('d, 'k) bsession -> float
Returns the the value of the integration step size used on the first step.
val get_current_time : ('d, 'k) bsession -> float
Returns the the current internal time reached by the solver.
val get_num_stab_lim_order_reds : ('d, 'k) bsession -> int
Returns the number of order reductions dictated by the BDF stability limit detection algorithm.
val get_tol_scale_factor : ('d, 'k) bsession -> float
Returns 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) bsession -> ('d, 'k) Nvector.t -> unit
Returns the solution error weights at the current time.
val get_est_local_errors : ('d, 'k) bsession -> ('d, 'k) Nvector.t -> unit
Returns the vector of estimated local errors.
val get_integrator_stats : ('d, 'k) bsession -> Cvode.integrator_stats
Returns the integrator statistics as a group.
val print_integrator_stats : ('d, 'k) bsession -> Stdlib.out_channel -> unit
Prints the integrator statistics on the given channel.
val get_num_nonlin_solv_iters : ('d, 'k) bsession -> int
Returns the cumulative number of nonlinear (functional or Newton) iterations.
val get_num_nonlin_solv_conv_fails : ('d, 'k) bsession -> int
Returns the cumulative number of nonlinear convergence failures.
val get_nonlin_solv_stats : ('d, 'k) bsession -> int * int
Returns both the numbers of nonlinear iterations performed nniters
and
nonlinear convergence failures nncfails
.
nniters
, nncfails
)exception AdjointNotInitialized
Adjoint sensitivity analysis was not initialized.
exception NoForwardCall
Neither Cvodes.Adjoint.forward_normal
nor Cvodes.Adjoint.forward_one_step
has been called.
exception ForwardReinitFailure
Reinitialization of the forward problem failed at the first checkpoint (corresponding to the initial time of the forward problem).
exception ForwardFailure
An error occured when integrating the forward problem from a checkpoint.
exception NoBackwardProblem
No backward problem has been created.
exception BadFinalTime
The final time was outside the interval over which the forward problem was solved.