Документация Engee

API

Страница в процессе перевода.

Guide

SVector

The simplest static array is the type SVector{N,T}, which provides an immutable vector of fixed length N and type T.

SVector defines a series of convenience constructors, so you can just type e.g. SVector(1,2,3). Alternatively there is an intelligent @SVector macro where you can use native Julia array literals syntax, comprehensions, and the zeros(), ones(), fill(), rand() and randn() functions, such as @SVector [1,2,3], @SVector Float64[1,2,3], @SVector [f(i) for i = 1:10], @SVector zeros(3), @SVector randn(Float32, 4), etc (Note: the range of a comprehension is evaluated at global scope by the macro, and must be made of combinations of literal values, functions, or global variables, but is not limited to just simple ranges. Extending this to (hopefully statically known by type-inference) local-scope variables is hoped for the future. The zeros(), ones(), fill(), rand(), randn(), and randexp() functions do not have this limitation.)

SMatrix

Statically sized N×M matrices are provided by SMatrix{N,M,T,L}.

Here L is the length of the matrix, such that N × M = L. However, convenience constructors are provided, so that L, T and even M are unnecessary. At minimum, you can type SMatrix{2}(1,2,3,4) to create a 2×2 matrix (the total number of elements must divide evenly into N). A convenience macro @SMatrix [1 2; 3 4] is provided (which also accepts comprehensions and the zeros(), ones(), fill(), rand(), randn(), and randexp() functions).

SArray

A container with arbitrarily many dimensions is defined as struct SArray{Size,T,N,L} <: StaticArray{Size,T,N}, where Size = Tuple{S1, S2, ...} is a tuple of Ints. You can easily construct one with the @SArray macro, supporting all the features of @SVector and @SMatrix (but with arbitrary dimension).

The main reason SVector and SMatrix are defined is to make it easier to define the types without the extra tuple characters (compare SVector{3} to SArray{Tuple{3}}).

Scalar

Sometimes you want to broadcast an operation, but not over one of your inputs. A classic example is attempting to displace a collection of vectors by the same vector. We can now do this with the Scalar type:

[[1,2,3], [4,5,6]] .+ Scalar([1,0,-1]) # [[2,2,2], [5,5,5]]

Scalar is simply an implementation of an immutable, 0-dimensional StaticArray.

The Size trait

The size of a statically sized array is a static parameter associated with the type of the array. The Size trait is provided as an abstract representation of the dimensions of a static array. An array sa::SA of size (dims...) is associated with Size{(dims...)}(). The following are equivalent (@pure) constructors:

Size{(dims...,)}()
Size(dims...)
Size(sa::StaticArray)
Size(SA) # SA <: StaticArray

This is extremely useful for (a) performing dispatch depending on the size of an array, and (b) passing array dimensions that the compiler can reason about.

An example of size-based dispatch for the determinant of a matrix would be:

det(x::StaticMatrix) = _det(Size(x), x)
_det(::Size{(1,1)}, x::StaticMatrix) = x[1,1]
_det(::Size{(2,2)}, x::StaticMatrix) = x[1,1]*x[2,2] - x[1,2]*x[2,1]
# and other definitions as necessary

Examples of using Size as a compile-time constant include

reshape(svector, Size(2,2))  # Convert SVector{4} to SMatrix{2,2}
SizedMatrix{3,3}(rand(3,3))  # Construct a random 3×3 SizedArray (see below)

Indexing

Statically sized indexing can be realized by indexing each dimension by a scalar, a StaticVector or :. Indexing in this way will result a statically sized array (even if the input was dynamically sized, in the case of StaticVector indices) of the closest type (as defined by similar_type).

Conversely, indexing a statically sized array with a dynamically sized index (such as a Vector{Integer} or UnitRange{Integer}) will result in a standard (dynamically sized) Array.

similar_type()

Since immutable arrays need to be constructed "all-at-once", we need a way of obtaining an appropriate constructor if the element type or dimensions of the output array differs from the input. To this end, similar_type is introduced, behaving just like similar, except that it returns a type. Relevant methods are:

similar_type(::Type{A}) where {A <: StaticArray} # defaults to A
similar_type(::Type{A}, ::Type{ElType}) where {A <: StaticArray, ElType} # Change element type
similar_type(::Type{A}, size::Size) where {A <: AbstractArray} # Change size
similar_type(::Type{A}, ::Type{ElType}, size::Size) where {A <: AbstractArray, ElType} # Change both

These setting will affect everything, from indexing, to matrix multiplication and broadcast. Users wanting introduce a new array type should only overload the last method in the above.

Use of similar will fall back to a mutable container, such as a MVector (see below), and it requires use of the Size trait if you wish to set a new static size (or else a dynamically sized Array will be generated when specifying the size as plain integers).

Collecting directly into static arrays

