Arrays

Basic functions

ndims(A::AbstractArray) → Integer

Returns the number of dimensions of A.

julia>A=ones(3,4,5);julia>ndims(A)3
size(A::AbstractArray[, dim...])

Returns a tuple containing the dimensions of A. Optionally you can specify the dimension(s) you want the length of, and get the length of that dimension, or a tuple of the lengths of dimensions you asked for.

julia>A=ones(2,3,4);julia>size(A,2)3julia>size(A,3,2)(4,3)
indices(A)

Returns the tuple of valid indices for array A.

indices(A, d)

Returns the valid range of indices for array A along dimension d.

length(A::AbstractArray) → Integer

Returns the number of elements in A.

julia>A=ones(3,4,5);julia>length(A)60
eachindex(A...)

Creates an iterable object for visiting each index of an AbstractArray A in an efficient manner. For array types that have opted into fast linear indexing (like Array), this is simply the range 1:length(A). For other array types, this returns a specialized Cartesian range to efficiently index into the array with indices specified for every dimension. For other iterables, including strings and dictionaries, this returns an iterator object supporting arbitrary index types (e.g. unevenly spaced or non-integer indices).

Example for a sparse 2-d array:

julia>A=sparse([1,1,2],[1,3,1],[1,2,-5])2×3sparsematrixwith3Int64nonzeroentries:[1,1]=1[2,1]=-5[1,3]=2julia>foriterineachindex(A)@showiter.I[1],iter.I[2]@showA[iter]end(iter.I[1],iter.I[2])=(1,1)A[iter]=1(iter.I[1],iter.I[2])=(2,1)A[iter]=-5(iter.I[1],iter.I[2])=(1,2)A[iter]=0(iter.I[1],iter.I[2])=(2,2)A[iter]=0(iter.I[1],iter.I[2])=(1,3)A[iter]=2(iter.I[1],iter.I[2])=(2,3)A[iter]=0

If you supply more than one AbstractArray argument, eachindex will create an iterable object that is fast for all arguments (a UnitRange if all inputs have fast linear indexing, a CartesianRange otherwise). If the arrays have different sizes and/or dimensionalities, eachindex returns an iterable that spans the largest range along each dimension.

linearindices(A)

Returns a UnitRange specifying the valid range of indices for A[i] where i is an Int. For arrays with conventional indexing (indices start at 1), or any multidimensional array, this is 1:length(A); however, for one-dimensional arrays with unconventional indices, this is indices(A,1).

Calling this function is the “safe” way to write algorithms that exploit linear indexing.

Base.linearindexing(A)

linearindexing defines how an AbstractArray most efficiently accesses its elements. If Base.linearindexing(A) returns Base.LinearFast(), this means that linear indexing with only one index is an efficient operation. If it instead returns Base.LinearSlow() (by default), this means that the array intrinsically accesses its elements with indices specified for every dimension. Since converting a linear index to multiple indexing subscripts is typically very expensive, this provides a traits-based mechanism to enable efficient generic code for all array types.

An abstract array subtype MyArray that wishes to opt into fast linear indexing behaviors should define linearindexing in the type-domain:

Base.linearindexing{T<:MyArray}(::Type{T})=Base.LinearFast()
countnz(A)

Counts the number of nonzero values in array A (dense or sparse). Note that this is not a constant-time operation. For sparse matrices, one should usually use nnz, which returns the number of stored values.

conj!(A)

Convert an array to its complex conjugate in-place.

stride(A, k::Integer)

Returns the distance in memory (in number of elements) between adjacent elements in dimension k.

julia>A=ones(3,4,5);julia>stride(A,2)3julia>stride(A,3)12
strides(A)

Returns a tuple of the memory strides in each dimension.

julia>A=ones(3,4,5);julia>strides(A)(1,3,12)
ind2sub(dims, index) → subscripts

Returns a tuple of subscripts into an array with dimensions dims, corresponding to the linear index index.

Example:

i,j,...=ind2sub(size(A),indmax(A))

provides the indices of the maximum element.

ind2sub(a, index) → subscripts

Returns a tuple of subscripts into array a corresponding to the linear index index.

sub2ind(dims, i, j, k...) → index

The inverse of ind2sub, returns the linear index corresponding to the provided subscripts.

LinAlg.checksquare(A)

Check that a matrix is square, then return its common dimension. For multiple arguments, return a vector.

Constructors

Array(dims)

Array{T}(dims) constructs an uninitialized dense array with element type T. dims may be a tuple or a series of integer arguments. The syntax Array(T,dims) is also available, but deprecated.

getindex(type[, elements...])

Construct a 1-d array of the specified type. This is usually called with the syntax Type[]. Element values can be specified using Type[a,b,c,...].

zeros(type, dims)

