The AbstractArray interface in Base Julia is still relatively young. The purpose of this library is to solidify extensions to the current AbstractArray interface, which are put to use in package ecosystems like DifferentialEquations.jl. Since these libraries are live, this package will serve as a staging ground for ideas before they are merged into Base Julia. For this reason, no functionality is exported so that if such functions are added and exported in a future Base Julia, there will be no issues with the upgrade.
Returns the parent array that x
wraps.
Returns true
if the size of T
can change, in which case operations
such as pop!
and popfirst!
are available for collections of type T
.
Given an array x
, this returns the indices along dimension d
. If x
is a tuple
of arrays, then the indices corresponding to dimension d
of all arrays in x
are
returned. If any indices are not equal along dimension d
, an error is thrown. A
tuple may be used to specify a different dimension for each array. If d
is not
specified, then indices for visiting each index of x
are returned.
A trait function for whether x
is a mutable or an immutable array. Used for
dispatching to in-place and out-of-place versions of functions.
Converts an array of structs formulation to a struct of arrays.
A trait function for whether a matrix x
is a sparse structured matrix.
A trait function for whether an array x
can use setindex!
.
Determine whether findstructralnz
accepts the parameter x
.
Returns iterators (I,J)
of the non-zeros in the structure of the matrix x
.
The same as the first two elements of findnz(::SparseMatrixCSC)
.
A trait function for whether matrix_colors(A)
is a fast algorithm or a slow
(graphically-based) method.
Returns an array for the sparsity colors of a matrix type A
. Also includes
an abstract type ColoringAlgorithm
for matrix_colors(A,alg::ColoringAlgorithm)
of non-structured matrices.
A trait function for whether scalar indexing is fast on a given array type.
A getindex
which is always allowed.
A setindex!
which is always allowed.
Returns an instance of the LU factorization object with the correct type cheaply.
Returns an instance of the LU factorization object with the correct type cheaply.
Is a form of vec
which is safe for all values in vector spaces, i.e., if it
is already a vector, like an AbstractVector or Number, it will return said
AbstractVector or Number.
Creates the zero'd matrix version of u
. Note that this is unique because
similar(u,length(u),length(u))
returns a mutable type, so is not type-matching,
while fill(zero(eltype(u)),length(u),length(u))
doesn't match the array type,
i.e., you'll get a CPU array from a GPU array. The generic fallback is
u .* u' .* false
, which works on a surprising number of types, but can be broken
with weird (recursive) broadcast overloads. For higher-order tensors, this
returns the matrix linear operator type, which acts on the vec
of the array.
Restructures the object y
into a shape of x
, keeping its values intact. For
simple objects, like an Array
, this simply amounts to a reshape. However, for
more complex objects such as an ArrayPartition
, not all of the structural
information is adequately contained in the type for standard tools to work. In
these cases, restructure
gives a way to convert, for example, an Array
into
a matching ArrayPartition
.
If first
of instances of type T
are known at compile time, return that first
element. Otherwise, return nothing
. For example, known_first(Base.OneTo{Int})
returns one(Int)
.
If last
of instances of type T
are known at compile time, return that
last element. Otherwise, return nothing
. For example,
known_last(StaticArrays.SOneTo{4})
returns 4.
If step
of instances of type T
are known at compile time, return that step.
Otherwise, returns nothing
. For example, known_step(UnitRange{Int})
returns
one(Int)
.
If length
of an instance of type T
is known at compile time, return it.
Otherwise, return nothing
.
Indicates the most efficient way to access elements from the collection in low-level code.
For GPUArrays
, returns ArrayInterface.GPU()
.
For AbstractArray
supporting a pointer
method, returns ArrayInterface.CPUPointer()
.
For other AbstractArray
s and Tuple
s, returns ArrayInterface.CPUIndex()
.
Otherwise, returns nothing
.
Returns the axis of an array of type T
containing contiguous data.
If no axis is contiguous, it returns Contiguous{-1}
.
If unknown, it returns nothing
.
Returns a tuple of boolean Val
s indicating whether that axis is contiguous.
Returns the size of contiguous batches if !isone(stride_rank(T, contiguous_axis(T)))
.
If isone(stride_rank(T, contiguous_axis(T)))
, then returns ContiguousBatch{0}()
.
If contiguous_axis(T) == -1
, returns ContiguousBatch{-1}()
.
If unknown, returns nothing
.
Returns the rank of each stride.
Returns a Val{true}()
if A
is column major, and a Val{false}()
otherwise.`
Returns a tuple of indicators for whether each axis is dense.
An axis i
of array A
is dense if stride(A, i) * size(A, i) == stride(A, j)
, where j
is the axis (if it exists) such that stride_rank(A)[i] + 1 == stride_rank(A)[j]
.
Returns the size of A
. If the sizes of any axes are known at compile time,
these should be returned as StaticInt
s. For example:
julia> using StaticArrays, ArrayInterface
julia> A = @SMatrix rand(3,4);
julia> ArrayInterface.size(A)
(static(3), static(4))
Returns the strides of array A
. If any strides are known at compile time,
these should be returned as StaticInt
s. For example:
julia> using ArrayInterface
julia> A = rand(3,4);
julia> ArrayInterface.strides(A)
(static(1), 3)
Returns offsets of indices with respect to 0. If values are known at compile time,
it should return them as StaticInt
s.
For example, if A isa Base.Matrix
, offsets(A) === (StaticInt(1), StaticInt(1))
.
Is the function f
whitelisted for LoopVectorization.@avx
?
Creates a static integer with value known at compile time. It is a number,
supporting basic arithmetic. Many operations with two StaticInt
integers
will produce another StaticInt
integer. If one of the arguments to a
function call isn't static (e.g., StaticInt(4) + 3
), then the StaticInt
number will promote to a dynamic value.
- JuliaLang/julia#22216
- JuliaLang/julia#22218
- JuliaLang/julia#22622
- JuliaLang/julia#25760
- JuliaLang/julia#25107
The following common array types are being understood and tested as part of this development.
- Array
- Various versions of sparse arrays
- SArray
- MArray
- FieldVector
- ArrayPartition
- VectorOfArray
- DistributedArrays
- GPUArrays (CLArrays and CuArrays)
- AFArrays
- MultiScaleArrays
- LabelledArrays
2.0: Changed the default of ismutable(array::AbstractArray) = true
. We previously defaulted to
Base.@pure ismutable(array::AbstractArray) = typeof(array).mutable
, but there are a lot of cases
where this tends to not work out in a way one would expect. For example, if you put a normal array
into an immutable struct that adds more information to it, this is considered immutable, even if
all of the setindex!
methods work (by forwarding to the mutable array). Thus, it seems safer to just
always assume mutability is standard for an array, and allow arrays to opt-out.