You can collect iterators into static arrays directly with StaticArrays.sacollect. The size needs to be specified, but the element type is optional.

Mutable arrays: MVector, MMatrix and MArray

These statically sized arrays are identical to the above, but are defined as mutable structs, instead of immutable structs. Because they are mutable, they allow setindex! to be defined (achieved through pointer manipulation, into a tuple).

As a consequence of Julia’s internal implementation, these mutable containers live on the heap, not the stack. Their memory must be allocated and tracked by the garbage collector. Nevertheless, there is opportunity for speed improvements relative to Base.Array because (a) there may be one less pointer indirection, (b) their (typically small) static size allows for additional loop unrolling and inlining, and consequentially (c) their mutating methods like map! are extremely fast. Benchmarking shows that operations such as addition and matrix multiplication are faster for MMatrix than Matrix, at least for sizes up to 14 × 14, though keep in mind that optimal speed will be obtained by using mutating functions (like map! or A_mul_B!) where possible, rather than reallocating new memory.

Mutable static arrays also happen to be very useful containers that can be constructed on the heap (with the ability to use setindex!, etc), and later copied as e.g. an immutable SVector to the stack for use, or into e.g. an Array{SVector} for storage.

Convenience macros @MVector, @MMatrix and @MArray are provided.

SizedArray: a decorate size wrapper for Array

Another convenient mutable type is the SizedArray, which is just a wrapper-type about a standard Julia Array which declares its known size. For example, if we knew that a was a 2×2 Matrix, then we can type sa = SizedArray{Tuple{2,2}}(a) to construct a new object which knows the type (the size will be verified automatically). For one and two dimensions, a more convenient syntax for obtaining a SizedArray is by using the SizedMatrix and SizedVector aliases, e.g. sa = SizedMatrix{2,2}(a).

Then, methods on sa will use the specialized code provided by the StaticArrays package, which in many cases will be much, much faster. For example, calling eigen(sa) will be significantly faster than eigen(a) since it will perform a specialized 2×2 matrix diagonalization rather than a general algorithm provided by Julia and LAPACK.

In some cases it will make more sense to use a SizedArray, and in other cases an MArray might be preferable.

FieldVector

Sometimes it is useful to give your own struct types the properties of a vector. StaticArrays can take care of this for you by allowing you to inherit from FieldVector{N, T}. For example, consider:

struct Point3D <: FieldVector{3, Float64}
    x::Float64
    y::Float64
    z::Float64
end

With this type, users can easily access fields to p = Point3D(x,y,z) using p.x, p.y or p.z, or alternatively via p[1], p[2], or p[3]. You may even permute the coordinates with p[SVector(3,2,1)]). Furthermore, Point3D is a complete AbstractVector implementation where you can add, subtract or scale vectors, multiply them by matrices, etc.

the three components of an ordinary v::SVector{3} can also be accessed as v.x, v.y, and v.z, so there is no need for a FieldVector to use this convention.

It is also worth noting that FieldVectors may be mutable or immutable, and that setindex! is defined for use on mutable types. For immutable containers, you may want to define a method for similar_type so that operations leave the type constant (otherwise they may fall back to SVector). For mutable containers, you may want to define a default constructor (no inputs) and an appropriate method for similar,

Implementing your own types

You can easily create your own StaticArray type, by defining linear getindex (and optionally setindex! for mutable types --- see setindex!(::MArray, val, i) in MArray.jl for an example of how to achieve this through pointer manipulation). Your type should define a constructor that takes a tuple of the data (and mutable containers may want to define a default constructor).

Other useful functions to overload may be similar_type (and similar for mutable containers).

Conversions from Array

In order to convert from a dynamically sized AbstractArray to one of the statically sized array types, you must specify the size explicitly. For example,

v = [1,2]

m = [1 2;
     3 4]

# ... a lot of intervening code

sv = SVector{2}(v)
sm = SMatrix{2,2}(m)
sa = SArray{Tuple{2,2}}(m)

sized_v = SizedVector{2}(v)
sized_m = SizedMatrix{2,2}(m)

We have avoided adding SVector(v::AbstractVector) as a valid constructor to help users avoid the type instability (and potential performance disaster, if used without care) of this innocuous looking expression.

Arrays of static arrays

Storing a large number of static arrays is convenient as an array of static arrays. For example, a collection of positions (3D coordinates --- SVector{3,Float64}) could be represented as a Vector{SVector{3,Float64}}.

Another common way of storing the same data is as a 3×N Matrix{Float64}. Rather conveniently, such types have exactly the same binary layout in memory, and therefore we can use reinterpret to convert between the two formats

function svectors(x::Matrix{T}, ::Val{N}) where {T,N}
    size(x,1) == N || error("sizes mismatch")
    isbitstype(T) || error("use for bitstypes only")
    reinterpret(SVector{N,T}, vec(x))