Create an array of all zeros of specified type. The type defaults to Float64 if not specified.

zeros(A)

Create an array of all zeros with the same element type and shape as A.

ones(type, dims)

Create an array of all ones of specified type. The type defaults to Float64 if not specified.

ones(A)

Create an array of all ones with the same element type and shape as A.

trues(dims)

Create a BitArray with all values set to true.

trues(A)

Create a BitArray with all values set to true of the same shape as A.

falses(dims)

Create a BitArray with all values set to false.

falses(A)

Create a BitArray with all values set to false of the same shape as A.

fill(x, dims)

Create an array filled with the value x. For example, fill(1.0,(10,10)) returns a 10×10 array of floats, with each element initialized to 1.0.

If x is an object reference, all elements will refer to the same object. fill(Foo(),dims) will return an array filled with the result of evaluating Foo() once.

fill!(A, x)

Fill array A with the value x. If x is an object reference, all elements will refer to the same object. fill!(A,Foo()) will return A filled with the result of evaluating Foo() once.

reshape(A, dims)

Create an array with the same data as the given array, but with different dimensions.

similar(array[, element_type=eltype(array)][, dims=size(array)])

Create an uninitialized mutable array with the given element type and size, based upon the given source array. The second and third arguments are both optional, defaulting to the given array’s eltype and size. The dimensions may be specified either as a single tuple argument or as a series of integer arguments.

Custom AbstractArray subtypes may choose which specific array type is best-suited to return for the given element type and dimensionality. If they do not specialize this method, the default is an Array{element_type}(dims...).

For example, similar(1:10,1,4) returns an uninitialized Array{Int,2} since ranges are neither mutable nor support 2 dimensions:

julia>similar(1:10,1,4)1×4Array{Int64,2}:4419743872437441387244197438880

Conversely, similar(trues(10,10),2) returns an uninitialized BitVector with two elements since BitArrays are both mutable and can support 1-dimensional arrays:

julia>similar(trues(10,10),2)2-elementBitArray{1}:falsefalse

Since BitArrays can only store elements of type Bool, however, if you request a different element type it will create a regular Array instead:

julia>similar(falses(10),Float64,2,4)2×4Array{Float64,2}:2.18425e-3142.18425e-3142.18425e-3142.18425e-3142.18425e-3142.18425e-3142.18425e-3142.18425e-314
similar(storagetype, indices)

Create an uninitialized mutable array analogous to that specified by storagetype, but with indices specified by the last argument. storagetype might be a type or a function.

Examples:

similar(Array{Int},indices(A))

creates an array that “acts like” an Array{Int} (and might indeed be backed by one), but which is indexed identically to A. If A has conventional indexing, this will be identical to Array{Int}(size(A)), but if A has unconventional indexing then the indices of the result will match A.

similar(BitArray,(indices(A,2),))

would create a 1-dimensional logical array whose indices match those of the columns of A.

similar(dims->zeros(Int,dims),indices(A))

would create an array of Int, initialized to zero, matching the indices of A.

reinterpret(type, A)

Change the type-interpretation of a block of memory. For example, reinterpret(Float32,UInt32(7)) interprets the 4 bytes corresponding to UInt32(7) as a Float32. For arrays, this constructs an array with the same binary data as the given array, but with the specified element type.

eye([T::Type=Float64, ]n::Integer)

n-by-n identity matrix. The default element type is Float64.

eye([T::Type=Float64, ]m::Integer, n::Integer)

m-by-n identity matrix. The default element type is Float64.

eye(A)

Constructs an identity matrix of the same dimensions and type as A.

julia>A=[123;456;789]3×3Array{Int64,2}:123456789julia>eye(A)3×3Array{Int64,2}:100010001

Note the difference from ones().

linspace(start, stop, n=50)

Construct a range of n linearly spaced elements from start to stop.

logspace(start, stop, n=50)

Construct a vector of n logarithmically spaced numbers from 10^start to 10^stop.

Mathematical operators and functions

All mathematical operations and functions are supported for arrays

broadcast(f, As...)

Broadcasts the arrays As to a common size by expanding singleton dimensions, and returns an array of the results f(as...) for each position.

broadcast!(f, dest, As...)

Like broadcast, but store the result of broadcast(f,As...) in the dest array. Note that dest is only used to store the result, and does not supply arguments to f unless it is also listed in the As, as in broadcast!(f,A,A,B) to perform A[:]=broadcast(f,A,B).

bitbroadcast(f, As...)

Like broadcast, but allocates a BitArray to store the result, rather then an Array.

Indexing, Assignment, and Concatenation

getindex(A, inds...)

Returns a subset of array A as specified by inds, where each ind may be an Int, a Range, or a Vector. See the manual section on array indexing for details.

view(A, inds...)

