DynamicBoundspODEsDiscrete.AM2Type
struct AM2 <: DynamicBoundspODEsDiscrete.Wilhelm2019Type

Use an second-order Adam's Moulton method style of relaxation.

DynamicBoundspODEsDiscrete.BDF2Type
struct BDF2 <: DynamicBoundspODEsDiscrete.Wilhelm2019Type

Use an second-order Backward Difference Formula method style of relaxation.

DynamicBoundspODEsDiscrete.CallbackHType
mutable struct CallbackH{V, F, T<:DynamicBoundspODEsDiscrete.Wilhelm2019Type} <: Function

A callback function used for the Wilhelm2019 integrator.

DynamicBoundspODEsDiscrete.DiscretizeRelaxType

DiscretizeRelax

An integrator for discretize and relaxation techniques.

  • x0f::Any

    Initial Condition for pODEs

  • Jx!::Any

    Jacobian w.r.t x

  • Jp!::Any

    Jacobian w.r.t p

  • p::Vector{Float64}

    Parameter value for pODEs

  • pL::Vector{Float64}

    Lower Parameter Bounds for pODEs

  • pU::Vector{Float64}

    Upper Parameter Bounds for pODEs

  • nx::Int64

    Number of state variables

  • np::Int64

    Number of decision variables

  • tspan::Tuple{Float64, Float64}

    Time span to integrate over

  • tsupports::Vector{Float64}

    Individual time points to evaluate

  • next_support_i::Int64

  • next_support::Float64

  • step_limit::Int64

    Maximum number of integration steps

  • step_count::Int64

    Steps taken

  • storage::Array{Vector{T}, 1} where T<:Number

    Stores solution X (from step 2) for each time

  • storage_apriori::Array{Vector{T}, 1} where T<:Number

    Stores solution X (from step 1) for each time

  • time::Vector{Float64}

    Stores each time t

  • support_dict::Dict{Int64, Int64}

    Support index to storage dictory

  • error_code::DynamicBoundsBase.TerminationStatusCode

    Holds data for numeric error encountered in integration step

  • P::Vector{T} where T<:Number

    Storage for bounds/relaxation of P

  • rP::Vector{T} where T<:Number

    Storage for bounds/relaxation of P - p

  • style::Number

    Relaxation Type

  • skip_step2::Bool

    Flag indicating that only apriori bounds should be computed

  • storage_buffer_size::Int64

  • print_relax_time::Bool

  • set_tf!::DynamicBoundspODEsDiscrete.TaylorFunctor!{F, K, S, T} where {T<:Number, S<:Real, F, K}

    Functor for evaluating Taylor coefficients over a set

  • method_f!::DynamicBoundspODEsDiscrete.AbstractStateContractor

  • exist_result::ExistStorage{F, K, S, T} where {T<:Number, S<:Real, F, K}

  • contractor_result::ContractorStorage{T} where T<:Number

  • step_result::StepResult{T} where T<:Number

  • step_params::StepParams

  • new_decision_pnt::Bool

  • new_decision_box::Bool

  • relax_t_dict_indx::Dict{Int64, Int64}

  • relax_t_dict_flt::Dict{Float64, Int64}

  • calculate_local_sensitivity::Bool

  • local_problem_storage::Any

  • constant_state_bounds::Union{Nothing, DynamicBoundsBase.ConstantStateBounds}

  • polyhedral_constraint::Union{Nothing, DynamicBoundsBase.PolyhedralConstraint}

  • prob::Any

DynamicBoundspODEsDiscrete.HermiteObreschkoffType

HermiteObreschkoff

A structure that stores the cofficient of the (P,Q)-Hermite-Obreschkoff method. (Offset due to method being zero indexed and Julia begin one indexed). The constructor HermiteObreschkoff(p::Val{P}, q::Val{Q}) where {P, Q} or HermiteObreschkoff(p::Int, q::Int) are used for the (P,Q)-method.

  • cpq::Vector{Float64}

    Cpq[i=1:p] index starting at i = 1 rather than 0

  • cqp::Vector{Float64}

    Cqp[i=1:q] index starting at i = 1 rather than 0

  • γ::Float64

    gamma for method

  • p::Int64

    Explicit order Hermite-Obreschkoff

  • q::Int64

    Implicit order Hermite-Obreschkoff

  • k::Int64

    Total order Hermite-Obreschkoff

  • skip_contractor::Bool

    Skips the contractor step of the Hermite Obreshkoff Contractor if set to true

DynamicBoundspODEsDiscrete.HermiteObreschkoffFunctorType

HermiteObreschkoffFunctor

