-
Notifications
You must be signed in to change notification settings - Fork 32
/
ConvexHull.jl
222 lines (156 loc) · 4.64 KB
/
ConvexHull.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
export ConvexHull, CH,
convex_hull!,
ConvexHull!,
swap
"""
ConvexHull{N, S1<:LazySet{N}, S2<:LazySet{N}} <: ConvexSet{N}
Type that represents the convex hull of the union of two sets, i.e., the set
```math
Z = \\{z ∈ ℝ^n : z = λx + (1-λ)y,\\qquad x ∈ X, y ∈ Y,λ ∈ [0, 1] \\}.
```
### Fields
- `X` -- set
- `Y` -- set
### Notes
The `EmptySet` is the neutral element for `ConvexHull`.
This type is always convex.
### Examples
The convex hull of two 100-dimensional Euclidean balls:
```jldoctest
julia> b1, b2 = Ball2(zeros(100), 0.1), Ball2(4*ones(100), 0.2);
julia> c = ConvexHull(b1, b2);
julia> typeof(c)
ConvexHull{Float64, Ball2{Float64, Vector{Float64}}, Ball2{Float64, Vector{Float64}}}
```
"""
struct ConvexHull{N,S1<:LazySet{N},S2<:LazySet{N}} <: ConvexSet{N}
X::S1
Y::S2
# default constructor with dimension check
function ConvexHull(X::LazySet{N}, Y::LazySet{N}) where {N}
@assert dim(X) == dim(Y) "sets in a convex hull must have the same " *
"dimension"
return new{N,typeof(X),typeof(Y)}(X, Y)
end
end
isoperationtype(::Type{<:ConvexHull}) = true
concrete_function(::Type{<:ConvexHull}) = convex_hull
isconvextype(::Type{<:ConvexHull}) = true
is_polyhedral(ch::ConvexHull) = is_polyhedral(ch.X) && is_polyhedral(ch.Y)
# EmptySet is the neutral element for ConvexHull
@neutral(ConvexHull, EmptySet)
# Universe is the absorbing element for ConvexHull
@absorbing(ConvexHull, Universe)
# interface for binary set operations
Base.first(ch::ConvexHull) = ch.X
second(ch::ConvexHull) = ch.Y
@declare_binary_operation(ConvexHull)
"""
CH
Alias for `ConvexHull`.
"""
const CH = ConvexHull
"""
swap(ch::ConvexHull)
Return a new `ConvexHull` object with the arguments swapped.
### Input
- `ch` -- convex hull of two sets
### Output
A new `ConvexHull` object with the arguments swapped.
"""
function swap(ch::ConvexHull)
return ConvexHull(ch.Y, ch.X)
end
"""
dim(ch::ConvexHull)
Return the dimension of a convex hull of two sets.
### Input
- `ch` -- convex hull of two sets
### Output
The ambient dimension of the convex hull of two sets.
"""
function dim(ch::ConvexHull)
return dim(ch.X)
end
"""
σ(d::AbstractVector, ch::ConvexHull)
Return a support vector of the convex hull of two sets in a given direction.
### Input
- `d` -- direction
- `ch` -- convex hull of two sets
### Output
A support vector of the convex hull in the given direction.
"""
function σ(d::AbstractVector, ch::ConvexHull)
σ1 = σ(d, ch.X)
σ2 = σ(d, ch.Y)
ρ1 = dot(d, σ1)
ρ2 = dot(d, σ2)
return ρ1 >= ρ2 ? σ1 : σ2
end
"""
ρ(d::AbstractVector, ch::ConvexHull)
Evaluate the support function of the convex hull of two sets in a given
direction.
### Input
- `d` -- direction
- `ch` -- convex hull of two sets
### Output
The evaluation of the support function of the convex hull in the given
direction.
"""
function ρ(d::AbstractVector, ch::ConvexHull)
return max(ρ(d, ch.X), ρ(d, ch.Y))
end
"""
isbounded(ch::ConvexHull)
Check whether the convex hull of two sets is bounded.
### Input
- `ch` -- convex hull of two sets
### Output
`true` iff both wrapped sets are bounded.
"""
function isbounded(ch::ConvexHull)
return isbounded(ch.X) && isbounded(ch.Y)
end
function isboundedtype(::Type{<:ConvexHull{N,S1,S2}}) where {N,S1,S2}
return isboundedtype(S1) && isboundedtype(S2)
end
"""
isempty(ch::ConvexHull)
Check whether the convex hull of two sets is empty.
### Input
- `ch` -- convex hull
### Output
`true` iff both wrapped sets are empty.
"""
function isempty(ch::ConvexHull)
return isempty(ch.X) && isempty(ch.Y)
end
"""
vertices_list(ch::ConvexHull; [apply_convex_hull]::Bool=true,
[backend]=nothing)
Return a list of vertices of the convex hull of two sets.
### Input
- `ch` -- convex hull of two 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`)
### Output
A list of vertices.
"""
function vertices_list(ch::ConvexHull;
apply_convex_hull::Bool=true,
backend=nothing)
vlist = vcat(vertices_list(ch.X), vertices_list(ch.Y))
if apply_convex_hull
convex_hull!(vlist; backend=backend)
end
return vlist
end
function translate(ch::ConvexHull, x::AbstractVector)
X = translate(first(ch), x)
Y = translate(second(ch), x)
return ConvexHull(X, Y)
end