Like getindex(), but returns a view into the parent array A with the given indices instead of making a copy. Calling getindex() or setindex!() on the returned SubArray computes the indices to the parent array on the fly without checking bounds.

@view A[inds...]

Creates a SubArray from an indexing expression. This can only be applied directly to a reference expression (e.g. @viewA[1,2:end]), and should not be used as the target of an assignment (e.g. @view(A[1,2:end])=...).

parent(A)

Returns the “parent array” of an array view type (e.g., SubArray), or the array itself if it is not a view.

parentindexes(A)

From an array view A, returns the corresponding indexes in the parent.

slicedim(A, d, i)

Return all the data of A where the index for dimension d equals i. Equivalent to A[:,:,...,i,:,:,...] where i is in position d.

setindex!(A, X, inds...)

Store values from array X within some subset of A as specified by inds.

broadcast_getindex(A, inds...)

Broadcasts the inds arrays to a common size like broadcast, and returns an array of the results A[ks...], where ks goes over the positions in the broadcast.

broadcast_setindex!(A, X, inds...)

Broadcasts the X and inds arrays to a common size and stores the value from each position in X at the indices given by the same positions in inds.

isassigned(array, i) → Bool

Tests whether the given array has a value associated with index i. Returns false if the index is out of bounds, or has an undefined reference.

cat(dims, A...)

Concatenate the input arrays along the specified dimensions in the iterable dims. For dimensions not in dims, all input arrays should have the same size, which will also be the size of the output array along that dimension. For dimensions in dims, the size of the output array is the sum of the sizes of the input arrays along that dimension. If dims is a single number, the different arrays are tightly stacked along that dimension. If dims is an iterable containing several dimensions, this allows one to construct block diagonal matrices and their higher-dimensional analogues by simultaneously increasing several dimensions for every new input array and putting zero blocks elsewhere. For example, cat([1,2],matrices...) builds a block diagonal matrix, i.e. a block matrix with matrices[1], matrices[2], ... as diagonal blocks and matching zero blocks away from the diagonal.

vcat(A...)

Concatenate along dimension 1.

julia>a=[12345]1×5Array{Int64,2}:12345julia>b=[678910;1112131415]2×5Array{Int64,2}:6789101112131415julia>vcat(a,b)3×5Array{Int64,2}:123456789101112131415julia>c=([123],[456])([123],[456])julia>vcat(c...)2×3Array{Int64,2}:123456
hcat(A...)

Concatenate along dimension 2.

julia>a=[1;2;3;4;5]5-elementArray{Int64,1}:12345julia>b=[67;89;1011;1213;1415]5×2Array{Int64,2}:6789101112131415julia>hcat(a,b)5×3Array{Int64,2}:167289310114121351415julia>c=([1;2;3],[4;5;6])([1,2,3],[4,5,6])julia>hcat(c...)3×2Array{Int64,2}:142536
hvcat(rows::Tuple{Vararg{Int}}, values...)

Horizontal and vertical concatenation in one call. This function is called for block matrix syntax. The first argument specifies the number of arguments to concatenate in each block row.

julia>a,b,c,d,e,f=1,2,3,4,5,6(1,2,3,4,5,6)julia>[abc;def]2×3Array{Int64,2}:123456julia>hvcat((3,3),a,b,c,d,e,f)2×3Array{Int64,2}:123456julia>[ab;cd;ef]3×2Array{Int64,2}:123456julia>hvcat((2,2,2),a,b,c,d,e,f)3×2Array{Int64,2}:123456

If the first argument is a single integer n, then all block rows are assumed to have n block columns.

flipdim(A, d)

Reverse A in dimension d.

julia>b=[12;34]2×2Array{Int64,2}:1234julia>flipdim(b,2)2×2Array{Int64,2}:2143
circshift(A, shifts)

Circularly shift the data in an array. The second argument is a vector giving the amount to shift in each dimension.

julia>b=reshape(collect(1:16),(4,4))4×4Array{Int64,2}:15913261014371115481216julia>circshift(b,(0,2))4×4Array{Int64,2}:91315101426111537121648julia>circshift(b,(-1,0))4×4Array{Int64,2}:26101437111548121615913
find(A)

Return a vector of the linear indexes of the non-zeros in A (determined by A[i]!=0). A common use of this is to convert a boolean array to an array of indexes of the true elements. If there are no non-zero elements of A, find returns an empty array.

julia>A=[truefalse;falsetrue]2×2Array{Bool,2}:truefalsefalsetruejulia>find(A)2-elementArray{Int64,1}:14
find(f::Function, A)

Return a vector I of the linear indexes of A where f(A[I]) returns true. If there are no such elements of A, find returns an empty array.

julia>A=[12;34]2×2Array{Int64,2}:1234julia>find(isodd,A)2-elementArray{Int64,1}:12
findn(A)

