TaylorModels.jl

TaylorModels.RTaylorModel1Type
RTaylorModel1{T,S}

Relative Taylor model in 1 variable, providing a rigurous polynomial approximation given by a Taylor polynomial pol (around x0) and a relative remainder rem for a function f(x) in one variable, valid in the interval dom. This corresponds to definition 2.3.2 of Mioara Joldes' thesis.

Fields:

  • pol: polynomial approximation, represented as TaylorSeries.Taylor1
  • rem: the interval bound
  • x0 : expansion point
  • dom: domain, interval over which the Taylor model is defined / valid

The approximation $f(x) = p(x) + \delta (x - x_0)^{n+1}$ is satisfied for all $x\in \mathcal{D}$; n is the order (degree) of the polynomial $p(x)=\sum_{i=0}^n p_i (x - x_0)^i$.

source
TaylorModels.TMSolType
TMSol{N,T,V1,V2,M}

Structure containing the solution of a validated integration.

Fields

`time :: AbstractVector{T}`  Vector containing the expansion time of the `TaylorModel1` solutions

`fp   :: AbstractVector{IntervalBox{N,T}}`  IntervalBox vector representing the flowpipe

`xTMv :: AbstractMatrix{TaylorModel1{TaylorN{T},T}}`  Matrix whose entry `xTMv[i,t]` represents
the `TaylorModel1` of the i-th dependent variable, obtained at time time[t].
source
TaylorModels.TaylorModel1Type
TaylorModel1{T,S}

Absolute Taylor model in 1 variable, providing a rigurous polynomial approximation given by a Taylor polynomial pol (around x0) and an absolute remainder rem for a function f(x) in one variable, valid in the interval dom. This corresponds to definition 2.1.3 of Mioara Joldes' thesis.

Fields:

  • pol: polynomial approximation, represented as TaylorSeries.Taylor1
  • rem: the interval bound
  • x0 : expansion point
  • dom: domain, interval over which the Taylor model is defined / valid

The approximation $f(x) = p(x) + \Delta$ is satisfied for all $x\in \mathcal{D}$ ($0\in \Delta$); n is the order (degree) of the polynomial $p(x)=\sum_{i=0}^n p_i (x - x_0)^i$.

source
TaylorModels._monot_bound_remainderMethod
_monot_bound_remainder(::Type{RTaylorModel1}, ::Val{true}, f::Function, polf::Taylor1, polfI::Taylor1, x0, I::Interval)

Computes the remainder exploiting monotonicity; see Prop 2.3.7 in Mioara Joldes' PhD thesis (pp 67).

source
TaylorModels._monot_bound_remainderMethod
_monot_bound_remainder(::Type{TaylorModel1}, ::Val{false}, f::Function, polf::Taylor1, polfI::Taylor1, x0, I::Interval)

Computes the remainder using Lagrange bound

source
TaylorModels._monot_bound_remainderMethod
_monot_bound_remainder(::Type{TaylorModel1}, ::Val{true}, f::Function, polf::Taylor1, polfI::Taylor1, x0, I::Interval)

Computes the remainder exploiting monotonicity; see Prop 2.2.1 in Mioara Joldes' PhD thesis (pp 52).

source
TaylorModels._picardMethod
picard(dx, x0, box)

Computes the picard (integral) operator for the initial condition x0. dx must be the rhs of the differential equation.

source
TaylorModels._rpaFunction
   _rpa(::Type{TaylorModel1}, f::Function, x0::Interval, I::Interval, _order::Integer)
   _rpa(::Type{RTaylorModel1}, f::Function, x0::Interval, I::Interval, _order::Integer)

Rigurous polynomial approximation (RPA) with absolute/relative remainder for the function f on the interval I, using a Taylor expansion around the interval x0 of order _order. The bound is computed by bound_remainder(@ref) exploiting monotonicity if possible, otherwise, it uses Lagrange bound.

source
TaylorModels._validate_step!Method
_validate_step(xTM1K, f!, dx, x0, params, t, box, dof; ε=1e-10, maxsteps=20, extrasteps=50)

Validate the Taylor Model solution for the current integration time step. This function implements the epsilon inflation algorithm proposed by Bünger with some custom adaptations.

Ref: Florian B"unger, A Taylor model toolbox for solving ODEs implemented in MATLAB/INTLAB, J. Comput. Appl. Math. 368, 112511, https://doi.org/10.1016/j.cam.2019.112511

source
TaylorModels.absorb_remainderMethod
absorb_remainder(a::TaylorModelN)

Returns a TaylorModelN, equivalent to a, such that the remainder is mostly absorbed in the constant and linear coefficients. The linear shift assumes that a is normalized to the IntervalBox(-1..1, Val(N)).

Ref: Xin Chen, Erika Abraham, and Sriram Sankaranarayanan, "Taylor Model Flowpipe Construction for Non-linear Hybrid Systems", in Real Time Systems Symposium (RTSS), pp. 183-192 (2012), IEEE Press.