end

Such a conversion does not copy the data, rather it refers to the same memory. Arguably, a Vector of SVectors is often preferable to a Matrix because it provides a better abstraction of the objects contained in the array and it allows the fast StaticArrays methods to act on elements.

However, the resulting object is a Base.ReinterpretArray, not an Array, which carries some runtime penalty on every single access. If you can afford the memory for a copy and can live with the non-shared mutation semantics, then it is better to pull a copy by e.g.

function svectorscopy(x::Matrix{T}, ::Val{N}) where {T,N}
    size(x,1) == N || error("sizes mismatch")
    isbitstype(T) || error("use for bitstypes only")
    copy(reinterpret(SVector{N,T}, vec(x)))
end

For example:

julia> M=reshape(collect(1:6), (2,3))
2×3 Array{Int64,2}:
 1  3  5
 2  4  6

julia> svectors(M, Val{2}())
3-element reinterpret(SArray{Tuple{2},Int64,1,2}, ::Array{Int64,1}):
 [1, 2]
 [3, 4]
 [5, 6]

julia> svectorscopy(M, Val{2}())
3-element Array{SArray{Tuple{2},Int64,1,2},1}:
 [1, 2]
 [3, 4]
 [5, 6]

Working with mutable and immutable arrays

Generally, it is performant to rebind an immutable array, such as

function average_position(positions::Vector{SVector{3,Float64}})
    x = zeros(SVector{3,Float64})
    for pos ∈ positions
        x = x + pos
    end
    return x / length(positions)
end

so long as the Type of the rebound variable (x, above) does not change.

On the other hand, the above code for mutable containers like Array, MArray or SizedArray is not very efficient. Mutable containers must be allocated and later garbage collected, and for small, fixed-size arrays this can be a leading contribution to the cost. In the above code, a new array will be instantiated and allocated on each iteration of the loop. In order to avoid unnecessary allocations, it is best to allocate an array only once and apply mutating functions to it:

function average_position(positions::Vector{SVector{3,Float64}})
    x = zeros(MVector{3,Float64})
    for pos ∈ positions
        x .+= pos
    end
    x ./= length(positions)
    return x
end

The functions setindex, push, pop, pushfirst, popfirst, insert and deleteat are provided for performing certain specific operations on static arrays, in analogy with the standard functions setindex!, push!, pop!, etc. (Note that if the size of the static array changes, the type of the output will differ from the input.)

When building static arrays iteratively, it is usually efficient to build up an MArray first and then convert. The allocation will be elided by recent Julia compilers, resulting in very efficient code:

function standard_basis_vector(T, ::Val{I}, ::Val{N}) where {I,N}
    v = zero(MVector{N,T})
    v[I] = one(T)
    SVector(v)
end

SIMD optimizations

It seems Julia and LLVM are smart enough to use processor vectorization extensions like SSE and AVX - however they are currently partially disabled by default. Run Julia with julia -O or julia -O3 to enable these optimizations, and many of your (immutable) StaticArray methods should become significantly faster!

Docstrings

# StaticArrays.StaticMatMulLikeType

StaticMatMulLike

Static wrappers used for multiplication dispatch.

# StaticArrays.ArgsType

Args

A help wrapper to distinguish SA(x...) and SA((x...,))

# StaticArrays.SAType

SA[ elements ]
SA{T}[ elements ]

Create SArray literals using array construction syntax. The element type is inferred by promoting elements to a common type or set to T when T is provided explicitly.

Examples:

  • SA[1.0, 2.0] creates a length-2 SVector of Float64 elements.

  • SA[1 2; 3 4] creates a 2×2 SMatrix of Ints.

  • SA[1 2] creates a 1×2 SMatrix of Ints.

  • SA{Float32}[1, 2] creates a length-2 SVector of Float32 elements.

A couple of helpful type aliases are also provided:

  • SA_F64[1, 2] creates a length-2 SVector of Float64 elements

  • SA_F32[1, 2] creates a length-2 SVector of Float32 elements

# StaticArrays.SHermitianCompactType

SHermitianCompact{N, T, L} <: StaticMatrix{N, N, T}

A StaticArray subtype that represents a Hermitian matrix. Unlike LinearAlgebra.Hermitian, SHermitianCompact stores only the lower triangle of the matrix (as an SVector). The lower triangle is stored in column-major order. For example, for an SHermitianCompact{3}, the indices of the stored elements can be visualized as follows:

┌ 1 ⋅ ⋅ ┐
| 2 4 ⋅ |
└ 3 5 6 ┘

Type parameters:

  • N: matrix dimension;

  • T: element type for lower triangle;

  • L: length of the SVector storing the lower triangular elements.

Note that L is always the Nth triangular number.