A functor used in computing bounds and relaxations via Hermite-Obreschkoff's method. The implementation of the parametric Hermite-Obreschkoff's method based on the non-parametric version given in (1).

  1. [Nedialkov NS, and Jackson KR. "An interval Hermite-Obreschkoff method for

computing rigorous bounds on the solution of an initial value problem for an ordinary differential equation." Reliable Computing 5.3 (1999): 289-310.](https://link.springer.com/article/10.1023/A:1009936607335)

  1. [Nedialkov NS. "Computing rigorous bounds on the solution of an

initial value problem for an ordinary differential equation." University of Toronto. 2000.](http://citeseerx.ist.psu.edu/viewdoc/download?doi=10.1.1.633.9654&rep=rep1&type=pdf)

DynamicBoundspODEsDiscrete.JacTaylorFunctor!Type

JacTaylorFunctor!

A callable structure used to evaluate the Jacobian of Taylor cofficients. This also contains some addition fields to be used as inplace storage when computing and preconditioning paralleliped based methods to representing enclosure of the pODEs (Lohner's QR, Hermite-Obreschkoff, etc.). The constructor given by JacTaylorFunctor!(g!, nx::Int, np::Int, k::Val{K}, t::T, q::Q) may be used were type T should use type Q for internal computations. The order of the TaylorSeries is k, the right-hand side function is g!, nx is the number of state variables, np is the number of parameters.

  • g!::Function

    Right-hand side function for pODE which operates in place as g!(dx,x,p,t)

  • nx::Int64

    Dimensionality of x

  • np::Int64

    Dimensionality of p

  • s::Int64

    Order of TaylorSeries

  • out::Vector{S} where S<:Real

    In-place temporary storage for Taylor coefficient calculation

  • y::Vector{S} where S<:Real

    Variables y = (x,p)

  • x::Array{ForwardDiff.Dual{Nothing, S, NY}, 1} where {S<:Real, NY}

    State variables x

  • p::Array{ForwardDiff.Dual{Nothing, S, NY}, 1} where {S<:Real, NY}

    Decision variables p

  • Jxsto::Matrix{S} where S<:Real

    Storage for sum of Jacobian w.r.t x

  • Jpsto::Matrix{S} where S<:Real

    Storage for sum of Jacobian w.r.t p

  • tjac::Matrix{S} where S<:Real

    Temporary for transpose of Jacobian w.r.t y

  • Jx::Array{Matrix{S}, 1} where S<:Real

    Storage for vector of Jacobian w.r.t x

  • Jp::Array{Matrix{S}, 1} where S<:Real

    Storage for vector of Jacobian w.r.t p

  • result::DiffResults.MutableDiffResult{1, Vector{S}, Tuple{Matrix{S}}} where S<:Real

    Jacobian Result from DiffResults

  • cfg::ForwardDiff.JacobianConfig{Nothing, S, NY, Tuple{Array{ForwardDiff.Dual{Nothing, S, NY}, 1}, Array{ForwardDiff.Dual{Nothing, S, NY}, 1}}} where {S<:Real, NY}

    Jacobian Configuration for ForwardDiff

  • xtaylor::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, ForwardDiff.Dual{Nothing, S, NY}}, 1} where {N, S<:Real, NY}

    Store temporary STaylor1 vector for calculations

  • xaux::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, ForwardDiff.Dual{Nothing, S, NY}}, 1} where {N, S<:Real, NY}

    Store temporary STaylor1 vector for calculations

  • dx::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, ForwardDiff.Dual{Nothing, S, NY}}, 1} where {N, S<:Real, NY}

    Store temporary STaylor1 vector for calculations

  • taux::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, T}, 1} where {N, T<:Real}

  • t::Float64

  • vnxt::Vector{Int64}

    Intermediate storage to avoid allocations in Taylor coefficient calc

  • fnxt::Vector{Float64}

    Intermediate storage to avoid allocations in Taylor coefficient calc

DynamicBoundspODEsDiscrete.LohnersFunctorType

LohnersFunctor

