Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Overload TM functions for Taylor model reachsets #190

Merged
merged 5 commits into from
Oct 24, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
23 changes: 20 additions & 3 deletions src/Flowpipes/reachsets.jl
Original file line number Diff line number Diff line change
Expand Up @@ -590,9 +590,10 @@ using TaylorModels: TaylorModel1, TaylorN
"""
TaylorModelReachSet{N} <: AbstractTaylorModelReachSet{N}

Taylor model reach-set represented as a of vector taylor models in one variable
(the "time" variable) whose coefficients are multivariate polynomials
(the "space" variables).
Taylor model reach-set represented as a vector of taylor models in one variable
(namely, the "time" variable) whose coefficients are multivariate polynomials
(namely in the "space" variables). It is assumed that the time domain is the same
for all components.

### Notes

Expand All @@ -615,6 +616,22 @@ end
@inline dim(R::TaylorModelReachSet) = get_numvars()
@inline vars(R::TaylorModelReachSet) = Tuple(Base.OneTo(length(R.X)),)

# overload getter functions for the taylor model
# we assume that the first element is representative
domain(R::TaylorModelReachSet) = domain(first(R.X))
remainder(R::TaylorModelReachSet) = remainder.(R.X)
polynomial(R::TaylorModelReachSet) = polynomial.(R.X)
get_order(R::TaylorModelReachSet) = get_order.(R.X)
expansion_point(R::TaylorModelReachSet) = [Xi.x0 for Xi in R.X]

# useful constants
@inline zeroBox(m) = IntervalBox(zeroI, m)
@inline unitBox(m) = IntervalBox(IA.Interval(0.0, 1.0), m)
@inline symBox(n::Integer) = IntervalBox(symI, n)
const zeroI = IA.Interval(0.0) # TODO use number type
const oneI = IA.Interval(1.0)
const symI = IA.Interval(-1.0, 1.0)

function shift(R::TaylorModelReachSet, t0::Number)
return TaylorModelReachSet(set(R), tspan(R) + t0)
end
Expand Down
6 changes: 6 additions & 0 deletions src/Initialization/exports.jl
Original file line number Diff line number Diff line change
Expand Up @@ -46,6 +46,12 @@ export
sup_func, # TODO keep?
setrep,
rsetrep,
reset_map,
guard,
source_invariant,
target_invariant,
# getter functions for Taylor model reach-sets
domain, remainder, polynomial, get_order, expansion_point,
numrsets,

# Concrete operations
Expand Down
3 changes: 3 additions & 0 deletions src/Initialization/init.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,9 @@ const TimeInterval = IA.Interval{Float64}
import TaylorModels
const TM = TaylorModels

# method extensions for Taylor model reach-sets
import TaylorModels: domain, remainder, polynomial, get_order

@inline function _isapprox(Δt::TimeInterval, Δs::TimeInterval)
return (inf(Δt) ≈ inf(Δs)) && (sup(Δt) ≈ sup(Δs))
end
Expand Down
16 changes: 12 additions & 4 deletions test/algorithms/TMJets.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,13 +34,21 @@
end

@testset "TMJets algorithm: linear IVPs" begin
prob, tspan = exponential_1d()
sol = solve(prob, tspan=tspan, TMJets())
prob, dt = exponential_1d()
sol = solve(prob, tspan=dt, TMJets())
@test sol.alg isa TMJets

# getter functions for a taylor model reach-set
R = sol[1]
@test domain(R) == tspan(R)
@test diam(remainder(R)[1]) < 1e-13
@test get_order(R) == [8]
@test polynomial(R) isa Vector{Taylor1{TaylorN{Float64}}}
@test expansion_point(R) ≈ [IA.Interval(0.0)]

# test intersection with invariant
prob, tspan = exponential_1d(invariant=HalfSpace([-1.0], -0.3)) # x >= 0.3
sol_inv = solve(prob, tspan=tspan, TMJets())
prob, dt = exponential_1d(invariant=HalfSpace([-1.0], -0.3)) # x >= 0.3
sol_inv = solve(prob, tspan=dt, TMJets())
@test [0.3] ∈ overapproximate(sol_inv[end], Zonotope)
m = length(sol_inv)
# check that the following reach-set escapes the invariant
Expand Down
2 changes: 2 additions & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@ using ReachabilityAnalysis: _isapprox, setrep, rsetrep,
import IntervalArithmetic
const IA = IntervalArithmetic

const IA = IntervalArithmetic

# load test models
include("models/exponential1D.jl")
include("models/motor.jl")
Expand Down