A Tour of LazySets

LazySets is a library for set-based computations in Euclidean space. The library offers both concrete and lazy set representations, where the latter stands for the ability to delay set computations until they are needed. The choice of the programming language Julia and the accompanying documentation allow researchers to easily translate set-based algorithms from mathematics to software in a platform-independent way, while achieving runtime performance that is comparable to statically compiled languages. Combining lazy operations in high dimensions and explicit computations in low dimensions, LazySets can be applied to solve complex, large-scale problems. The library can handle the most common operations between sets, including Minkowski sum, Cartesian product, convex hull, intersection and union.

The purpose of this section is to give a high-level overview of the library and several examples to illustrate basic functionality. Further sections delve into more advanced uses.

Creating sets

Sets can be created using the constructor for each set type. Here we define the two-dimensional half-space

\[H = \{ (x, y) ∈ ℝ^2 : x ≤ y\} = \{ (x, y) ∈ ℝ^2 : x - y ≤ 0\}.\]

using LazySets

H = HalfSpace([1.0, -1.0], 0.0)
HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], 0.0)

The first argument is the normal vector to the half-space and the second argument the displacement, written in the standard form $a⋅x ≤ b$. We can define the same set using symbolic variables with the help of the Julia package Symbolics.jl:

using Symbolics

var = @variables x y

H′ = HalfSpace(x ≤ y, var)  # for `′`, type H\prime<tab>
HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0)

We can use LazySets to check that H and H′ represent the same sets (mathematically the set equivalence $H = H'$ corresponds to a double inclusion check $H ⊆ H' ∧ H' ⊆ H$):

isequivalent(H, H′)
true

High-dimensional sets can be conveniently defined using array notation. Here we create the nine-dimensional half-space

\[H9 = \left\{ x ∈ ℝ^{9} : ∑_{i=1}^{9} i x_i ≤ 10 \right\}.\]

var = @variables x[1:9]

expr = sum(i*x[i] for i in 1:9) ≤ 10.0
H9 = HalfSpace(expr, var)
HalfSpace{Float64, Vector{Float64}}([1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0], 10.0)

Next we ask LazySets to create two intervals $[0, 1]$ and $[2, 3]$ and then take their Cartesian product.

A = Interval(0, 1)
B = Interval(2, 3)

C = A × B
CartesianProduct{Float64, Interval{Float64}, Interval{Float64}}(Interval{Float64}([0, 1]), Interval{Float64}([2, 3]))

In LazySets, intervals are implemented as thin wrappers around those available in IntervalArithmetic.jl. Use the convert function to represent the set as a Hyperrectangle, which is a dedicated set type that represents a hyperrectangle with a center and a radius vector:

H = convert(Hyperrectangle, C)
Hyperrectangle{Float64, Vector{Float64}, Vector{Float64}}([0.5, 2.5], [0.5, 0.5])

IntervalArithmetic.jl represents multi-dimensional intervals with the type IntervalBox, to which we can also convert:

using IntervalArithmetic

Hbox = convert(IntervalBox, H)

typeof(Hbox)
IntervalArithmetic.IntervalBox{2, Float64}

Hyperrectangles are a special subclass of zonotopes.

Z = convert(Zonotope, C)
Zonotope{Float64, Vector{Float64}, Matrix{Float64}}([0.5, 2.5], [0.5 0.0; 0.0 0.5])

There are many other set representations available in this library. They are introduced in further sections below.

Operating with sets

The result of intersecting two or more half-spaces is a polyhedron:

var = @variables x y

H1 = HalfSpace(y ≥ x, var)
H2 = HalfSpace(y ≥ -x, var)

P = H1 ∩ H2  # for `∩`, type `\cap<tab>`
HPolyhedron{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0)])

Naturally, LazySets can plot the resulting set using the Plots.jl library:

using Plots

plot(H1, lab="H1")
plot!(H2, lab="H2")
plot!(P, lab="P = H1 ∩ H2", linewidth=2.0, linestyle=:dash)
Example block output

