Arrays
This section of the manual describes the Arrays module.
ReachabilityBase.Arrays — ModuleArraysThis module provides machinery for vectors and matrices.
ReachabilityBase.Arrays.argmaxabs — Functionargmaxabs(x::AbstractArray)Return the index with the absolute-wise maximum entry.
Input
x– array
Output
The index i such that |x[i]| >= |x[j]| for all j.
ReachabilityBase.Arrays.cross_product — Methodcross_product(M::AbstractMatrix{N}) where {N<:Real}Compute the high-dimensional cross product of $n-1$ $n$-dimensional vectors.
Input
M– $n × n - 1$-dimensional matrix
Output
A vector.
Algorithm
The cross product is defined as follows:
\[\left[ \dots, (-1)^{n+1} \det(M^{[i]}), \dots \right]^T\]
where $M^{[i]}$ is defined as $M$ with the $i$-th row removed. See Althoff, Stursberg, Buss: Computing Reachable Sets of Hybrid Systems Using a Combination of Zonotopes and Polytopes. 2009.
ReachabilityBase.Arrays.distance — Methoddistance(x::AbstractVector{N}, y::AbstractVector{N}; [p]::Real=N(2)) where {N}Compute the distance between two vectors with respect to the given p-norm, computed as
\[ \|x - y\|_p = \left( \sum_{i=1}^n | x_i - y_i |^p \right)^{1/p}\]
Input
x– vectory– vectorp– (optional, default:2.0) thep-norm used;p = 2.0corresponds to the usual Euclidean norm
Output
A scalar representing $‖ x - y ‖_p$.
ReachabilityBase.Arrays.dot_zero — Functiondot_zero(x::AbstractVector{N}, y::AbstractVector{N}) where{N<:Real}Dot product with preference for zero value in the presence of infinity values.
Input
x– first vectory– second vector
Output
The dot product of x and y, but with the rule that 0 * Inf == 0.
ReachabilityBase.Arrays.extend — Functionextend(M::AbstractMatrix; check_rank=true)Return an invertible extension of M whose first n columns span the column space of M, assuming that size(M) = (m, n), m > n and the rank of M is n.
Input
M– rectangularm × nmatrix withm > nand full rank (i.e. its rank isn)check_rank– (optional, default:true) iftrue, check the rank assumption, otherwise do not perform this check
Output
The tuple (Mext, inv_Mext), where Mext is a square m × m invertible matrix that extends M, i.e. in the sense that Mext = [M | Q2], and the rank of Mext is m. Here, inv_Mext is the inverse of Mext.
Algorithm
First we compute the QR decomposition of M to extract a suitable subspace of column vectors (Q2) that are orthogonal to the column span of M. Then we observe that the inverse of the extended matrix Mext = [M | Q2] is [R⁻¹Qᵀ; Q2ᵀ].
ReachabilityBase.Arrays.extend_with_zeros — Functionextend_with_zeros(x::AbstractVector, indices::AbstractVector{<:Int})Extend a vector with zeros in the given dimensions.
Input
x– vectorindices– indices to extend, from the interval1:(length(x) + length(indices))
Output
A new vector.
Notes
The indices in the extension list are interpreted on the output vector. This is best understood with an example. Let x = [1, 2, 3] and indices = [3, 5]. Then the output vector is y = [1, 2, 0, 3, 0]. Indeed, y[3] == y[5] == 0.
ReachabilityBase.Arrays.hasfullrowrank — Functionhasfullrowrank(M::AbstractMatrix)Check whether a matrix has full row rank.
Input
M– matrix
Output
true iff the matrix has full row rank.
ReachabilityBase.Arrays.inner — Functioninner(x::AbstractVector{N}, A::AbstractMatrix{N}, y::AbstractVector{N}
     ) where {N}Compute the inner product $xᵀ A y$.