An SHermitianCompact may be constructed either:

  • from an AbstractVector containing the lower triangular elements; or

  • from a Tuple containing both upper and lower triangular elements in column major order; or

  • from another StaticMatrix.

For the latter two cases, only the lower triangular elements are used; the upper triangular elements are ignored.

# StaticArrays.SOneToType

SOneTo(n)

Return a statically-sized AbstractUnitRange starting at 1, functioning as the axes of a StaticArray.

# StaticArrays.ScalarType

Scalar{T}(x::T)

Construct a statically-sized 0-dimensional array that contains a single element, x. This type is particularly useful for influencing broadcasting operations.

# StaticArrays.TSizeType

Size that stores whether a Matrix is a Transpose Useful when selecting multiplication methods, and avoiding allocations when dealing with the Transpose type by passing around the original matrix. Should pair with parent.

# StaticArrays._InitialValue — _Type

_InitialValue

A singleton type for representing "universal" initial value (identity element).

The idea is that, given op for mapfoldl, virtually, we define an "extended" version of it by

op′(::_InitialValue, x) = x
op′(acc, x) = op(acc, x)

This is just a conceptually useful model to have in mind and we don’t actually define op′ here (yet?). But see Base.BottomRF for how it might work in action.

(It is related to that you can always turn a semigroup without an identity into a monoid by "adjoining" an element that acts as the identity.)

# Base.setindexMethod

setindex(vec::StaticArray, x, index::Int)

Return a new array with the item at index replaced by x.

Examples

julia> setindex(@SVector[1,2,3], 4, 2)
3-element SVector{3, Int64} with indices SOneTo(3):
 1
 4
 3

julia> setindex(@SMatrix[2 4; 6 8], 1, 2)
2×2 SMatrix{2, 2, Int64, 4} with indices SOneTo(2)×SOneTo(2):
 2  4
 1  8

# Base.similarMethod

similar(static_array)
similar(static_array, T)
similar(array, ::Size)
similar(array, T, ::Size)

Constructs and returns a mutable but statically-sized array (i.e. a StaticArray). If the input array is not a StaticArray, then the Size is required to determine the output size (or else a dynamically sized array will be returned).

# LinearAlgebra.qrMethod

qr(A::StaticMatrix,
   pivot::Union{Val{true}, Val{false}, LinearAlgebra.PivotingStrategy} = Val(false))

Compute the QR factorization of A. The factors can be obtained by iteration:

julia> A = @SMatrix rand(3,4);

julia> Q, R = qr(A);

julia> Q * R ≈ A
true

or by using getfield:

julia> F = qr(A);

julia> F.Q * F.R ≈ A
true

# StaticArrays._construct_similar — _Method

_construct_similar(a, ::Size, elements::NTuple)

Construct a static array of similar type to a with the given elements.

When a is an instance or a concrete type the element type eltype(a) is used. However, when a is a UnionAll type such as SMatrix{2,2}, the promoted type of elements is used instead.

# StaticArrays._lind — _Method

Obtain an expression for the linear index of var[k,j], taking transposes into account

# StaticArrays._muladd_expr — _Method

Combine left and right sides of an assignment expression, short-cutting lhs = α * rhs + β * lhs, element-wise. If α = 1, the multiplication by α is removed. If β = 0, the second rhs term is removed.

# StaticArrays._size — _Method

Return either the statically known Size() or runtime size()

# StaticArrays.arithmetic_closureMethod

arithmetic_closure(T)

Return the type which values of type T will promote to under a combination of the arithmetic operations +, -, * and /.

julia> import StaticArrays.arithmetic_closure

julia> arithmetic_closure(Bool)
Float64

julia> arithmetic_closure(Int32)
Float64

julia> arithmetic_closure(BigFloat)
BigFloat

julia> arithmetic_closure(BigInt)
BigFloat

# StaticArrays.check_dimsMethod

Validate the dimensions of a matrix multiplication, including matrix-vector products

# StaticArrays.construct_typeMethod

SA′ = construct_type(::Type{SA}, x) where {SA<:StaticArray}

Pick a proper constructor SA′ based on x if SA(x)/SA(x...) has no specific definition. The default returned SA′ is SA itself for user defined StaticArrays. This differs from similar_type() in that SA′ should always be a subtype of SA.

To distinguish SA(x...) and SA(x::Tuple), the former calls construct_type(SA, StaticArrays.Args(x)) instead of construct_type(SA, x).

Please make sure SA'(x) has a specific definition if the default behavior is overloaded. Otherwise construction might fall into infinite recursion.


The adaption rules for official StaticArrays could be summarized as:

SA <: FieldArray: eltype adaptable

FieldArrays are always static-sized. We only derive SA′'s eltype using type promotion if needed.

SA <: Union{SArray, MArray, SHermitianCompact, SizedArray}: size/eltype adaptable

  • SA(x::Tuple)