A polyhedron consists of a finite intersection of half-spaces, which are also called linear constraints:

constraints_list(P)
2-element Vector{HalfSpace{Float64, Vector{Float64}}}:
 HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0)
 HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0)

The intersection of $P$ with the triangle $R$ with vertices

$\{(-1.0, -1.0), (1.0, -1.0), (0.0, 1.0)\}$ can be constructed with:

R = VPolygon([[-1.0, -1.0], [1.0, -1.0], [0.0, 1.0]])

Q = P ∩ R
Intersection{Float64, HPolyhedron{Float64, Vector{Float64}}, VPolygon{Float64, Vector{Float64}}}(HPolyhedron{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0)]), VPolygon{Float64, Vector{Float64}}([[-1.0, -1.0], [1.0, -1.0], [0.0, 1.0]]), LazySets.IntersectionCache(-1))
plot!(Q, lab="Q = P ∩ R")
plot!(R, lab="R", linestyle=:dashdot, alpha=.2)
Example block output

The resulting set Q is a lazy (binary) intersection: by design, all operations in LazySets are lazy by default. Concrete operations can be performed by calling corresponding lower-case functions:

Q = intersection(P, R)
HPolytope{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0), HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0), HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)])

Alternatively, we can use the concretize function to transform a lazy set representation into a suitable concrete set representation.

concretize(P ∩ R)
HPolytope{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0), HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0), HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)])

The result is an HPolytope (instead of an HPolyhedron), since polytopes are special cases of bounded polyhedra.

isbounded(P)
false
isbounded(Q)
true

LazySets can compute various useful operations. For example, it can compute the area of a two-dimensional bounded set:

area(Q)
0.3333333333333333

When using rational numbers instead of floating-point numbers, we get an exact value for the area:

area(rationalize(Q))
1//3

The vertices can be retrieved using vertices_list,

vertices_list(Q)
4-element Vector{Vector{Float64}}:
 [0.0, 1.0]
 [-0.3333333333333333, 0.3333333333333333]
 [-0.0, 0.0]
 [0.3333333333333333, 0.3333333333333333]

or in iterator-fashion with vertices:

vertices(Q)
ReachabilityBase.Iteration.VectorIterator{Vector{Vector{Float64}}}([[0.0, 1.0], [-0.3333333333333333, 0.3333333333333333], [-0.0, 0.0], [0.3333333333333333, 0.3333333333333333]])

The constraint representation of $Q$ is obtained similarly.

constraints_list(Q)
4-element Vector{HalfSpace{Float64, Vector{Float64}}}:
 HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0)
 HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0)
 HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0)
 HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)

The constraint iterator is called with constraints.

constraints(Q)
ReachabilityBase.Iteration.VectorIterator{Vector{HalfSpace{Float64, Vector{Float64}}}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0), HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0), HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)])

The Minkowski sum of two sets can be computed using (type \oplus<tab>) (or alternatively just + for convenience):

# ball in the 2-norm of given center and radius
B = Ball2([0.8, 1.0], 0.3)

E = Q ⊕ B
MinkowskiSum{Float64, HPolytope{Float64, Vector{Float64}}, Ball2{Float64, Vector{Float64}}}(HPolytope{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0), HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0), HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)]), Ball2{Float64, Vector{Float64}}([0.8, 1.0], 0.3))

Approximating and plotting such sets relies on techniques based on the support function, discussed at length in other sections of this manual.

plot(E, lab="E = Q ⊕ B", legend=:bottomright)
plot!(Q, lab="Q")
plot!(B, lab="B", ratio=1., xlims=(-0.5, 1.5))
Example block output

Linear and affine maps are defined combining the * and + operators.

M = [0 1; 1 0.]
b = [0.5, -0.5]