Input
x– vector on the leftA– matrixy– vector on the right
Output
The (scalar) result of the multiplication.
ReachabilityBase.Arrays.is_cyclic_permutation — Functionis_cyclic_permutation(candidate::AbstractVector, paragon::AbstractVector)Checks if the elements in candidate are a cyclic permutation of the elements in paragon.
Input
candidate– candidate vectorparagon– paragon vector
Output
A boolean indicating if the elements of candidate are in the same order as in paragon or any of its cyclic permutations.
ReachabilityBase.Arrays.is_right_turn — Functionis_right_turn([O::AbstractVector{<:Real}=[0, 0]], u::AbstractVector{<:Real},
              v::AbstractVector{<:Real})Determine whether the acute angle defined by three 2D points O, u, v in the plane is a right turn (< 180° counter-clockwise) with respect to the center O. Determine if the acute angle defined by two 2D vectors is a right turn (< 180° counter-clockwise) with respect to the center O.
Input
O– (optional; default:[0, 0]) 2D center pointu– first 2D directionv– second 2D direction
Output
true iff the vectors constitute a right turn.
ReachabilityBase.Arrays.isabove — Functionisabove(u::AbstractVector, v1::AbstractVector, v2::AbstractVector)Checks whether the difference v1 - v2 points toward the given direction u.
Input
u– directionv1– first vectorv2– second vector
Output
A Boolean indicating whether the difference of the given vectors points toward the given direction.
Algorithm
The result is equivalent to dot(u, v1 - v2) > 0, but the implementation avoids the allocation of the difference vector.
ReachabilityBase.Arrays.isinvertible — Functionisinvertible(M::Matrix; [cond_tol]::Number=1.0e6)A sufficient check of a matrix being invertible (or nonsingular).
Input
M– matrixcond_tol– (optional, default:1.0e6) tolerance of matrix condition
Output
If the result is true, M is invertible. If the result is false, the matrix is non-square or not well-conditioned.
Algorithm
We check whether the matrix is square and well-conditioned (via iswellconditioned).
ReachabilityBase.Arrays.ismultiple — Functionismultiple(u::AbstractVector{<:Real}, v::AbstractVector{<:Real})Check whether two vectors are linearly dependent.
Input
u– first vectorv– second vector
Output
(true, k) iff the vectors are identical up to a scaling factor k ≠ 0 such that u = k * v, and (false, 0) otherwise.
Examples
julia> ismultiple([1, 2, 3], [2, 4, 6])
(true, 0.5)
julia> ismultiple([1, 2, 3], [3, 2, 1])
(false, 0)
julia> ismultiple([1, 2, 3], [-1, -2, -3])
(true, -1.0)
ReachabilityBase.Arrays.ispermutation — Functionispermutation(u::AbstractVector{T}, v::AbstractVector) where {T}Check that two vectors contain the same elements up to reordering.
Input
u– first vectorv– second vector
Output
true iff the vectors are identical up to reordering.
Examples
julia> ispermutation([1, 2, 2], [2, 2, 1])
true
julia> ispermutation([1, 2, 2], [1, 1, 2])
falseNotes
Containment check is performed using Comparison._in(e, v), so in the case of floating point numbers, the precision to which the check is made is determined by the type of elements in v. See _in and _isapprox for more information.
Note that approximate equality is not an equivalence relation. Hence the result may depend on the order of the elements.
ReachabilityBase.Arrays.issquare — Functionissquare(M::AbstractMatrix)Check whether a matrix is square.
Input
M– matrix
Output
true iff the matrix is square.
ReachabilityBase.Arrays.iswellconditioned — Functioniswellconditioned(M::Matrix; [cond_tol]::Number=1.0e6)A check of a matrix being sufficiently well-conditioned.
Input
M– matrixcond_tol– (optional, default:1.0e6) tolerance of matrix condition
Output
true iff M is well-conditioned subject to the cond_tol parameter.
Algorithm
We check whether the matrix condition number cond(M) is below the prescribed tolerance cond_tol.
ReachabilityBase.Arrays.matrix_type — Functionmatrix_type(T)Return a corresponding matrix type with respect to type T.
Input
T– vector or matrix type
Output
A matrix type that corresponds in some sense (see Notes below) to T.
Notes
- If 
Tis a sparse vector or a sparse matrix, the corresponding type is also a sparse matrix. - If 
Tis a regular vector (i.e.Vector) or a regular matrix (i.e.Matrix), the corresponding type is also a regular matrix. - Otherwise, the corresponding type is a regular matrix.
 