source
TaylorModels.bound_integrationFunction
bound_integration(xTM::TaylorModel1{T,S}, δ)
bound_integration(xTM::Vector{TaylorModel1{T,S}}, δ)

Bound the remainder of the integration of a xTM::TaylorModel1, where δ is the domain used to bound the integration. The remainder corresponds to $δ * remainder(a) + a.pol[order] * δ^(order+1) / (order+1)$. This is tighter that the one used by Berz+Makino, which corresponds to $Δ = δ * remainder(a) + a.pol[order] * δ^(order+1)$.


bound_integration(xTM::RTaylorModel1{T,S}, δ)
bound_integration(xTM::Vector{RTaylorModel1{T,S}, δ)

Remainder bound for the integration of a xTM::RTaylorModel1, where δ is the domain used to bound the integration. The remainder corresponds to $Δ = δ * remainder(a)/(order+2) + getcoeff(polynomial(a), order)/(order+1)$.

source
TaylorModels.bound_remainderMethod

bound_remainder(::Type{RTaylorModel1}, f::Function, polf::Taylor1, polfI::Taylor1, x0::Interval, I::Interval)

Bound the relative remainder of the polynomial approximation of f given by the Taylor polynomial polf around x0 on the interval I. It requires an the interval extension polfI of a polynomial that approximates f for the whole interval I, in order to compute the Lagrange remainder.

If polfI[end] has a definite sign, then it is monotonic in the interval I, which is exploited; otherwise, the last coefficients bounds the relative remainder. This corresponds to Prop 2.3.7 in Mioara Joldes' PhD thesis (pp 67).

source
TaylorModels.bound_remainderMethod

bound_remainder(::Type{TaylorModel1}f::Function, polf::Taylor1, polfI::Taylor1, x0::Interval, I::Interval)

Bound the absolute remainder of the polynomial approximation of f given by the Taylor polynomial polf around x0 on the interval I. It requires the interval extension polfI of the polynomial that approximates f for the whole interval I, in order to compute the Lagrange remainder.

If polfI[end] has a definite sign, then it is monotonic in the intervals [I.lo, x0] and [x0.hi, I.hi], which is exploited; otherwise, it is used to compute the Lagrange remainder. This corresponds to Prop 2.2.1 in Mioara Joldes PhD thesis (pp 52).

source
TaylorModels.bound_taylor1Function
bound_taylor1(fT::TaylorModel1, I=domain(fT)::Interval)

Compute a tight polynomial bound for the Taylor model fT in the interval I, considering whether its derivative ftd has a definite sign.

source
TaylorModels.bound_taylor1Method
bound_taylor1(fT::Taylor1, I::Interval)

Compute a tight polynomial bound for the Taylor polynomial fT in the interval I.

Note: Algorithm 2.1.1 corresponds to evaluate(fT, I) or simply fT(I). This function uses the roots of the derivative offt` to obtain a tighter bound.

source
TaylorModels.bound_taylor1Method
bound_taylor1(fT::Taylor1, fTd::Taylor1, I::Interval)

Compute a tight polynomial bound for the Taylor polynomial fT in the interval I, considering whether its derivative ftd has a definite sign.

source
TaylorModels.fp_rpaFunction
fp_rpa(tm::TaylorModel1{Interval{T},T})
fp_rpa(tm::RTaylorModel1{Interval{T},T})

Convert a tm TaylorModel1 to a TaylorModel1 whose polynomial coefficients are Float64. The accumulated error is added to the remainder. The mid point of the expansion interval is preferentially rounded down if it is not an exactly representable value.

source
TaylorModels.initialize!Method
initialize!(X0::Vector{TaylorModel1{TaylorN{T}, T}}, orderQ, orderT, x, xI)
initialize!(X0::Vector{TaylorModel1{TaylorN{T}, T}}, orderQ, orderT, x)

Initialize the auxiliary integration variables assuming that the given vector of taylor models X0 is normalized to the domain [-1, 1]^n in space.

source
TaylorModels.initialize!Method
initialize!(X0::IntervalBox, orderQ, orderT, x, xI)
initialize!(X0::IntervalBox, orderQ, orderT, x)

Initialize the internal integration variables and normalize the given interval box to the domain [-1, 1]^n.

source
TaylorModels.iscontractiveMethod
iscontractive(Δ, Δx)

Checks if Δ .⊂ Δx is satisfied. If `Δ ⊆ Δx is satisfied, it returns true if all cases where == holds corresponds to the zero Interval.

source
TaylorModels.linear_dominated_bounderMethod
linear_dominated_bounder(fT::TaylorModel1, ϵ=1e-3::Float64, max_iter=5::Int)

Compute a tighter polynomial bound for the Taylor model fT by the linear dominated bounder algorithm. The linear dominated algorithm is applied until the bound of fT gets tighter than ϵ or the number of steps reachs max_iter. The returned bound corresponds to the improved polynomial bound with the remainder of the TaylorModel1 included.

source
TaylorModels.linear_dominated_bounderMethod
linear_dominated_bounder(fT::TaylorModelN, ϵ=1e-3::Float64, max_iter=5::Int)

Compute a tighter polynomial bound for the Taylor model fT by the linear dominated bounder algorithm. The linear dominated algorithm is applied until the bound of fT gets tighter than ϵ or the number of steps reachs max_iter. The returned bound corresponds to the improved polynomial bound with the remainder of the TaylorModelN included.

source
TaylorModels.mince_in_timeMethod
mince_in_time(sol::TMSol; var=0, timediv=1) --> ::Vector{Interval}

For var=0, this function divides the time domain of each entry of sol in timediv parts (timediv==1 is the initial domain), and returns the time intervals where the solution holds. This is useful for plotting or finding specific events. For var ≥ 1, this function evaluates the flowpipes at the split domain intervals, which is useful to decrease the overapproximations associated to the whole time domain.

source
TaylorModels.picard_iterationMethod
picard_iteration(f!, dx, xTM1K, params, t, x0, box)

Computes the picard (integral) operator for the set of equations f! and the initial conditionx0. This method returns the remainder of the resulting Taylor Model without the remainder of the initial condition.

source
TaylorModels.picard_iterationMethod
picard_iteration(f!, dx, xTM1K, params, t, x0, box, ::Val{true})

Computes the picard (integral) operator for the set of equations f! and the initial conditionx0. The Val parameter enables the selection of the desired method. This one returns the remainder of the resulting Taylor Model with the remainder of the initial condition included.

source
TaylorModels.quadratic_fast_bounderMethod
quadratic_fast_bounder(fT::TaylorModel1)

Compute a tighter polynomial bound by the quadratic fast bounder. The returned bound corresponds to the "improved" polynomial bound with the remainder of the TaylorModel1 included. This "improved" bound can be one of the following two: 1) An improved bound: if the domain of fT has a local minimizer, then an improved bound is returned. 2) Original TaylorModel bound: if the local minimizer criteria is not satisfied, then the original bound of fT is returned.