Return a vector of indexes for each dimension giving the locations of the non-zeros in A (determined by A[i]!=0). If there are no non-zero elements of A, findn returns a 2-tuple of empty arrays.

julia>A=[120;003;040]3×3Array{Int64,2}:120003040julia>findn(A)([1,1,3,2],[1,2,2,3])julia>A=zeros(2,2)2×2Array{Float64,2}:0.00.00.00.0julia>findn(A)(Int64[],Int64[])
findnz(A)

Return a tuple (I,J,V) where I and J are the row and column indexes of the non-zero values in matrix A, and V is a vector of the non-zero values.

julia>A=[120;003;040]3×3Array{Int64,2}:120003040julia>findnz(A)([1,1,3,2],[1,2,2,3],[1,2,4,3])
findfirst(A)

Return the linear index of the first non-zero value in A (determined by A[i]!=0). Returns 0 if no such value is found.

julia>A=[00;10]2×2Array{Int64,2}:0010julia>findfirst(A)2
findfirst(A, v)

Return the linear index of the first element equal to v in A. Returns 0 if v is not found.

julia>A=[46;22]2×2Array{Int64,2}:4622julia>findfirst(A,2)2julia>findfirst(A,3)0
findfirst(predicate::Function, A)

Return the linear index of the first element of A for which predicate returns true. Returns 0 if there is no such element.

julia>A=[14;22]2×2Array{Int64,2}:1422julia>findfirst(iseven,A)2julia>findfirst(x->x>10,A)0
findlast(A)

Return the linear index of the last non-zero value in A (determined by A[i]!=0). Returns 0 if there is no non-zero value in A.

julia>A=[10;10]2×2Array{Int64,2}:1010julia>findlast(A)2julia>A=zeros(2,2)2×2Array{Float64,2}:0.00.00.00.0julia>findlast(A)0
findlast(A, v)

Return the linear index of the last element equal to v in A. Returns 0 if there is no element of A equal to v.

julia>A=[12;21]2×2Array{Int64,2}:1221julia>findlast(A,1)4julia>findlast(A,2)3julia>findlast(A,3)0
findlast(predicate::Function, A)

Return the linear index of the last element of A for which predicate returns true. Returns 0 if there is no such element.

julia>A=[12;34]2×2Array{Int64,2}:1234julia>findlast(isodd,A)2julia>findlast(x->x>5,A)0
findnext(A, i::Integer)

Find the next linear index >= i of a non-zero element of A, or 0 if not found.

julia>A=[00;10]2×2Array{Int64,2}:0010julia>findnext(A,1)2julia>findnext(A,3)0
findnext(predicate::Function, A, i::Integer)

Find the next linear index >= i of an element of A for which predicate returns true, or 0 if not found.

julia>A=[14;22]2×2Array{Int64,2}:1422julia>findnext(isodd,A,1)1julia>findnext(isodd,A,2)0
findnext(A, v, i::Integer)

Find the next linear index >= i of an element of A equal to v (using ==), or 0 if not found.

julia>A=[14;22]2×2Array{Int64,2}:1422julia>findnext(A,4,4)0julia>findnext(A,4,3)3
findprev(A, i::Integer)

Find the previous linear index <= i of a non-zero element of A, or 0 if not found.

julia>A=[00;12]2×2Array{Int64,2}:0012julia>findprev(A,2)2julia>findprev(A,1)0
findprev(predicate::Function, A, i::Integer)

Find the previous linear index <= i of an element of A for which predicate returns true, or 0 if not found.

julia>A=[46;12]2×2Array{Int64,2}:4612julia>findprev(isodd,A,1)0julia>findprev(isodd,A,3)2
findprev(A, v, i::Integer)

Find the previous linear index <= i of an element of A equal to v (using ==), or 0 if not found.

julia>A=[00;12]2×2Array{Int64,2}:0012julia>findprev(A,1,4)2julia>findprev(A,1,1)0
permutedims(A, perm)

Permute the dimensions of array A. perm is a vector specifying a permutation of length ndims(A). This is a generalization of transpose for multi-dimensional arrays. Transpose is equivalent to permutedims(A,[2,1]).

ipermutedims(A, perm)

Like permutedims(), except the inverse of the given permutation is applied.

permutedims!(dest, src, perm)

Permute the dimensions of array src and store the result in the array dest. perm is a vector specifying a permutation of length ndims(src). The preallocated array dest should have size(dest)==size(src)[perm] and is completely overwritten. No in-place permutation is supported and unexpected results will happen if src and dest have overlapping memory regions.

squeeze(A, dims)

Remove the dimensions specified by dims from array A. Elements of dims must be unique and within the range 1:ndims(A). size(A,i) must equal 1 for all i in dims.