ReachabilityBase.Arrays.nonzero_columns — Functionnonzero_columns(A::AbstractMatrix; [comparison]=isapproxzero)Return all columns that have at least one non-zero entry.
Input
A– matrixcomparison– (optional; default:isapproxzero) function to check for equality with zero
Output
A vector of indices.
ReachabilityBase.Arrays.nonzero_indices — Functionnonzero_indices(v::AbstractVector{N}) where {N<:Real}Return the indices in which a vector is non-zero.
Input
v– vector
Output
A vector of ascending indices i such that the vector is non-zero in dimension i.
ReachabilityBase.Arrays.projection_matrix — Functionprojection_matrix(block::AbstractVector{Int}, n::Int, [N]::Type{<:Number}=Float64)Return the projection matrix associated to the given block of variables.
Input
block– integer vector with the variables of interestn– integer representing the ambient dimensionN– (optional, default:Float64) number type
Output
A sparse matrix that corresponds to the projection onto the variables in block.
Examples
julia> proj = projection_matrix([1, 3], 4)
2×4 SparseArrays.SparseMatrixCSC{Float64, Int64} with 2 stored entries:
 1.0   ⋅    ⋅    ⋅
  ⋅    ⋅   1.0   ⋅
julia> Matrix(proj)
2×4 Matrix{Float64}:
 1.0  0.0  0.0  0.0
 0.0  0.0  1.0  0.0ReachabilityBase.Arrays.rand_pos_neg_zerosum_vector — Functionrand_pos_neg_zerosum_vector(n::Int; [N]::Type{<:Real}=Float64,
                                    [rng]::AbstractRNG=GLOBAL_RNG)Create a vector of random numbers such that the total sum is (approximately) zero, no duplicates exist, all positive entries come first, and all negative entries come last.
Input
n– length of the vectorN– (optional; default:Float64) numeric typerng– (optional; default:GLOBAL_RNG) random number generator
Output
A vector as described above.
Algorithm
This is the first phase of the algorithm described here.
Base.rationalize — Functionrationalize(::Type{T}, x::AbstractArray{N}, tol::Real) where {T<:Integer, N<:AbstractFloat}Approximate an array of floating point numbers as a rational vector with entries of the given integer type.
Input
T– (optional, default:Int) integer type to represent the rationalsx– vector of floating point numberstol– (optional, default:eps(N)) tolerance the result at entryiwill differ fromx[i]by no more thantol
Output
An array of type Rational{T} where the i-th entry is the rationalization of the i-th component of x.
Notes
See also Base.rationalize.
ReachabilityBase.Arrays.rectify — Functionrectify(x::AbstractArray{N}) where {N<:Real}Rectify an array, i.e., take the element-wise maximum with zero.
Input
x– array
Output
A copy of the array where each negative entry is replaced by zero.
ReachabilityBase.Arrays.remove_duplicates_sorted! — Functionremove_duplicates_sorted!(v::AbstractVector)Remove duplicate entries in a sorted vector.
Input
v– sorted vector
Output
The input vector without duplicates.
ReachabilityBase.Arrays.remove_zero_columns — Functionremove_zero_columns(A::AbstractMatrix)Return a matrix with all columns containing only zero entries removed.
Input
A– matrix
Output
The original matrix A if it contains no zero columns or otherwise a new matrix where those columns have been removed.
ReachabilityBase.Arrays.right_turn — Functionright_turn([O::AbstractVector{<:Real}=[0, 0]], u::AbstractVector{<:Real},
           v::AbstractVector{<:Real})Compute a scalar that determines whether the acute angle defined by three 2D points O, u, v in the plane is a right turn (< 180° counter-clockwise) with respect to the center O.
