-
Notifications
You must be signed in to change notification settings - Fork 32
/
ConvexHullArray.jl
244 lines (176 loc) · 5.73 KB
/
ConvexHullArray.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
export ConvexHullArray, CHArray,
array
"""
ConvexHullArray{N, S<:LazySet{N}} <: ConvexSet{N}
Type that represents the symbolic convex hull of a finite number of sets.
### Fields
- `array` -- array of sets
### Notes
The `EmptySet` is the neutral element for `ConvexHullArray`.
A `ConvexHullArray` is always convex.
### Examples
Convex hull of 100 two-dimensional balls whose centers follow a sinusoidal:
```jldoctest
julia> b = [Ball2([2*pi*i/100, sin(2*pi*i/100)], 0.05) for i in 1:100];
julia> c = ConvexHullArray(b);
```
"""
struct ConvexHullArray{N,S<:LazySet{N}} <: ConvexSet{N}
array::Vector{S}
end
isoperationtype(::Type{<:ConvexHullArray}) = true
isconvextype(::Type{<:ConvexHullArray}) = true
# constructor for an empty hull with optional size hint and numeric type
function ConvexHullArray(n::Int=0, N::Type=Float64)
a = Vector{LazySet{N}}()
sizehint!(a, n)
return ConvexHullArray(a)
end
# EmptySet is the neutral element for ConvexHullArray
@neutral(ConvexHullArray, EmptySet)
# Universe is the absorbing element for ConvexHullArray
@absorbing(ConvexHullArray, Universe)
# add functions connecting ConvexHull and ConvexHullArray
@declare_array_version(ConvexHull, ConvexHullArray)
"""
CHArray
Alias for `ConvexHullArray`.
"""
const CHArray = ConvexHullArray
"""
array(cha::ConvexHullArray)
Return the array of a convex hull of a finite number of sets.
### Input
- `cha` -- convex hull of a finite number of sets
### Output
The array of a convex hull of a finite number of sets.
"""
function array(cha::ConvexHullArray)
return cha.array
end
"""
dim(cha::ConvexHullArray)
Return the dimension of the convex hull of a finite number of sets.
### Input
- `cha` -- convex hull of a finite number of sets
### Output
The ambient dimension of the convex hull of a finite number of sets, or `0` if
there is no set in the array.
"""
function dim(cha::ConvexHullArray)
return length(cha.array) == 0 ? 0 : dim(cha.array[1])
end
"""
σ(d::AbstractVector, cha::ConvexHullArray)
Return a support vector of a convex hull of a finite number of sets in a given
direction.
### Input
- `d` -- direction
- `cha` -- convex hull of a finite number of sets
### Output
A support vector in the given direction.
"""
function σ(d::AbstractVector, cha::ConvexHullArray)
@assert !isempty(cha.array) "an empty convex hull is not allowed"
return _σ_union(d, array(cha))
end
"""
ρ(d::AbstractVector, cha::ConvexHullArray)
Evaluate the support function of a convex hull of a finite number of sets in a
given direction.
### Input
- `d` -- direction
- `cha` -- convex hull of a finite number of sets
### Output
The evaluation of the support function of the convex hull of a finite number of
sets in the given direction.
### Algorithm
This algorithm calculates the maximum over all ``ρ(d, X_i)``, where the
``X_1, …, X_k`` are the sets in the array of `cha`.
"""
function ρ(d::AbstractVector, cha::ConvexHullArray)
return maximum(ρ(d, Xi) for Xi in cha)
end
"""
isbounded(cha::ConvexHullArray)
Check whether a convex hull of a finite number of sets is bounded.
### Input
- `cha` -- convex hull of a finite number of sets
### Output
`true` iff all wrapped sets are bounded.
"""
function isbounded(cha::ConvexHullArray)
return all(isbounded, cha.array)
end
function isboundedtype(::Type{<:ConvexHullArray{N,S}}) where {N,S}
return isboundedtype(S)
end
"""
isempty(cha::ConvexHullArray)
Check whether a convex hull of a finite number of sets is empty.
### Input
- `cha` -- convex hull of a finite number of sets
### Output
`true` iff all wrapped sets are empty.
"""
function isempty(cha::ConvexHullArray)
return all(isempty, array(cha))
end
"""
vertices_list(cha::ConvexHullArray; [apply_convex_hull]::Bool=true,
[backend]=nothing, [prune]::Bool=apply_convex_hull)
Return a list of vertices of the convex hull of a finite number of sets.
### Input
- `cha` -- convex hull of a finite number of sets
- `apply_convex_hull` -- (optional, default: `true`) if `true`, post-process the
vertices using a convex-hull algorithm
- `backend` -- (optional, default: `nothing`) backend for computing
the convex hull (see argument `apply_convex_hull`)
- `prune` -- (optional, default: `apply_convex_hull`) alias for
`apply_convex_hull`
### Output
A list of vertices.
"""
function vertices_list(cha::ConvexHullArray;
apply_convex_hull::Bool=true,
backend=nothing,
prune::Bool=apply_convex_hull)
vlist = vcat([vertices_list(Xi) for Xi in cha]...)
if apply_convex_hull || prune
convex_hull!(vlist; backend=backend)
end
return vlist
end
# list of constraints of a convex-hull array of singletons
function constraints_list(X::ConvexHullArray{N,Singleton{N,VT}}) where {N,VT}
n = dim(X)
ST = n == 2 ? VPolygon : VPolytope
V = convert(ST, X)
return constraints_list(V)
end
# membership in a convex-hull array of singletons
function ∈(x::AbstractVector, X::ConvexHullArray)
n = length(x)
ST = n == 2 ? VPolygon : VPolytope
V = convert(ST, X)
return x ∈ V
end
function concretize(cha::ConvexHullArray)
a = array(cha)
@assert !isempty(a) "an empty convex hull is not allowed"
if all(is_polyhedral, a)
return _convex_hull_polytopes(cha)
end
X = cha
@inbounds for (i, Y) in enumerate(a)
if i == 1
X = concretize(Y)
else
X = convex_hull(X, concretize(Y))
end
end
return X
end
function translate(cha::ConvexHullArray, x::AbstractVector)
return ConvexHullArray([translate(X, x) for X in array(cha)])
end