X = M * E + b
AffineMap{Float64, MinkowskiSum{Float64, HPolytope{Float64, Vector{Float64}}, Ball2{Float64, Vector{Float64}}}, Float64, Matrix{Float64}, Vector{Float64}}([0.0 1.0; 1.0 0.0], MinkowskiSum{Float64, HPolytope{Float64, Vector{Float64}}, Ball2{Float64, Vector{Float64}}}(HPolytope{Float64, Vector{Float64}}(HalfSpace{Float64, Vector{Float64}}[HalfSpace{Float64, Vector{Float64}}([1.0, -1.0], -0.0), HalfSpace{Float64, Vector{Float64}}([-1.0, -1.0], 0.0), HalfSpace{Float64, Vector{Float64}}([2.0, 1.0], 1.0), HalfSpace{Float64, Vector{Float64}}([-2.0, 1.0], 1.0)]), Ball2{Float64, Vector{Float64}}([0.8, 1.0], 0.3)), [0.5, -0.5])

The set $X$, and in general, the composition of sets and operations, is again a set (in this case an operation).

isoperation(X)
true

Therefore, all set operations (and plotting) still work in the same way.

plot!(X, lab="X = M*E + b", legend=:topright)
Example block output

Convex hulls can be computed with convex_hull (concrete), ConvexHull (lazy, binary), or ConvexHullArray (lazy, n-ary).

C = ConvexHullArray([E, Q, X])

plot!(C, c=:grey, lab="C = CH(E, Q, X)")
Example block output

The Cartesian product of sets is defined with × (type \times<tab>) (or alternatively just * for convenience).

X = Interval(0, 1) × Ball1([1.0, 2.0, 3.0], 1.0) × Interval(2, 4)
CartesianProduct{Float64, CartesianProduct{Float64, Interval{Float64}, Ball1{Float64, Vector{Float64}}}, Interval{Float64}}(CartesianProduct{Float64, Interval{Float64}, Ball1{Float64, Vector{Float64}}}(Interval{Float64}([0, 1]), Ball1{Float64, Vector{Float64}}([1.0, 2.0, 3.0], 1.0)), Interval{Float64}([2, 4]))

Set unions are defined with (type \cup<tab>).

Y = Ball2([1.0, 2.0, 3.0], 1.0) × LineSegment([0.0, 0.0], [1.0, 1.0])

U = X ∪ Y
UnionSet{Float64, CartesianProduct{Float64, CartesianProduct{Float64, Interval{Float64}, Ball1{Float64, Vector{Float64}}}, Interval{Float64}}, CartesianProduct{Float64, Ball2{Float64, Vector{Float64}}, LineSegment{Float64, Vector{Float64}}}}(CartesianProduct{Float64, CartesianProduct{Float64, Interval{Float64}, Ball1{Float64, Vector{Float64}}}, Interval{Float64}}(CartesianProduct{Float64, Interval{Float64}, Ball1{Float64, Vector{Float64}}}(Interval{Float64}([0, 1]), Ball1{Float64, Vector{Float64}}([1.0, 2.0, 3.0], 1.0)), Interval{Float64}([2, 4])), CartesianProduct{Float64, Ball2{Float64, Vector{Float64}}, LineSegment{Float64, Vector{Float64}}}(Ball2{Float64, Vector{Float64}}([1.0, 2.0, 3.0], 1.0), LineSegment{Float64, Vector{Float64}}([0.0, 0.0], [1.0, 1.0])))

There are other set operations not mentioned in this section. See the remaining sections of this manual for further examples.

Exploring the type hierarchy

Every set type in the library inherits from the parametric abstract type LazySet{N}, where N is a parameter for the numeric type (typically, double-precision floating point numbers, Float64). This way one can easily choose between, e.g., floating point (Float64) and exact (Rational) precision with no additional performance penalty: At runtime, Julia uses multiple dispatch on N and JIT-compiles into type-specific code.

Since the LazySets type hierarchy is rather involved, a way to visualize it is to use AbstractTrees.jl. The package provides several utilities for working with tree-like data structures. The function print_tree provides a very detailed answer.

using AbstractTrees

AbstractTrees.children(x::Type) = LazySets.subtypes(x, false)