julia>a=reshape(collect(1:4),(2,2,1,1))2×2×1×1Array{Int64,4}:[:,:,1,1]=1324julia>squeeze(a,3)2×2×1Array{Int64,3}:[:,:,1]=1324
vec(a::AbstractArray) → Vector

Reshape array a as a one-dimensional column vector.

julia>a=[123;456]2×3Array{Int64,2}:123456julia>vec(a)6-elementArray{Int64,1}:142536
promote_shape(s1, s2)

Check two array shapes for compatibility, allowing trailing singleton dimensions, and return whichever shape has more dimensions.

checkbounds(A, I...)

Throw an error if the specified indices I are not in bounds for the given array A.

checkbounds(Bool, A, I...)

Return true if the specified indices I are in bounds for the given array A. Subtypes of AbstractArray should specialize this method if they need to provide custom bounds checking behaviors; however, in many cases one can rely on A‘s indices and checkindex.

See also checkindex.

checkindex(Bool, inds::AbstractUnitRange, index)

Return true if the given index is within the bounds of inds. Custom types that would like to behave as indices for all arrays can extend this method in order to provide a specialized bounds checking implementation.

randsubseq(A, p) → Vector

Return a vector consisting of a random subsequence of the given array A, where each element of A is included (in order) with independent probability p. (Complexity is linear in p*length(A), so this function is efficient even if p is small and A is large.) Technically, this process is known as “Bernoulli sampling” of A.

randsubseq!(S, A, p)

Like randsubseq, but the results are stored in S (which is resized as needed).

Array functions

cumprod(A, dim=1)

Cumulative product along a dimension dim (defaults to 1). See also cumprod!() to use a preallocated output array, both for performance and to control the precision of the output (e.g. to avoid overflow).

julia>a=[123;456]2×3Array{Int64,2}:123456julia>cumprod(a,1)2×3Array{Int64,2}:12341018julia>cumprod(a,2)2×3Array{Int64,2}:126420120
cumprod!(B, A[, dim])

Cumulative product of A along a dimension, storing the result in B. The dimension defaults to 1.

cumsum(A, dim=1)

Cumulative sum along a dimension dim (defaults to 1). See also cumsum!() to use a preallocated output array, both for performance and to control the precision of the output (e.g. to avoid overflow).

julia>a=[123;456]2×3Array{Int64,2}:123456julia>cumsum(a,1)2×3Array{Int64,2}:123579julia>cumsum(a,2)2×3Array{Int64,2}:1364915
cumsum!(B, A[, dim])

Cumulative sum of A along a dimension, storing the result in B. The dimension defaults to 1.

cumsum_kbn(A[, dim])

Cumulative sum along a dimension, using the Kahan-Babuska-Neumaier compensated summation algorithm for additional accuracy. The dimension defaults to 1.

cummin(A[, dim])

Cumulative minimum along a dimension. The dimension defaults to 1.

cummax(A[, dim])

Cumulative maximum along a dimension. The dimension defaults to 1.

diff(A[, dim])

Finite difference operator of matrix or vector.

gradient(F[, h])

Compute differences along vector F, using h as the spacing between points. The default spacing is one.

rot180(A)

Rotate matrix A 180 degrees.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rot180(a)2×2Array{Int64,2}:4321
rot180(A, k)

Rotate matrix A 180 degrees an integer k number of times. If k is even, this is equivalent to a copy.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rot180(a,1)2×2Array{Int64,2}:4321julia>rot180(a,2)2×2Array{Int64,2}:1234
rotl90(A)

Rotate matrix A left 90 degrees.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rotl90(a)2×2Array{Int64,2}:2413
rotl90(A, k)

Rotate matrix A left 90 degrees an integer k number of times. If k is zero or a multiple of four, this is equivalent to a copy.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rotl90(a,1)2×2Array{Int64,2}:2413julia>rotl90(a,2)2×2Array{Int64,2}:4321julia>rotl90(a,3)2×2Array{Int64,2}:3142julia>rotl90(a,4)2×2Array{Int64,2}:1234
rotr90(A)

Rotate matrix A right 90 degrees.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rotr90(a)2×2Array{Int64,2}:3142
rotr90(A, k)

Rotate matrix A right 90 degrees an integer k number of times. If k is zero or a multiple of four, this is equivalent to a copy.

julia>a=[12;34]2×2Array{Int64,2}:1234julia>rotr90(a,1)2×2Array{Int64,2}:3142julia>rotr90(a,2)2×2Array{Int64,2}:4321julia>rotr90(a,3)2×2Array{Int64,2}:2413julia>rotr90(a,4)2×2Array{Int64,2}:1234
reducedim(f, A, region[, v0])

