Arrays

This section of the manual describes the Arrays module.

ReachabilityBase.Arrays.argmaxabsFunction
argmaxabs(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.

source
ReachabilityBase.Arrays.cross_productMethod
cross_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.

source
ReachabilityBase.Arrays.distanceMethod
distance(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 – vector
  • y – vector
  • p – (optional, default: 2.0) the p-norm used; p = 2.0 corresponds to the usual Euclidean norm

Output

A scalar representing $‖ x - y ‖_p$.

source
ReachabilityBase.Arrays.dot_zeroFunction
dot_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 vector
  • y – second vector

Output

The dot product of x and y, but with the rule that 0 * Inf == 0.

source
ReachabilityBase.Arrays.extendFunction
extend(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 – rectangular m × n matrix with m > n and full rank (i.e. its rank is n)
  • check_rank – (optional, default: true) if true, 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ᵀ].

source
ReachabilityBase.Arrays.extend_with_zerosFunction
extend_with_zeros(x::AbstractVector, indices::AbstractVector{<:Int})

Extend a vector with zeros in the given dimensions.

Input

  • x – vector
  • indices – indices to extend, from the interval 1:(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.

source
ReachabilityBase.Arrays.innerFunction
inner(x::AbstractVector{N}, A::AbstractMatrix{N}, y::AbstractVector{N}
     ) where {N}

Compute the inner product $xᵀ A y$.

Input

  • x – vector on the left
  • A – matrix
  • y – vector on the right

Output

The (scalar) result of the multiplication.

source
ReachabilityBase.Arrays.is_cyclic_permutationFunction
is_cyclic_permutation(candidate::AbstractVector, paragon::AbstractVector)

Checks if the elements in candidate are a cyclic permutation of the elements in paragon.

Input

  • candidate – candidate vector
  • paragon – 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.

source
ReachabilityBase.Arrays.is_right_turnFunction
is_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 point
  • u – first 2D direction
  • v – second 2D direction

Output

true iff the vectors constitute a right turn.

source
ReachabilityBase.Arrays.isaboveFunction
isabove(u::AbstractVector, v1::AbstractVector, v2::AbstractVector)

Checks whether the difference v1 - v2 points toward the given direction u.

Input

  • u – direction
  • v1 – first vector
  • v2 – 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.

source
ReachabilityBase.Arrays.isinvertibleFunction
isinvertible(M::Matrix; [cond_tol]::Number=1.0e6)

A sufficient check of a matrix being invertible (or nonsingular).

Input

  • M – matrix
  • cond_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).

source
ReachabilityBase.Arrays.ismultipleFunction
ismultiple(u::AbstractVector{<:Real}, v::AbstractVector{<:Real})

Check whether two vectors are linearly dependent.

Input

  • u – first vector
  • v – 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)
source
ReachabilityBase.Arrays.ispermutationFunction
ispermutation(u::AbstractVector{T}, v::AbstractVector) where {T}

Check that two vectors contain the same elements up to reordering.

Input

  • u – first vector
  • v – 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])
false

Notes

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.

source
ReachabilityBase.Arrays.iswellconditionedFunction
iswellconditioned(M::Matrix; [cond_tol]::Number=1.0e6)

A check of a matrix being sufficiently well-conditioned.

Input

  • M – matrix
  • cond_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.

source
ReachabilityBase.Arrays.matrix_typeFunction
matrix_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 T is a sparse vector or a sparse matrix, the corresponding type is also a sparse matrix.
  • If T is 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.
source
ReachabilityBase.Arrays.nonzero_columnsFunction
nonzero_columns(A::AbstractMatrix; [comparison]=isapproxzero)

Return all columns that have at least one non-zero entry.

Input

  • A – matrix
  • comparison – (optional; default: isapproxzero) function to check for equality with zero

Output

A vector of indices.

source
ReachabilityBase.Arrays.nonzero_indicesFunction
nonzero_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.

source
ReachabilityBase.Arrays.projection_matrixFunction
projection_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 interest
  • n – integer representing the ambient dimension
  • N – (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.0
source
ReachabilityBase.Arrays.rand_pos_neg_zerosum_vectorFunction
rand_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 vector
  • N – (optional; default: Float64) numeric type
  • rng – (optional; default: GLOBAL_RNG) random number generator

Output

A vector as described above.

Algorithm

This is the first phase of the algorithm described here.

source
Base.rationalizeFunction
rationalize(::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 rationals
  • x – vector of floating point numbers
  • tol – (optional, default: eps(N)) tolerance the result at entry i will differ from x[i] by no more than tol

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.

source
ReachabilityBase.Arrays.rectifyFunction
rectify(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.

source
ReachabilityBase.Arrays.remove_zero_columnsFunction
remove_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.

source
ReachabilityBase.Arrays.right_turnFunction
right_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 point
  • u – first 2D point
  • v – 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.

source
ReachabilityBase.Arrays.same_signFunction
same_sign(A::AbstractArray{N}; [optimistic]::Bool=false) where {N}

Check whether all elements of the given array have the same sign.

Input

  • A – array
  • optimistic – (optional; default: false) flag for expressing that the expected result is true

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|\]

source
ReachabilityBase.Arrays.samedirFunction
samedir(u::AbstractVector{<:Real}, v::AbstractVector{<:Real})

Check whether two vectors point in the same direction.

Input

  • u – first vector
  • v – 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)
source
ReachabilityBase.Arrays.substituteFunction
substitute(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.

source
ReachabilityBase.Arrays.substitute!Function
substitute!(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.

source
ReachabilityBase.Arrays.to_matrixFunction
to_matrix(vectors::AbstractVector{VN},
          [m]=length(vectors[1]),
          [MT]=matrix_type(VN)) where {VN}

Input

  • vectors – list of vectors
  • m – (optional; default: length(vectors[1])) number of rows
  • MT – (optional; default: matrix_type(VN)) type of target matrix

Output

A matrix with the column vectors from vectors in the same order.

source
ReachabilityBase.Arrays.uniform_partitionFunction
 uniform_partition(n::Int, block_size::Int)

Compute a uniform block partition of the given size.

Input

  • n – number of dimensions of the partition
  • block_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:2

If 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:10
source
ReachabilityBase.Arrays.vector_typeFunction
vector_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 T is a sparse vector or a sparse matrix, the corresponding type is also a sparse vector.
  • If T is 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.
source