print_tree(LazySet)
LazySet
├─ AbstractAffineMap
│  ├─ AffineMap
│  ├─ ExponentialMap
│  ├─ ExponentialProjectionMap
│  ├─ InverseLinearMap
│  ├─ LinearMap
│  ├─ ResetMap
│  └─ Translation
├─ AbstractPolynomialZonotope
│  ├─ DensePolynomialZonotope
│  ├─ SimpleSparsePolynomialZonotope
│  └─ SparsePolynomialZonotope
├─ Bloating
├─ CachedMinkowskiSumArray
├─ CartesianProduct
├─ CartesianProductArray
├─ Complement
├─ ConvexSet
│  ├─ AbstractCentrallySymmetric
│  │  ├─ AbstractBallp
│  │  │  ├─ Ball2
│  │  │  └─ Ballp
│  │  └─ Ellipsoid
│  ├─ AbstractPolyhedron
│  │  ├─ AbstractPolytope
│  │  │  ├─ AbstractCentrallySymmetricPolytope
│  │  │  │  ├─ AbstractZonotope
│  │  │  │  │  ⋮
│  │  │  │  │
│  │  │  │  └─ Ball1
│  │  │  ├─ AbstractPolygon
│  │  │  │  ├─ AbstractHPolygon
│  │  │  │  │  ⋮
│  │  │  │  │
│  │  │  │  └─ VPolygon
│  │  │  ├─ HPolytope
│  │  │  ├─ Tetrahedron
│  │  │  └─ VPolytope
│  │  ├─ HPolyhedron
│  │  ├─ HalfSpace
│  │  ├─ Hyperplane
│  │  ├─ Line
│  │  ├─ Line2D
│  │  ├─ Star
│  │  └─ Universe
│  ├─ ConvexHull
│  ├─ ConvexHullArray
│  └─ EmptySet
├─ Intersection
├─ IntersectionArray
├─ MinkowskiSum
├─ MinkowskiSumArray
├─ Polygon
├─ QuadraticMap
├─ Rectification
├─ UnionSet
└─ UnionSetArray

Since the list does not fit into the default size, some types are hidden. LazySets has several set interfaces whose names start with Abstract. Visualizing the subtypes of a specific set interface can be done similarly. Consider the class of zonotopic sets, AbstractZonotope, which are those that can be represented as

\[Z = \left\{ x ∈ ℝ^n : x = c + ∑_{i=1}^p ξ_i g_i,~~ ξ_i ∈ [-1, 1]~~ ∀ i = 1,…, p \right\},\]

where $c ∈ ℝ^n$ is called the center and the vectors $\{g_i\}_{i=1}^p$, $g_i ∈ ℝ^n$, are called the generators.

print_tree(AbstractZonotope)
AbstractZonotope
├─ AbstractHyperrectangle
│  ├─ AbstractSingleton
│  │  ├─ Singleton
│  │  └─ ZeroSet
│  ├─ BallInf
│  ├─ Hyperrectangle
│  ├─ Interval
│  └─ SymmetricIntervalHull
├─ HParallelotope
├─ LineSegment
├─ RotatedHyperrectangle
└─ Zonotope

All the sets belonging to the abstract zonotope interface have in common that they can be characterized by a center and a generator matrix (equivalently, as the finite Minkowski sum of line segments, or as the image of a unit infinity-norm ball in $ℝ^n$ by an affine transformation). Hence, new set types implementing a few interface functions can make use of all the available functionality that are implemented by generic algorithms working at the interface level.

Last but not least, we remark that one of the key design choices of LazySets (which, admittedly, may be confusing at first) is that both set representations and set operations subtype LazySet. In other words, an "operation between sets" and "a set" are on the same footing in terms of belonging to the same type hierarchy. One of the advantages of these two mathematical concepts to be "just types" is to conveniently compose representations and operations (existing and user-created ones).

An overview of only the set operations can be obtained like so:

filter(isoperationtype, LazySets.subtypes(LazySet, true))
23-element Vector{Type}:
 AffineMap
 Bloating
 CachedMinkowskiSumArray
 CartesianProduct
 CartesianProductArray
 Complement
 ConvexHull
 ConvexHullArray
 ExponentialMap
 ExponentialProjectionMap
 ⋮
 MinkowskiSum
 MinkowskiSumArray
 QuadraticMap
 Rectification
 ResetMap
 SymmetricIntervalHull
 Translation
 UnionSet
 UnionSetArray

The total amount of available representations is larger (if we had used subtypes(LazySet), that would only show the direct subtypes of LazySet).

length(LazySets.subtypes(LazySet, true))
54

Random sampling and splitting

Another useful operation is to perform random sampling from a given set.

# returns a vector of samples (vectors)
S = sample(Q, 300, include_vertices=true)

plot(Q, lab="Q")

# transform into a vector of singletons (for plotting)
plot!(Singleton.(S), lab="S", c=:magenta)
Example block output

Set membership is computed with the operator (type \in<tab>):

all(s -> s ∈ Q, S)
true

Below we show a straightforward way to enclose Q with smaller, simpler sets (hyperrectangles in this case).

# compute an enclosing hyperrectangle
Qbox = overapproximate(Q, Hyperrectangle)

# split it into 6 pieces horizontally and 10 pieces vertically
Bs = split(Qbox, [6, 10])

# find all the boxes that intersect Q
idx = findall(Bi -> !isdisjoint(Bi, Q), Bs)

plot(Bs[idx])
plot!(Q, lab="Q")
plot!(Singleton.(S), lab="S", c=:magenta)
Example block output

With a bit of extra work we can also split the set into non axis-aligned parallelotopes:

# compute an enclosing zonotope
Z = overapproximate(Q, Zonotope, OctDirections)

# reduce to order 1
P = overapproximate(Z, HParallelotope)

# split along the first (resp. second) generator three (resp. four) times
Zs = split(convert(Zonotope, P), [1, 2], [3, 4])

# find all the zonotopes that intersect Q
idx = findall(Zi -> !isdisjoint(Zi, Q), Zs)

plot(Zs[idx])
plot!(Q, lab="Q")
plot!(Singleton.(S), lab="S", c=:magenta)
Example block output

Polyhedral approximations

In this library we mainly consider representations of closed convex sets in the usual sense from convex geometry: A set is closed if it contains all its limit points. A set $S$ is convex if for any $m$ points $v_j ∈ S$ and $m$ non-negative numbers $λ_j$ that sum up to $1$ we have that $∑_{j=1}^m λ_j v_j ∈ S$ as well. Alternatively, a closed convex set is an intersection of (possibly infinitely many) closed half-spaces.

One of the central functions are support_function(d, X) (with the alias ρ(d, X) (type \rho<tab>)) and support_vector (with the alias σ (type \sigma<tab>)). The support function of a set $X$ along direction $d ∈ ℝ^n$ corresponds to the (signed) distance of the maximum element in $X$ along direction $d$. The support vector is a (among possibly infinitely many) maximizer of the support function.

The support function also characterizes the half-space that tightly covers a set in a given direction. Hence it can be used to obtain a polyhedral (over-)approximation of any set. This topic is further covered in other sections of the manual.

Specialization

One of the key features of LazySets is specialization. On a more technical level, the library aims at a high level of optimization to get the best possible performance given the following two restrictions: the set type and the set dimension.

The submodule LazySets.Arrays exports the "one-hot vector" SingleEntryVector that can be used, e.g., to represent polyhedra with axis-aligned half-spaces as information on the type. If H is a (non-flat) 100-dimensional hyperrectangle, it has $2^{100}$ vertices. Computing the support function of M*X for any square matrix M along the canonical direction e_{50} = [0, …, 0, 1, 0, …, 0] takes around 10us on a modern computer.

using LazySets.Arrays: SingleEntryVector

using BenchmarkTools

d = SingleEntryVector(50, 100, 1.0)
M = rand(100, 100)  # a random 100×100 matrix
H = rand(Hyperrectangle, dim=100)  # a random 100-dimensional hyperrectangle