Reduce 2-argument function f along dimensions of A. region is a vector specifying the dimensions to reduce, and v0 is the initial value to use in the reductions. For +, *, max and min the v0 argument is optional.

The associativity of the reduction is implementation-dependent; if you need a particular associativity, e.g. left-to-right, you should write your own loop. See documentation for reduce().

julia>a=reshape(collect(1:16),(4,4))4×4Array{Int64,2}:15913261014371115481216julia>reducedim(max,a,2)4×1Array{Int64,2}:13141516julia>reducedim(max,a,1)1×4Array{Int64,2}:481216
mapreducedim(f, op, A, region[, v0])

Evaluates to the same as reducedim(op,map(f,A),region,f(v0)), but is generally faster because the intermediate array is avoided.

julia>a=reshape(collect(1:16),(4,4))4×4Array{Int64,2}:15913261014371115481216julia>mapreducedim(isodd,*,a,1)1×4Array{Bool,2}:falsefalsefalsefalsejulia>mapreducedim(isodd,|,a,1,true)1×4Array{Bool,2}:truetruetruetrue
mapslices(f, A, dims)

Transform the given dimensions of array A using function f. f is called on each slice of A of the form A[...,:,...,:,...]. dims is an integer vector specifying where the colons go in this expression. The results are concatenated along the remaining dimensions. For example, if dims is [1,2] and A is 4-dimensional, f is called on A[:,:,i,j] for all i and j.

julia>a=reshape(collect(1:16),(2,2,2,2))2×2×2×2Array{Int64,4}:[:,:,1,1]=1324[:,:,2,1]=5768[:,:,1,2]=9111012[:,:,2,2]=13151416julia>mapslices(sum,a,[1,2])1×1×2×2Array{Int64,4}:[:,:,1,1]=10[:,:,2,1]=26[:,:,1,2]=42[:,:,2,2]=58
sum_kbn(A)

Returns the sum of all array elements, using the Kahan-Babuska-Neumaier compensated summation algorithm for additional accuracy.

Combinatorics

randperm([rng, ]n)

Construct a random permutation of length n. The optional rng argument specifies a random number generator (see Random Numbers). To randomly permute a arbitrary vector, see shuffle() or shuffle!().

invperm(v)

Return the inverse permutation of v

isperm(v) → Bool

Returns true if v is a valid permutation.

permute!(v, p)

Permute vector v in-place, according to permutation p. No checking is done to verify that p is a permutation.

To return a new permutation, use v[p]. Note that this is generally faster than permute!(v,p) for large vectors.

ipermute!(v, p)

Like permute!, but the inverse of the given permutation is applied.

randcycle([rng, ]n)

Construct a random cyclic permutation of length n. The optional rng argument specifies a random number generator, see Random Numbers.

shuffle([rng, ]v)

Return a randomly permuted copy of v. The optional rng argument specifies a random number generator (see Random Numbers). To permute v in-place, see shuffle!(). To obtain randomly permuted indices, see randperm().

shuffle!([rng, ]v)

In-place version of shuffle(): randomly permute the array v in-place, optionally supplying the random-number generator rng.

reverse(v[, start=1[, stop=length(v)]])

Return a copy of v reversed from start to stop.

reverseind(v, i)

Given an index i in reverse(v), return the corresponding index in v so that v[reverseind(v,i)]==reverse(v)[i]. (This can be nontrivial in the case where v is a Unicode string.)

reverse!(v[, start=1[, stop=length(v)]]) → v

In-place version of reverse().

BitArrays

BitArrays are space-efficient “packed” boolean arrays, which store one bit per boolean value. They can be used similarly to Array{Bool} arrays (which store one byte per boolean value), and can be converted to/from the latter via Array(bitarray) and BitArray(array), respectively.

flipbits!(B::BitArray{N}) → BitArray{N}

Performs a bitwise not operation on B. See ~ operator.

julia>A=trues(2,2)2×2BitArray{2}:truetruetruetruejulia>flipbits!(A)2×2BitArray{2}:falsefalsefalsefalse
rol!(dest::BitVector, src::BitVector, i::Integer) → BitVector

Performs a left rotation operation on src and puts the result into dest. i controls how far to rotate the bits.

rol!(B::BitVector, i::Integer) → BitVector

Performs a left rotation operation in-place on B. i controls how far to rotate the bits.

rol(B::BitVector, i::Integer) → BitVector

Performs a left rotation operation, returning a new BitVector. i controls how far to rotate the bits. See also rol!().

julia>A=BitArray([true,true,false,false,true])5-elementBitArray{1}:truetruefalsefalsetruejulia>rol(A,1)5-elementBitArray{1}:truefalsefalsetruetruejulia>rol(A,2)5-elementBitArray{1}:falsefalsetruetruetruejulia>rol(A,5)5-elementBitArray{1}:truetruefalsefalsetrue
ror!(dest::BitVector, src::BitVector, i::Integer) → BitVector