A functor used in computing bounds and relaxations via Lohner's method. The implementation of the parametric Lohner's method described in the paper in (1) based on the non-parametric version given in (2).

  1. [Sahlodin, Ali M., and Benoit Chachuat. "Discretize-then-relax approach for

convex/concave relaxations of the solutions of parametric ODEs." Applied Numerical Mathematics 61.7 (2011): 803-820.](https://www.sciencedirect.com/science/article/abs/pii/S0168927411000316)

  1. [R.J. Lohner, Computation of guaranteed enclosures for the solutions of

ordinary initial and boundary value problems, in: J.R. Cash, I. Gladwell (Eds.), Computational Ordinary Differential Equations, vol. 1, Clarendon Press, 1992, pp. 425–436.](http://www.goldsztejn.com/old-papers/Lohner-1992.pdf)

DynamicBoundspODEsDiscrete.PICallbackType
mutable struct PICallback{FH, FJ} <: DynamicBoundspODEsDiscrete.AbstractIntervalCallback

Functor object d that computes h(xmid, P) and hj(X,P) in-place using X, P information stored in the fields when d() is run.

DynamicBoundspODEsDiscrete.StepParamsType

StepParams

LEPUS and Integration parameters.

  • atol::Float64

    Absolute error tolerance of integrator

  • rtol::Float64

    Relative error tolerance of integrator

  • hmin::Float64

    Minimum stepsize

  • repeat_limit::Int64

    Number of repetitions allowed for refinement

  • is_adaptive::Bool

    Indicates an adaptive stepsize is used

  • skip_step2::Bool

    Indicates the contractor step should be skipped

DynamicBoundspODEsDiscrete.StepResultType

StepResult{S}

Results passed to the next step.

  • xⱼ::Vector{Float64}

    nominal value of the state variables

  • Xⱼ::Vector

    relaxations/bounds of the state variables

  • A_Q::DynamicBoundspODEsDiscrete.FixedCircularBuffer{Matrix{Float64}}

    storage for parallelepid enclosure of xⱼ

  • A_inv::DynamicBoundspODEsDiscrete.FixedCircularBuffer{Matrix{Float64}}

  • Δ::DynamicBoundspODEsDiscrete.FixedCircularBuffer{Vector{S}} where S

    storage for parallelepid enclosure of xⱼ

  • predicted_hj::Float64

    predicted step size for next step

  • time::Float64

    new time

DynamicBoundspODEsDiscrete.TaylorFunctor!Type

TaylorFunctor!

A function g!(out, y) that perfoms a Taylor coefficient calculation. Provides preallocated storage. Evaluating this function out is a vector of length nx(s+1) where 1:(s+1) are the Taylor coefficients of the first component, (s+2):nx(s+1) are the Taylor coefficients of the second component, and so on. This may be constructed using TaylorFunctor!(g!, nx::Int, np::Int, k::Val{K}, t::T, q::Q) were type T should use type Q for internal computations. The order of the TaylorSeries is k, the right-hand side function is g!, nx is the number of state variables, np is the number of parameters.

  • g!::Function

    Right-hand side function for pODE which operates in place as g!(dx,x,p,t)

  • nx::Int64

    Dimensionality of x

  • np::Int64

    Dimensionality of p

  • k::Int64

    Order of TaylorSeries, that is the first k terms are used in the approximation and N = k+1 term is bounded

  • x::Vector{S} where S<:Real

    State variables x

  • p::Vector{S} where S<:Real

    Decision variables p

  • xtaylor::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, S}, 1} where {N, S<:Real}

    Store temporary STaylor1 vector for calculations

  • xaux::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, S}, 1} where {N, S<:Real}

    Store temporary STaylor1 vector for calculations

  • dx::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, S}, 1} where {N, S<:Real}

    Store temporary STaylor1 vector for calculations

  • taux::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, T}, 1} where {N, T<:Real}

  • vnxt::Vector{Int64}

  • fnxt::Vector{Float64}

DynamicBoundspODEsDiscrete.Wilhelm2019Type
mutable struct Wilhelm2019{T<:DynamicBoundspODEsDiscrete.Wilhelm2019Type, ICB1<:DynamicBoundspODEsDiscrete.PICallback, ICB2<:DynamicBoundspODEsDiscrete.PICallback, PRE, CTR<:DynamicBoundspODEsDiscrete.AbstractContractor, IC<:Function, F, Z, J, PRB, N, C, AMAT} <: DynamicBoundsBase.AbstractODERelaxIntegrator