This algorithm is a slight modification to the Makino & Berz algorithm. For this algorithm the linear part is bounded by solving a simple set of linear equations (compared to the iterative procedure by Makino & Berz).

source
TaylorModels.quadratic_fast_bounderMethod
quadratic_fast_bounder(fT::TaylorModelN)

Compute a tighter polynomial bound by the quadratic fast bounder. The returned bound corresponds to the "improved" polynomial bound with the remainder of the TaylorModelN included. This "improved" bound can be one of the following two: 1) An improved bound: if the domain of fT has a local minimizer, then an improved bound is returned. 2) Original TaylorModel bound: if the local minimizer criteria is not satisfied, then the original bound of fT is returned.

This algorithm is a slight modification to the Makino & Berz algorithm. For this algorithm the linear part is bounded by solving a simple set of linear equations (compared to the iterative procedure by Makino & Berz).

source
TaylorModels.remainder_taylorstep!Method
remainder_taylorstep!(f!, t, x, dx, xI, dxI, δI, δt, params)

Returns a remainder for the integration step for the dependent variables (x) checking that the solution satisfies the criteria for existence and uniqueness.

source
TaylorModels.rpaFunction

rpa(g::Function, tmf::TaylorModel1) rpa(g::Function, tmf::RTaylorModel1) rpa(g::Function, tmf::TaylorModelN)

Rigurous polynomial approximation (RPA) for the function g using the Taylor Model with absolute/relative remainder tmf. The bound is computed exploiting monotonicity if possible, otherwise, it uses Lagrange bound.

source
TaylorModels.shrink_wrapping!Method
shrink_wrapping!(xTMN::TaylorModelN)

Returns a modified inplace xTMN, which has absorbed the remainder by the modified shrink-wrapping method of Florian Bünger. The domain of xTMN is the normalized interval box [-1,1]^N.

Ref: Florian B"unger, Shrink wrapping for Taylor models revisited, Numer Algor 78:1001–1017 (2018), https://doi.org/10.1007/s11075-017-0410-1

source
TaylorSeries.integrateFunction
integrate(a::TM{T,S}, c0)
integrate(a::TM{TaylorN{T},S}, c0, cc0)

Integrates the one-variable Taylor Model (TaylorModel1 or RTaylorModel1) with respect to the independent variable. c0 is the integration constant; if omitted it is taken as zero. When the coefficients of a are TaylorN variables, the domain is specified by cc0::IntervalBox.


integrate(fT, which)

Integrates a fT::TaylorModelN with respect to which variable. The returned TaylorModelN corresponds to the Taylor Model of the definite integral ∫f(x) - ∫f(expansion_point).

source