Performs a right rotation operation on src and puts the result into dest. i controls how far to rotate the bits.

ror!(B::BitVector, i::Integer) → BitVector

Performs a right rotation operation in-place on B. i controls how far to rotate the bits.

ror(B::BitVector, i::Integer) → BitVector

Performs a right rotation operation on B, returning a new BitVector. i controls how far to rotate the bits. See also ror!().

julia>A=BitArray([true,true,false,false,true])5-elementBitArray{1}:truetruefalsefalsetruejulia>ror(A,1)5-elementBitArray{1}:truetruetruefalsefalsejulia>ror(A,2)5-elementBitArray{1}:falsetruetruetruefalsejulia>ror(A,5)5-elementBitArray{1}:truetruefalsefalsetrue

Sparse Vectors and Matrices

Sparse vectors and matrices largely support the same set of operations as their dense counterparts. The following functions are specific to sparse arrays.

sparse(I, J, V[, m, n, combine])

Create a sparse matrix S of dimensions mxn such that S[I[k],J[k]]=V[k]. The combine function is used to combine duplicates. If m and n are not specified, they are set to maximum(I) and maximum(J) respectively. If the combine function is not supplied, combine defaults to + unless the elements of V are Booleans in which case combine defaults to |. All elements of I must satisfy 1<=I[k]<=m, and all elements of J must satisfy 1<=J[k]<=n. Numerical zeros in (I, J, V) are retained as structural nonzeros; to drop numerical zeros, use dropzeros!().

For additional documentation and an expert driver, see Base.SparseArrays.sparse!.

sparsevec(I, V[, m, combine])

Create a sparse vector S of length m such that S[I[k]]=V[k]. Duplicates are combined using the combine function, which defaults to + if no combine argument is provided, unless the elements of V are Booleans in which case combine defaults to |.

sparsevec(D::Dict[, m])

Create a sparse vector of length m where the nonzero indices are keys from the dictionary, and the nonzero values are the values from the dictionary.

issparse(S)

Returns true if S is sparse, and false otherwise.

sparse(A)

Convert an AbstractMatrix A into a sparse matrix.

sparsevec(A)

Convert a vector A into a sparse vector of length m.

full(S)

Convert a sparse matrix or vector S into a dense matrix or vector.

nnz(A)

Returns the number of stored (filled) elements in a sparse array.

spzeros([type, ]m[, n])

Create a sparse vector of length m or sparse matrix of size mxn. This sparse array will not contain any nonzero values. No storage will be allocated for nonzero values during construction. The type defaults to Float64 if not specified.

spones(S)

Create a sparse array with the same structure as that of S, but with every nonzero element having the value 1.0.

julia>A=sparse([1,2,3,4],[2,4,3,1],[5.,4.,3.,2.])4×4sparsematrixwith4Float64nonzeroentries:[4,1]=2.0[1,2]=5.0[3,3]=3.0[2,4]=4.0julia>spones(A)4×4sparsematrixwith4Float64nonzeroentries:[4,1]=1.0[1,2]=1.0[3,3]=1.0[2,4]=1.0

Note the difference from speye().

speye([type, ]m[, n])

Create a sparse identity matrix of size mxm. When n is supplied, create a sparse identity matrix of size mxn. The type defaults to Float64 if not specified.

speye(S)

Create a sparse identity matrix with the same size as S.

julia>A=sparse([1,2,3,4],[2,4,3,1],[5.,4.,3.,2.])4×4sparsematrixwith4Float64nonzeroentries:[4,1]=2.0[1,2]=5.0[3,3]=3.0[2,4]=4.0julia>speye(A)4×4sparsematrixwith4Float64nonzeroentries:[1,1]=1.0[2,2]=1.0[3,3]=1.0[4,4]=1.0

Note the difference from spones().

spdiagm(B, d[, m, n])

Construct a sparse diagonal matrix. B is a tuple of vectors containing the diagonals and d is a tuple containing the positions of the diagonals. In the case the input contains only one diagonal, B can be a vector (instead of a tuple) and d can be the diagonal position (instead of a tuple), defaulting to 0 (diagonal). Optionally, m and n specify the size of the resulting sparse matrix.

julia>spdiagm(([1,2,3,4],[4,3,2,1]),(-1,1))5×5sparsematrixwith8Int64nonzeroentries:[2,1]=1[1,2]=4[3,2]=2[2,3]=3[4,3]=3[3,4]=2[5,4]=4[4,5]=1
sprand([rng, ][type, ]m, [n, ]p::AbstractFloat[, rfn])