If SA is fully static-sized, then we first try to fill SA with x's elements. If failed and length(SA) == 1, then we try to fill SA with x itself.

If SA is not fully static-sized, then we always try to fill SA with x's elements, and the constructor’s Size is derived based on:

  1. If SA <: StaticVector, then we use length(x) as the output Length

  2. If SA <: StaticMatrix{M}, then we use (M, N) (N = length(x) ÷ M) as the output Size

  3. If SA <: StaticMatrix{M,M} where M, then we use (N, N) (N = sqrt(length(x)) as the output Size.

    • SA(x…​)

Similar to Tuple, but we never fill SA with x itself.

  • SA(x::StaticArray)

We treat x as Tuple whenever possible. If failed, then try to inherit x's Size.

  • SA(x::AbstractArray)

x is used to provide eltype. Thus SA must be static sized.

# StaticArrays.deleteatMethod

deleteat(vec::StaticVector, index::Integer)

Return a new vector with the item at the given index removed.

Examples

julia> deleteat(@SVector[6, 5, 4, 3, 2, 1], 2)
5-element SVector{5, Int64} with indices SOneTo(5):
 6
 4
 3
 2
 1

# StaticArrays.dimmatchFunction

dimmatch(x::StaticDimension, y::StaticDimension)

Return whether dimensions x and y match at compile time, that is:

  • if x and y are both Ints, check that they are equal

  • if x or y are Dynamic(), return true

# StaticArrays.gen_by_accessFunction

gen_by_access(expr_gen, a::Type{<:AbstractArray}, asym = :wrapped_a)

Statically generate outer code for fully unrolled multiplication loops. Returned code does wrapper-specific tests (for example if a symmetric matrix view is U or L) and the body of the if expression is then generated by function expr_gen. The function expr_gen receives access pattern description symbol as its argument and this symbol is then consumed by uplo_access to generate the right code for matrix element access.

The name of the matrix to test is indicated by asym.

# StaticArrays.gen_by_accessMethod

gen_by_access(expr_gen, a::Type{<:AbstractArray}, b::Type{<:AbstractArray})

Similar to genbyaccess with only one type argument. The difference is that tests for both arrays of type a and b are generated and expr_gen receives two access arguments, first for matrix a and the second for matrix b.

# StaticArrays.insertMethod

insert(vec::StaticVector, index::Integer, item)

Return a new vector with item inserted into vec at the given index.

Examples

julia> insert(@SVector[6, 5, 4, 2, 1], 4, 3)
6-element SVector{6, Int64} with indices SOneTo(6):
 6
 5
 4
 3
 2
 1

# StaticArrays.mul_result_structureMethod

mul_result_structure(a::Type, b::Type)

Get a structure wrapper that should be applied to the result of multiplication of matrices of given types (a*b).

# StaticArrays.multiplied_dimensionMethod

Calculate the product of the dimensions being multiplied. Useful as a heuristic for unrolling.

# StaticArrays.popMethod

pop(vec::StaticVector)

Return a new vector with the last item in vec removed.

Examples

julia> pop(@SVector[1,2,3])
2-element SVector{2, Int64} with indices SOneTo(2):
 1
 2

# StaticArrays.popfirstMethod

popfirst(vec::StaticVector)

Return a new vector with the first item in vec removed.

Examples

julia> popfirst(@SVector[1,2,3])
2-element SVector{2, Int64} with indices SOneTo(2):
 2
 3

# StaticArrays.pushMethod

push(vec::StaticVector, item)

Return a new StaticVector with item inserted on the end of vec.

Examples

julia> push(@SVector[1, 2, 3], 4)
4-element SVector{4, Int64} with indices SOneTo(4):
 1
 2
 3
 4

# StaticArrays.pushfirstMethod

pushfirst(vec::StaticVector, item)

Return a new StaticVector with item inserted at the beginning of vec.

Examples

julia> pushfirst(@SVector[1, 2, 3, 4], 5)
5-element SVector{5, Int64} with indices SOneTo(5):
 5
 1
 2
 3
 4

# StaticArrays.sacollectFunction

sacollect(SA, gen)

Construct a statically-sized vector of type SA.from a generator gen. SA needs to have a size parameter since the length of vec is unknown to the compiler. SA can optionally specify the element type as well.

Example:

sacollect(SVector{3, Int}, 2i+1 for i in 1:3)
sacollect(SMatrix{2, 3}, i+j for i in 1:2, j in 1:3)
sacollect(SArray{2, 3}, i+j for i in 1:2, j in 1:3)

This creates the same statically-sized vector as if the generator were collected in an array, but is more efficient since no array is allocated.

Equivalent:

SVector{3, Int}([2i+1 for i in 1:3])

# StaticArrays.same_sizeMethod

Returns the common Size of the inputs (or else throws a DimensionMismatch)

# StaticArrays.sizematchMethod

sizematch(::Size, ::Size)
sizematch(::Tuple, ::Tuple)

Determine whether two sizes match, in the sense that they have the same number of dimensions, and their dimensions match as determined by dimmatch.

# StaticArrays.sizematchMethod

sizematch(::Size, A::AbstractArray)

Determine whether array A matches the given size. If A is a StaticArray, the check is performed at compile time, otherwise, the check is performed at runtime.

# StaticArrays.uplo_accessMethod

uplo_access(sa, asym, k, j, uplo)

Generate code for matrix element access, for a matrix of size sa locally referred to as asym in the context where the result will be used. Both indices k and j need to be statically known for this function to work. uplo is the access pattern mode generated by the gen_by_access function.

# StaticArrays.@MArrayMacro

@MArray [a b; c d]
@MArray [[a, b];[c, d]]
@MArray [i+j for i in 1:2, j in 1:2]
@MArray ones(2, 2, 2)

A convenience macro to construct MArray with arbitrary dimension. See @SArray for detailed features.

# StaticArrays.@MMatrixMacro

@MMatrix [a b c d]
@MMatrix [[a, b];[c, d]]
@MMatrix [i+j for i in 1:2, j in 1:2]
@MMatrix ones(2, 2)

A convenience macro to construct MMatrix. See @SArray for detailed features.

# StaticArrays.@MVectorMacro

@MVector [a, b, c, d]
@MVector [i for i in 1:2]
@MVector ones(2)

A convenience macro to construct MVector. See @SArray for detailed features.

# StaticArrays.@SArrayMacro

@SArray [a b; c d]
@SArray [[a, b];[c, d]]
@SArray [i+j for i in 1:2, j in 1:2]
@SArray ones(2, 2, 2)

A convenience macro to construct SArray with arbitrary dimension. It supports:

  1. (typed) array literals.

Every argument inside the square brackets is treated as a scalar during expansion. Thus @SArray[a; b] always forms a SVector{2} and @SArray [a b; c] always throws an error.

  1. comprehensions

The range of a comprehension is evaluated at global scope by the macro, and must be made of combinations of literal values, functions, or global variables.

  1. initialization functions

Only support zeros(), ones(), fill(), rand(), randn(), and randexp()

# StaticArrays.@SMatrixMacro

@SMatrix [a b c d]
@SMatrix [[a, b];[c, d]]
@SMatrix [i+j for i in 1:2, j in 1:2]
@SMatrix ones(2, 2)

A convenience macro to construct SMatrix. See @SArray for detailed features.

# StaticArrays.@SVectorMacro

@SVector [a, b, c, d]
@SVector [i for i in 1:2]
@SVector ones(2)

A convenience macro to construct SVector. See @SArray for detailed features.

# StaticArraysCore.DynamicType

Dynamic()

Used to signify that a dimension of an array is not known statically.

# StaticArraysCore.FieldArrayType

abstract FieldArray{N, T, D} <: StaticArray{N, T, D}

Inheriting from this type will make it easy to create your own rank-D tensor types. A FieldArray will automatically define getindex and setindex! appropriately. An immutable FieldArray will be as performant as an SArray of similar length and element type, while a mutable FieldArray will behave similarly to an MArray.

Note that you must define the fields of any FieldArray subtype in column major order. If you want to use an alternative ordering you will need to pay special attention in providing your own definitions of getindex, setindex! and tuple conversion.

If you define a FieldArray which is parametric on the element type you should consider defining similar_type as in the FieldVector example.

Example

struct Stiffness <: FieldArray{Tuple{2,2,2,2}, Float64, 4}
    xxxx::Float64
    yxxx::Float64
    xyxx::Float64
    yyxx::Float64
    xxyx::Float64
    yxyx::Float64
    xyyx::Float64
    yyyx::Float64
    xxxy::Float64
    yxxy::Float64
    xyxy::Float64
    yyxy::Float64
    xxyy::Float64
    yxyy::Float64
    xyyy::Float64
    yyyy::Float64
end

# StaticArraysCore.FieldMatrixType

abstract FieldMatrix{N1, N2, T} <: FieldArray{Tuple{N1, N2}, 2}

Inheriting from this type will make it easy to create your own rank-two tensor types. A FieldMatrix will automatically define getindex and setindex! appropriately. An immutable FieldMatrix will be as performant as an SMatrix of similar length and element type, while a mutable FieldMatrix will behave similarly to an MMatrix.

Note that the fields of any subtype of FieldMatrix must be defined in column major order unless you are willing to implement your own getindex.

If you define a FieldMatrix which is parametric on the element type you should consider defining similar_type as in the FieldVector example.

Example

struct Stress <: FieldMatrix{3, 3, Float64}
    xx::Float64
    yx::Float64
    zx::Float64
    xy::Float64
    yy::Float64
    zy::Float64
    xz::Float64
    yz::Float64
    zz::Float64
end

Note that the fields of any subtype of FieldMatrix must be defined in column major order. This means that formatting of constructors for literal FieldMatrix can be confusing. For example

sigma = Stress(1.0, 2.0, 3.0,
               4.0, 5.0, 6.0,
               7.0, 8.0, 9.0)

3×3 Stress:
 1.0  4.0  7.0
 2.0  5.0  8.0
 3.0  6.0  9.0

will give you the transpose of what the multi-argument formatting suggests. For clarity, you may consider using the alternative

sigma = Stress(SA[1.0 2.0 3.0;
                  4.0 5.0 6.0;
                  7.0 8.0 9.0])

# StaticArraysCore.FieldVectorType

abstract FieldVector{N, T} <: FieldArray{Tuple{N}, 1}

Inheriting from this type will make it easy to create your own vector types. A FieldVector will automatically define getindex and setindex! appropriately. An immutable FieldVector will be as performant as an SVector of similar length and element type, while a mutable FieldVector will behave similarly to an MVector.

If you define a FieldVector which is parametric on the element type you should consider defining similar_type to preserve your array type through array operations as in the example below.

Example

struct Vec3D{T} <: FieldVector{3, T}
    x::T
    y::T
    z::T
end

StaticArrays.similar_type(::Type{<:Vec3D}, ::Type{T}, s::Size{(3,)}) where {T} = Vec3D{T}

# StaticArraysCore.MArrayType

MArray{S, T, N, L}(undef)
MArray{S, T, N, L}(x::NTuple{L})
MArray{S, T, N, L}(x1, x2, x3, ...)

Construct a statically-sized, mutable array MArray. The data may optionally be provided upon construction and can be mutated later. The S parameter is a Tuple-type specifying the dimensions, or size, of the array - such as Tuple{3,4,5} for a 3×4×5-sized array. The N parameter is the dimension of the array; the L parameter is the length of the array and is always equal to prod(S). Constructors may drop the L, N and T parameters if they are inferrable from the input (e.g. L is always inferrable from S).

MArray{S}(a::Array)

Construct a statically-sized, mutable array of dimensions S (expressed as a Tuple{...}) using the data from a. The S parameter is mandatory since the size of a is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.MMatrixType

MMatrix{S1, S2, T, L}(undef)
MMatrix{S1, S2, T, L}(x::NTuple{L, T})
MMatrix{S1, S2, T, L}(x1, x2, x3, ...)

Construct a statically-sized, mutable matrix MMatrix. The data may optionally be provided upon construction and can be mutated later. The L parameter is the length of the array and is always equal to S1 * S2. Constructors may drop the L, T and even S2 parameters if they are inferrable from the input (e.g. L is always inferrable from S1 and S2).

MMatrix{S1, S2}(mat::Matrix)

Construct a statically-sized, mutable matrix of dimensions S1 × S2 using the data from mat. The parameters S1 and S2 are mandatory since the size of mat is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.MVectorType

MVector{S,T}(undef)
MVector{S,T}(x::NTuple{S, T})
MVector{S,T}(x1, x2, x3, ...)

Construct a statically-sized, mutable vector MVector. Data may optionally be provided upon construction, and can be mutated later. Constructors may drop the T and S parameters if they are inferrable from the input (e.g. MVector(1,2,3) constructs an MVector{3, Int}).

MVector{S}(vec::Vector)

Construct a statically-sized, mutable vector of length S using the data from vec. The parameter S is mandatory since the length of vec is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.SArrayType

SArray{S, T, N, L}(x::NTuple{L})
SArray{S, T, N, L}(x1, x2, x3, ...)

Construct a statically-sized array SArray. Since this type is immutable, the data must be provided upon construction and cannot be mutated later. The S parameter is a Tuple-type specifying the dimensions, or size, of the array - such as Tuple{3,4,5} for a 3×4×5-sized array. The N parameter is the dimension of the array; the L parameter is the length of the array and is always equal to prod(S). Constructors may drop the L, N and T parameters if they are inferrable from the input (e.g. L is always inferrable from S).

SArray{S}(a::Array)

Construct a statically-sized array of dimensions S (expressed as a Tuple{...}) using the data from a. The S parameter is mandatory since the size of a is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.SMatrixType

SMatrix{S1, S2, T, L}(x::NTuple{L, T})
SMatrix{S1, S2, T, L}(x1, x2, x3, ...)

Construct a statically-sized matrix SMatrix. Since this type is immutable, the data must be provided upon construction and cannot be mutated later. The L parameter is the length of the array and is always equal to S1 * S2. Constructors may drop the L, T and even S2 parameters if they are inferrable from the input (e.g. L is always inferrable from S1 and S2).

SMatrix{S1, S2}(mat::Matrix)

Construct a statically-sized matrix of dimensions S1 × S2 using the data from mat. The parameters S1 and S2 are mandatory since the size of mat is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.SVectorType

SVector{S, T}(x::NTuple{S, T})
SVector{S, T}(x1, x2, x3, ...)

Construct a statically-sized vector SVector. Since this type is immutable, the data must be provided upon construction and cannot be mutated later. Constructors may drop the T and S parameters if they are inferrable from the input (e.g. SVector(1,2,3) constructs an SVector{3, Int}).

SVector{S}(vec::Vector)

Construct a statically-sized vector of length S using the data from vec. The parameter S is mandatory since the length of vec is unknown to the compiler (the element type may optionally also be specified).

# StaticArraysCore.SizeType

Size(dims::Int...)

Size is used extensively throughout the StaticArrays API to describe compile-time knowledge of the size of an array. The dimensions are stored as a type parameter and are statically propagated by the compiler, resulting in efficient, type-inferrable code. For example, to create a static matrix of zeros, use A = zeros(SMatrix{3,3}). The static size of A can be obtained by Size(A). (rather than size(zeros(3,3)), which returns Base.Tuple{2,Int}).

Note that if dimensions are not known statically (e.g., for standard Arrays), Dynamic() should be used instead of an Int.

Size(a::AbstractArray)
Size(::Type{T<:AbstractArray})

The Size constructor can be used to extract static dimension information from a given array. For example:

julia> Size(zeros(SMatrix{3, 4}))
Size(3, 4)

julia> Size(zeros(3, 4))
Size(StaticArrays.Dynamic(), StaticArrays.Dynamic())

This has multiple uses, including "trait"-based dispatch on the size of a statically-sized array. For example:

det(x::StaticMatrix) = _det(Size(x), x)
_det(::Size{(1,1)}, x::StaticMatrix) = x[1,1]
_det(::Size{(2,2)}, x::StaticMatrix) = x[1,1]*x[2,2] - x[1,2]*x[2,1]
# and other definitions as necessary

# StaticArraysCore.SizedArrayType

SizedArray{Tuple{dims...}}(array)

Wraps an AbstractArray with a static size, so to take advantage of the (faster) methods defined by the static array package. The size is checked once upon construction to determine if the number of elements (length) match, but the array may be reshaped.

The aliases SizedVector{N} and SizedMatrix{N,M} are provided as more convenient names for one and two dimensional SizedArrays. For example, to wrap a 2x3 array a in a SizedArray, use SizedMatrix{2,3}(a).

# StaticArraysCore.StaticArrayType

abstract type StaticArray{S, T, N} <: AbstractArray{T, N} end
StaticScalar{T}     = StaticArray{Tuple{}, T, 0}
StaticVector{N,T}   = StaticArray{Tuple{N}, T, 1}
StaticMatrix{N,M,T} = StaticArray{Tuple{N,M}, T, 2}

StaticArrays are Julia arrays with fixed, known size.

Dev docs

They must define the following methods:

  • Constructors that accept a flat tuple of data.

  • getindex() with an integer (linear indexing) (preferably @inline with @boundscheck).

  • Tuple(), returning the data in a flat Tuple.

It may be useful to implement:

  • similar_type(::Type{MyStaticArray}, ::Type{NewElType}, ::Size{NewSize}), returning a type (or type constructor) that accepts a flat tuple of data.

For mutable containers you may also need to define the following:

  • setindex! for a single element (linear indexing).

  • similar(::Type{MyStaticArray}, ::Type{NewElType}, ::Size{NewSize}).

  • In some cases, a zero-parameter constructor, MyStaticArray{...}() for unintialized data is assumed to exist.

(see also SVector, SMatrix, SArray, MVector, MMatrix, MArray, SizedArray, FieldVector, FieldMatrix and FieldArray)

# StaticArraysCore.similar_typeFunction

similar_type(static_array)
similar_type(static_array, T)
similar_type(array, ::Size)
similar_type(array, T, ::Size)

Returns a constructor for a statically-sized array similar to the input array (or type) static_array/array, optionally with different element type T or size Size. If the input array is not a StaticArray then the Size is mandatory.

This differs from similar() in that the resulting array type may not be mutable (or define setindex!()), and therefore the returned type may need to be constructed with its data.

Note that the (optional) size must be specified as a static Size object (so the compiler can infer the result statically).

New types should define the signature similar_type(::Type{A},::Type{T},::Size{S}) where {A<:MyType,T,S} if they wish to overload the default behavior.

# StaticArraysCore.size_to_tupleMethod

size_to_tuple(::Type{S}) where S<:Tuple

Converts a size given by Tuple{N, M, ...} into a tuple (N, M, ...).