Input
O– (optional; default:[0, 0]) 2D center pointu– first 2D pointv– second 2D point
Output
A scalar representing the rotation. If the result is 0, the points are collinear; if it is positive, the points constitute a positive angle of rotation around O from u to v; otherwise they constitute a negative angle.
Algorithm
The cross product is used to determine the sense of rotation.
ReachabilityBase.Arrays.same_sign — Functionsame_sign(A::AbstractArray{N}; [optimistic]::Bool=false) where {N}Check whether all elements of the given array have the same sign.
Input
A– arrayoptimistic– (optional; default:false) flag for expressing that the expected result istrue
Output
true if and only if all elements in M have the same sign.
Algorithm
If optimistic is false, we check the sign of the first element and compare to the sign of all elements.
If optimistic is true, we compare the absolute element sum with the sum of the absolute of the elements; this is faster if the result is true because there is no branching.
\[ |\sum_i A_i| = \sum_i |A_i|\]
ReachabilityBase.Arrays.samedir — Functionsamedir(u::AbstractVector{<:Real}, v::AbstractVector{<:Real})Check whether two vectors point in the same direction.
Input
u– first vectorv– second vector
Output
(true, k) iff the vectors are identical up to a positive scaling factor k such that u = k * v, and (false, 0) otherwise.
Examples
julia> samedir([1, 2, 3], [2, 4, 6])
(true, 0.5)
julia> samedir([1, 2, 3], [3, 2, 1])
(false, 0)
julia> samedir([1, 2, 3], [-1, -2, -3])
(false, 0)
ReachabilityBase.Arrays.SingleEntryVector — TypeSingleEntryVector{N} <: AbstractVector{N}A sparse unit vector with arbitrary one-element.
Fields
i– index of non-zero entryn– vector lengthv– non-zero entry
ReachabilityBase.Arrays.substitute — Functionsubstitute(substitution::Dict{Int, T}, x::AbstractVector{T}) where {T}Apply a substitution to a given vector.
Input
substitution– substitution (a mapping from an index to a new value)x– vector
Output
A fresh vector corresponding to x after substitution was applied.
ReachabilityBase.Arrays.substitute! — Functionsubstitute!(substitution::Dict{Int, T}, x::AbstractVector{T}) where {T}Apply a substitution to a given vector.
Input
substitution– substitution (a mapping from an index to a new value)x– vector (modified in this function)
Output
The same (but see the Notes below) vector x but after substitution was applied.
Notes
The vector x is modified in-place if it has type Vector or SparseVector. Otherwise, we first create a new Vector from it.
ReachabilityBase.Arrays.to_matrix — Functionto_matrix(vectors::AbstractVector{VN},
          [m]=length(vectors[1]),
          [MT]=matrix_type(VN)) where {VN}Input
vectors– list of vectorsm– (optional; default:length(vectors[1])) number of rowsMT– (optional; default:matrix_type(VN)) type of target matrix
Output
A matrix with the column vectors from vectors in the same order.
ReachabilityBase.Arrays.to_negative_vector — Functionto_negative_vector(v::AbstractVector{N}) where {N}Negate a vector and convert to type Vector.
Input
v– vector
Output
A Vector equivalent to $-v$.
ReachabilityBase.Arrays.uniform_partition — Function uniform_partition(n::Int, block_size::Int)Compute a uniform block partition of the given size.
Input
n– number of dimensions of the partitionblock_size– size of each block
Output
A vector of ranges, Vector{UnitRange{Int}}, such that the size of each block is the same, if possible.
Examples
If the number of dimensions n is 2, we have two options: either two blocks of size 1 or one block of size 2:
julia> uniform_partition(2, 1)
2-element Vector{UnitRange{Int64}}:
 1:1
 2:2
julia> uniform_partition(2, 2)
1-element Vector{UnitRange{Int64}}:
 1:2If the block_size argument is not compatible with (i.e., does not divide) n, the output is filled with one block of the size needed to reach n:
julia> uniform_partition(3, 1)
3-element Vector{UnitRange{Int64}}:
 1:1
 2:2
 3:3
julia> uniform_partition(3, 2)
2-element Vector{UnitRange{Int64}}:
 1:2
 3:3
julia> uniform_partition(10, 6)
2-element Vector{UnitRange{Int64}}:
 1:6
 7:10ReachabilityBase.Arrays.vector_type — Functionvector_type(T)Return a corresponding vector type with respect to type T.
Input
T– vector or matrix type
Output
A vector type that corresponds in some sense (see Notes below) to T.
Notes
- If 
Tis a sparse vector or a sparse matrix, the corresponding type is also a sparse vector. - If 
Tis a regular vector (i.e.Vector) or a regular matrix (i.e.Matrix), the corresponding type is also a regular vector. - Otherwise, the corresponding type is a regular vector.