Create a random length m sparse vector or m by n sparse matrix, in which the probability of any element being nonzero is independently given by p (and hence the mean density of nonzeros is also exactly p). Nonzero values are sampled from the distribution specified by rfn and have the type type. The uniform distribution is used in case rfn is not specified. The optional rng argument specifies a random number generator, see Random Numbers.

sprandn([rng, ]m, [n, ]p::AbstractFloat)

Create a random sparse vector of length m or sparse matrix of size m by n with the specified (independent) probability p of any entry being nonzero, where nonzero values are sampled from the normal distribution. The optional rng argument specifies a random number generator, see Random Numbers.

nonzeros(A)

Return a vector of the structural nonzero values in sparse array A. This includes zeros that are explicitly stored in the sparse array. The returned vector points directly to the internal nonzero storage of A, and any modifications to the returned vector will mutate A as well. See rowvals() and nzrange().

rowvals(A::SparseMatrixCSC)

Return a vector of the row indices of A. Any modifications to the returned vector will mutate A as well. Providing access to how the row indices are stored internally can be useful in conjunction with iterating over structural nonzero values. See also nonzeros() and nzrange().

nzrange(A::SparseMatrixCSC, col)

Return the range of indices to the structural nonzero values of a sparse matrix column. In conjunction with nonzeros() and rowvals(), this allows for convenient iterating over a sparse matrix :

A=sparse(I,J,V)rows=rowvals(A)vals=nonzeros(A)m,n=size(A)fori=1:nforjinnzrange(A,i)row=rows[j]val=vals[j]# perform sparse wizardry...endend
dropzeros!(A::SparseMatrixCSC, trim::Bool = true)

Removes stored numerical zeros from A, optionally trimming resulting excess space from A.rowval and A.nzval when trim is true.

For an out-of-place version, see dropzeros(). For algorithmic information, see Base.SparseArrays.fkeep!().

dropzeros(A::SparseMatrixCSC, trim::Bool = true)

Generates a copy of A and removes stored numerical zeros from that copy, optionally trimming excess space from the result’s rowval and nzval arrays when trim is true.

For an in-place version and algorithmic information, see dropzeros!().

dropzeros!(x::SparseVector, trim::Bool = true)

Removes stored numerical zeros from x, optionally trimming resulting excess space from x.nzind and x.nzval when trim is true.

For an out-of-place version, see Base.SparseArrays.dropzeros(). For algorithmic information, see Base.SparseArrays.fkeep!().

dropzeros(x::SparseVector, trim::Bool = true)

Generates a copy of x and removes numerical zeros from that copy, optionally trimming excess space from the result’s nzind and nzval arrays when trim is true.

For an in-place version and algorithmic information, see Base.SparseArrays.dropzeros!().

permute{Tv,Ti,Tp<:Integer,Tq<:Integer}(A::SparseMatrixCSC{Tv,Ti}, p::AbstractVector{Tp},
q::AbstractVector{Tq})()

Bilaterally permute A, returning PAQ (A[p,q]). Column-permutation q‘s length must match A‘s column count (length(q)==A.n). Row-permutation p‘s length must match A‘s row count (length(p)==A.m).

For expert drivers and additional information, see Base.SparseArrays.permute!().

permute!{Tv,Ti,Tp<:Integer,Tq<:Integer}(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti},
p::AbstractVector{Tp}, q::AbstractVector{Tq}[, C::SparseMatrixCSC{Tv,Ti}])

Bilaterally permute A, storing result PAQ (A[p,q]) in X. Stores intermediate result (AQ)^T (transpose(A[:,q])) in optional argument C if present. Requires that none of X, A, and, if present, C alias each other; to store result PAQ back into A, use the following method lacking X:

permute!{Tv,Ti,Tp<:Integer,Tq<:Integer}(A::SparseMatrixCSC{Tv,Ti},p::AbstractVector{Tp},q::AbstractVector{Tq}[,C::SparseMatrixCSC{Tv,Ti}[,workcolptr::Vector{Ti}]])

X‘s dimensions must match those of A (X.m==A.m and X.n==A.n), and X must have enough storage to accommodate all allocated entries in A (length(X.rowval)>=nnz(A) and length(X.nzval)>=nnz(A)). Column-permutation q‘s length must match A‘s column count (length(q)==A.n). Row-permutation p‘s length must match A‘s row count (length(p)==A.m).

C‘s dimensions must match those of transpose(A) (C.m==A.n and C.n==A.m), and C must have enough storage to accommodate all allocated entries in A (length(C.rowval) >= nnz(A)``and``length(C.nzval) >= nnz(A)`).

For additional (algorithmic) information, and for versions of these methods that forgo argument checking, see (unexported) parent methods Base.SparseArrays.unchecked_noalias_permute!() and Base.SparseArrays.unchecked_aliasing_permute!().

See also: Base.SparseArrays.permute()