Skip to content

Commit

Permalink
merge code for overapproximate of TaylorModelReachSet
Browse files Browse the repository at this point in the history
  • Loading branch information
schillic committed Jan 25, 2024
1 parent 296c76e commit bc0e41e
Showing 1 changed file with 23 additions and 67 deletions.
90 changes: 23 additions & 67 deletions src/ReachSets/TaylorModelReachSet.jl
Original file line number Diff line number Diff line change
Expand Up @@ -256,13 +256,12 @@ function overapproximate(R::TaylorModelReachSet{N}, ::ZonotopeEnclosure) where {
end

# tdom defined in the methods below is the same as Δt - t0, but the domain inclusion check
# in TM.evauate may fail, so we unpack the domain again here; also note that
# in TM.evaluate may fail, so we unpack the domain again here; also note that
# by construction the TMs in time are centered at zero

# overapproximate taylor model reachset with one hyperrectangle
function _overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Hyperrectangle};
Δt::TimeInterval=tspan(R), dom=symBox(dim(R))) where {N}

# overapproximate Taylor model reachset with a set
function overapproximate(R::TaylorModelReachSet, T::Type{<:LazySet};
Δt::TimeInterval=tspan(R), dom=symBox(dim(R)), kwargs...)
# dimension of the reachset
n = dim(R)

Expand All @@ -274,9 +273,6 @@ function _overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Hyperrectangle};
Δt tspan(R) || throw(ArgumentError("the given time range $Δt does not belong to the " *
"reach-set's time span, $(tspan(R))"))

# normalized time domain
tdom = domain(R)

# evaluate the Taylor model in time
X = set(R)
tdom = Δt - tstart(R) # normalize time (to TM-internal time)
Expand All @@ -295,13 +291,23 @@ function _overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Hyperrectangle};
# fX̂ is a TaylorModelN whose coefficients are floats
fX̂ = fp_rpa.(X̂)

# evaluate the spatial variables in the given domain
X̂b = IntervalBox([evaluate(fX̂[i], dom) for i in 1:n]...)
H = convert(Hyperrectangle, X̂b)
return ReachSet(H, Δt)
# evaluate the spatial variables in the given domain and overapproximate
Y = _overapproximate_TM(fX̂, T, dom; kwargs...)

return ReachSet(Y, Δt)
end

function _overapproximate_TM(fX̂, ::Type{<:Hyperrectangle}, dom; kwargs...)
# TODO outsource to LazySets
return convert(Hyperrectangle, IntervalBox([evaluate(fX̂[i], dom) for i in eachindex(fX̂)]...))
end

function _overapproximate_TM(fX̂, ::Type{<:Zonotope}, dom; kwargs...)
# overapproximate a Taylor model with a zonotope using LazySets
return overapproximate(fX̂, Zonotope; kwargs...)
end

LazySets.box_approximation(R::TaylorModelReachSet) = _overapproximate(R, Hyperrectangle)
LazySets.box_approximation(R::TaylorModelReachSet) = overapproximate(R, Hyperrectangle)

# overapproximate taylor model reachset with several hyperrectangles,
# by splitting in time (ntdiv), and/or space (with either nsdiv for uniform partitions
Expand All @@ -323,7 +329,7 @@ function overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Hyperrectangle};

# no splitting
if isnothing(partition) && nsdiv == 1 && ntdiv == 1
return _overapproximate(R, Hyperrectangle; Δt, dom)
return overapproximate(R, Hyperrectangle; Δt, dom)
end
@assert Δt == tspan(R) "time subdomain not implemented"
@assert dom == S "spatial subdomain not implemented"
Expand Down Expand Up @@ -378,44 +384,6 @@ function overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Hyperrectangle};
return R̂out
end

# overapproximate taylor model reachset with one zonotope
function overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Zonotope};
Δt::TimeInterval=tspan(R), dom=symBox(dim(R))) where {N}

# consistency checks
dim(R) == dim(dom) ||
throw(ArgumentError("the dimension of the reach-set should match the dimension of the domain, " *
"but they are $(dim(R)) and $(dim(dom)) respectively"))
dom symBox(dim(R)) || throw(ArgumentError("`dom` must be a subset of [-1, 1]^n"))
Δt tspan(R) || throw(ArgumentError("the given time range $Δt does not belong to the " *
"reach-set's time span, $(tspan(R))"))

# dimension of the reachset
n = dim(R)

# evaluate the Taylor model in time
X = set(R)
tdom = Δt - tstart(R) # normalize time (to TM-internal time)
tdom = tdom domain(R) # intersection handles round-off errors
# X_Δt is a vector of TaylorN (spatial variables) whose coefficients are intervals
X_Δt = evaluate(X, tdom)

# transform the domain if needed
X_Δt = _taylor_shift(X_Δt, dom)

# builds the associated taylor model for each coordinate j = 1...n
# X̂ is a TaylorModelN whose coefficients are intervals
= [TaylorModelN(X_Δt[j], zeroI, zeroBox(n), symBox(n)) for j in 1:n]

# compute floating point rigorous polynomial approximation
# fX̂ is a TaylorModelN whose coefficients are floats
fX̂ = fp_rpa.(X̂)

# overapproximate the Taylor model with a zonotope using LazySets
Z = overapproximate(fX̂, Zonotope)
return ReachSet(Z, Δt)
end

# overapproximate taylor model reachset with several zonotopes given the partition,
# that specifies the way in which the symmetric box [-1, 1]^n is split:
#
Expand Down Expand Up @@ -448,21 +416,9 @@ function overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Zonotope},
return ReachSet(ConvexHullArray(Z), Δt)
end

# evaluate at a given time and overapproximate the resulting set with a zonotope
function overapproximate(R::TaylorModelReachSet{N}, ::Type{<:Zonotope}, t::AbstractFloat;
remove_zero_generators=true) where {N}
@assert t tspan(R) "the given time point $t does not belong to the reach-set's time span, $(tspan(R))"

X = set(R)
tn = t - tstart(R) # normalize time (to TM-internal time)
tn = interval(tn) domain(R) # intersection handles round-off errors
X_Δt = evaluate(X, tn)
n = dim(R)
= [TaylorModelN(X_Δt[j], zeroI, zeroBox(n), symBox(n)) for j in 1:n]
fX̂ = fp_rpa.(X̂)
Zi = overapproximate(fX̂, Zonotope; remove_zero_generators=remove_zero_generators)
Δt = TimeInterval(t, t)
return ReachSet(Zi, Δt)
# evaluate at a given time point and overapproximate the resulting set
function overapproximate(R::TaylorModelReachSet, T::Type{<:LazySet}, t::AbstractFloat; kwargs...)
return overapproximate(R, T, TimeInterval(t, t); kwargs...)
end

# convert a hyperrectangular set to a taylor model reachset
Expand Down

0 comments on commit bc0e41e

Please sign in to comment.