An integrator that bounds the numerical solution of the pODEs system.

  • integrator_type::DynamicBoundspODEsDiscrete.Wilhelm2019Type

  • time::Vector{Float64}

  • steps::Int64

  • p::Vector{Float64}

  • pL::Vector{Float64}

  • pU::Vector{Float64}

  • nx::Int64

  • np::Int64

  • xL::Matrix{Float64}

  • xU::Matrix{Float64}

  • integrator_states::DynamicBoundsBase.IntegratorStates

  • evaluate_interval::Bool

  • evaluate_state::Bool

  • differentiable_flag::Bool

  • P::Vector{IntervalArithmetic.Interval{Float64}}

  • X::Matrix{IntervalArithmetic.Interval{Float64}}

  • X0P::Vector{IntervalArithmetic.Interval{Float64}}

  • pi_callback1::DynamicBoundspODEsDiscrete.PICallback

  • pi_callback2::DynamicBoundspODEsDiscrete.PICallback

  • pi_precond!::Any

  • pi_contractor::DynamicBoundspODEsDiscrete.AbstractContractor

  • inclusion_flag::Bool

  • exclusion_flag::Bool

  • extended_division_flag::Bool

  • ic::Function

  • h1::DynamicBoundspODEsDiscrete.CallbackH{Z, F, DynamicBoundspODEsDiscrete.ImpEuler} where {F, Z}

  • h2::DynamicBoundspODEsDiscrete.CallbackH{Z, F, T} where {T<:DynamicBoundspODEsDiscrete.Wilhelm2019Type, F, Z}

  • hj1::DynamicBoundspODEsDiscrete.CallbackHJ{J, DynamicBoundspODEsDiscrete.ImpEuler} where J

  • hj2::DynamicBoundspODEsDiscrete.CallbackHJ{J, T} where {T<:DynamicBoundspODEsDiscrete.Wilhelm2019Type, J}

  • mccallback1::McCormick.MCCallback

  • mccallback2::McCormick.MCCallback

  • IC_relax::Vector

  • state_relax::Matrix

  • var_relax::Vector

  • param::Array{Array{Vector{Z}, 1}, 1} where Z

  • kmax::Int64

  • calculate_local_sensitivity::Bool

  • local_problem_storage::Any

  • prob::Any

  • constant_state_bounds::Union{Nothing, DynamicBoundsBase.ConstantStateBounds}

  • relax_t_dict_flt::Dict{Float64, Int64}

  • relax_t_dict_indx::Dict{Int64, Int64}

DynamicBoundspODEsDiscrete.compute_X0!Method

compute_X0!(d::DiscretizeRelax)

Initializes the circular buffer that holds Δ with the out - mid(out) at index 1 and a zero vector at all other indices.

DynamicBoundspODEsDiscrete.existence_uniqueness!Method

existence_uniqueness!

Implements the adaptive higher-order enclosure approach detailed in Nedialkov's dissertation (Nedialko S. Nedialkov. Computing rigorous bounds on the solution of an initial value problem for an ordinary differential equation. 1999. Universisty of Toronto, PhD Dissertation, Algorithm 5.1, page 73-74). The arguments are s::ExistStorage{F,K,S,T}, params::StepParams, t::Float64, j::Int64.

DynamicBoundspODEsDiscrete.jacobian_taylor_coeffs!Method

jacobiantaylorcoeffs!

Computes the Jacobian of the Taylor coefficients w.r.t. y = (x,p) storing the output inplace to result. A JacobianConfig object without tag checking, cfg, is required input and is initialized from cfg = ForwardDiff.JacobianConfig(nothing, out, y). The JacTaylorFunctor! used for the evaluation is g and inputs are x and p.

DynamicBoundspODEsDiscrete.jetcoeffs!Method
jetcoeffs!(
    eqsdiff!,
    t::Number,
    x::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, U<:Number}, 1},
    xaux::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, U<:Number}, 1},
    dx::Array{DynamicBoundspODEsDiscrete.StaticTaylorSeries.STaylor1{N, U<:Number}, 1},
    order::Int64,
    params,
    vnxt::Vector{Int64},
    fnxt::Vector{Float64}
)

A variant of the jetcoeffs! function used in TaylorIntegration.jl (https://github.com/PerezHz/TaylorIntegration.jl/blob/master/src/explicitode.jl) which preallocates taux and updates taux coefficients to avoid further allocations.

The TaylorIntegration.jl package is licensed under the MIT "Expat" License: Copyright (c) 2016-2020: Jorge A. Perez and Luis Benet. Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions: The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

DynamicBoundspODEsDiscrete.set_JxJp!Method

set_JxJp!

Extracts the Jacobian of the Taylor coefficients w.r.t. x, Jx, and the Jacobian of the Taylor coefficients w.r.t. p, Jp, from result. The order of the Taylor series is s, the dimensionality of x is nx, the dimensionality of p is np, and tjac is preallocated storage for the transpose of the Jacobian w.r.t. y = (x,p).

DynamicBoundspODEsDiscrete.set_Δ!Method

set_Δ!

Initializes the circular buffer, Δ, that holds Δ_i with the out - mid(out) at index 1 and a zero vector at all other indices.

DynamicBoundspODEsDiscrete.αMethod

calc_alpha

Computes the stepsize for the adaptive step-routine via a golden section rootfinding method. The step size is rounded down.

DynamicBoundspODEsDiscrete.μ!Method

μ!(xⱼ,x̂ⱼ,η)

Used to compute the arguments of Jacobians (x̂ⱼ + η(xⱼ - x̂ⱼ)) used by the parametric Mean Value Theorem. The result is stored to out.

DynamicBoundspODEsDiscrete.ρ!Method

ρ!(out,p,p̂ⱼ,η)

Used to compute the arguments of Jacobians (p̂ⱼ + η(p - p̂ⱼ)) used by the parametric Mean Value Theorem. The result is stored to out.