out = @benchmark ρ($d, $M * $H)
out
BenchmarkTools.Trial: 10000 samples with 4 evaluations.
 Range (minmax):  7.158 μs 13.402 μs   GC (min … max): 0.00% … 0.00%
 Time  (median):     7.201 μs                GC (median):    0.00%
 Time  (mean ± σ):   7.300 μs ± 307.734 ns   GC (mean ± σ):  0.00% ± 0.00%

  ▇█▄▄▂          ▁▄▄▃▁                                       ▂
  ███████▇▅▅▅▄▃▁█████▇▆▅▆▅▆▅▄▄▁▄▅██▇▆▃▃▃▃▃▃▄▄▄▁▄▁▁▁▁▁▁▁▃▁▅▄ █
  7.16 μs      Histogram: log(frequency) by time      8.62 μs <

 Memory estimate: 896 bytes, allocs estimate: 1.

The Lazy paradigm

For the several set representations, such as polyhedra in constraint and in vertex representation, ellipsoids, balls in different norms, and specific classes of polyhedra (such as intervals, zonotopes, or hyperrectangles), LazySets implements specialized methods for many set operations. Applications that use LazySets can explore different approaches with minimal changes in their code. Conversion between set representations as well as overapproximation and underapproximation functionality are available.

Set operations can be performed in two possible (complementary) ways: concretely or lazily. A concrete operation returns a set of a dedicated type (such as Hyperrectangle). A lazy operation simply returns a wrapper object representing the result of the operation between the given sets. For instance, consider the transformation with a linear map, i.e., the set $Y = \{y : y = Ax \textrm{ for some } x ∈ X\}$. In LazySets, linear_map(A, X) returns a concrete set representation of $Y$. The algorithm that is actually used, as well as the type of $Y$, depend on the types of its arguments (this is the multiple dispatch paradigm) and the dimension of X.

For example, if X is a polygon in vertex representation (VPolygon), then linear_map(M, X) applies M to each vertex of X. However, if X is a 30-dimensional polyhedron in half-space representation (HPolyhedron), then the half-space representation is used if M is invertible – if not, and if MX does not fall into any special case, the vertex representation has to be computed. For polyhedra manipulation we use the library Polyhedra.jl (the cost of converting between representations increases exponentially with the dimension, so computing with concrete polyhedra is usually avoided in high dimensions). Observe that if we are interested in P(MX) for some projection matrix P, that operation can be performed efficiently using the support function, i.e., without actually computing MX. On the other hand, there are special cases of polyhedra for which concrete operations can be performed efficiently. For instance, if X is a hyperrectangle (or a zonotope), the resulting set is a zonotope and the computation is efficient.

To give an example of lazy operations, LinearMap(A, X), or simply A * X, computes the lazy linear map, which is just a new object that wraps the computation of the linear map until it is actually needed. In other words, LinearMap(A, X) can be used to reason about the linear map even if computing the result is expensive (e.g., if X is high-dimensional), since this command just builds an object representing the linear map of A and X. In LazySets, Unicode symbols such as A * X, X ⊕ Y, X ⊖ Y, X × Y, all default to lazy operations by design. Unicode is written using the LaTeX macro of the symbol followed by the TAB key, such as \oplus<tab> for the (lazy) Minkowski sum ⊕.

Finally, one can combine lazy set operations to build lazy expressions that represent several operations between sets, such as Q = (Z ⊕ A*X) × T. By means of the basic tools of (convex) geometry, useful information about Q can be obtained without actually computing the linear map, Minkowski sum and Cartesian product in the above computation. For example, if the computation only involves querying information about Q in a restricted number of directions, as it is often the case in applications, the lazy approach can be realized very efficiently!

How to contribute

Development happens on github. New contributors should follow the links provided in the About section of this documentation.

How to cite

When citing LazySets.jl, please use the entry in CITATION.bib.

Further publications using LazySets can be found in the Publications section of the README. If you would like to list your work, feel free to create a pull request.