Gridap.Arrays

Gridap.ArraysModule

This module provides:

  • An extension of the AbstractArray interface in order to properly deal with mutable caches.
  • A mechanism to generate lazy arrays resulting from operations between arrays.
  • A collection of concrete implementations of AbstractArray.

The exported names in this module are:

source
Gridap.Arrays.ArrayBlockType
struct ArrayBlock{A,N}
  array::Array{A,N}
  touched::Array{Bool,N}
end

Block-wise storage of arrays, where touched indicates which blocks are active. Accessing a non-touched block returns nothing.

ArrayBlock is mostly used for multi-field and skeleton computations and assembly, where each block corresponds to a field or plus/minus respectively. Blocks might be nested in the case of multi-field computations on skeletons.

source
Gridap.Arrays.BlockMapType
struct BlockMap{N} <: Map
    size::NTuple{N,Int}
    indices::Vector{CartesianIndex{N}}
end

A BlockMap maps M = length(indices) arrays to a single array of N-dimensinal blocks, where only the blocks indexed by indices are touched and contain the corresponding entry of the input arrays.

Constructors:

BlockMap(l::Integer,i::Integer) ≡ BlockMap((l,),[CartesianIndex((i,))])
BlockMap(l::Integer,inds::Vector{<:Integer}) ≡ BlockMap((l,),[CartesianIndex((i,)) for i in inds])
BlockMap(s::NTuple,i::Integer) ≡ BlockMap(s,[CartesianIndex(i)])
BlockMap(s::NTuple,inds::Vector{<:NTuple}) ≡ BlockMap(s,[CartesianIndex(i) for i in inds])

Usage:

  lazy_map(BlockMap(2,[1,2]),a,b)
source
Gridap.Arrays.BroadcastingType
Broadcasting(f)

Returns a mapping that represents the "broadcasted" version of the function f.

Example

using Gridap.Arrays

a = [3,2]
b = [2,1]

bm = Broadcasting(+)

c = evaluate(bm,a,b)

println(c)

# output
[5, 3]
source
Gridap.Arrays.CachedArrayType
mutable struct CachedArray{T, N, A<:AbstractArray{T, N}} <: AbstractArray{T, N}

Type providing a re-sizable array.

The size of a CachedArray is changed via the setsize! function.

A CachedArray can be build with the constructors

using Gridap.Arrays
# Create an empty CachedArray
a = CachedArray(Float64,2)
# Resize to new shape (2,3)
setsize!(a,(2,3))
size(a)
# output
(2, 3)
source
Gridap.Arrays.CompressedArrayType
struct CompressedArray{T,N,A,P} <: AbstractArray{T,N}
  values::A
  ptrs::P
end

Type representing an array with a reduced set of values. The array is represented by a short array of values, namely the field values, and a large array of indices, namely the field ptrs. The i-th component of the resulting array is defined as values[ptrs[i]]. The type parameters A, and P are restricted to be array types by the inner constructor of this struct.

source
Gridap.Arrays.CompressedArrayMethod
CompressedArray(values::AbstractArray,ptrs::AbstractArray)

Creates a CompressedArray object by the given arrays of values and ptrs.

source
Gridap.Arrays.LazyArrayType

Subtype of AbstractArray which is the result of lazy_map. It represents the result of lazy_mapping a Map to a set of arrays that contain the mapping arguments. This struct makes use of the cache provided by the mapping in order to compute its indices (thus allowing to prevent allocation). The array is lazy, i.e., the values are only computed on demand. It extends the AbstractArray API with two methods:

array_cache(a::AbstractArray) getindex!(cache,a::AbstractArray,i...)

source
Gridap.Arrays.LocalPenaltySolveMapType
LocalPenaltySolveMap(; factorize! = lu!, pivot = NoPivot())

A map for solving local constrained linear systems, relying on a factorization method.

Given a left-hand-side 2x2 block matrix matrixmat and a set of 2xN right-hand-side arrays lhs, returns an N-Tuple of arrays containing the solutions to the linear systems.

Each system is given by A*[x_i; λ_i] = b_i, where A = [App, Aλp; Apλ, 0] is the augmented matrix, and b_i = [Bp; Bλ] is the right-hand side vector. The solution is computed using a penalty method, as x_i = ldiv!(factorize!(C,pivot),d_i) with C = App + μT * Apλ * Aλp and d_i = Bp + μT * Apλ * Bλ, where μT is a penalty parameter. The penalty parameter μT is heuristically chosen as μT = norm(App)/norm(Apλ*Aλp).

source
Gridap.Arrays.LocalSolveMapType
LocalSolveMap(; factorize! = lu!, pivot = NoPivot())

A map for solving local linear systems, relying on a factorization method.

Given a left-hand-side matrix mat and a set of N right-hand-side arrays lhs, returns an N-Tuple of arrays containing the solutions to the linear systems defined by

Each system is given by A*x_i = b_i, and the solution is computed as x_i = ldiv!(factorize!(A,pivot),b_i)

source
Gridap.Arrays.MapType

Abstract type representing a function (mapping) that provides a cache and an in-place evaluation for performance. This is the type to be used in the lazy_map function.

Derived types must implement the following method:

and optionally these ones:

The mapping interface can be tested with the test_map function.

Note that most of the functionality implemented in terms of this interface relies in duck typing. That is, it is not strictly needed to work with types that inherit from Map. This is specially useful in order to accommodate existing types into this framework without the need to implement a wrapper type that inherits from Map. For instance, a default implementation is available for Function objects. However, we recommend that new types inherit from Map.

source
Gridap.Arrays.MergeBlockMapType
struct MergeBlockMap{N,M} <: Map
    size::NTuple{N,Int}
    indices::Vector{Vector{Tuple{CartesianIndex{M},CartesianIndex{N}}}}
end

A MergeBlockMap create a single array of N-dimensional blocks from L=length(indices) input arrays of M-dimensional blocks.

For the l-th input array a_l, the vector of tuples in indices[l] contains the mapping between the indices of the blocks in a_l and the indices of the blocks in the output array, i.e a_out[P[2]] = a_l[P[1]] ∀ P ∈ indices[l], ∀ l.

source
Gridap.Arrays.OperationType
Operation(op)

Returns the map that results after applying an operation f over a set of map(s) args. That is Operation(f)(args)(x...) is formally defined as f(map(a->a(x...),args)...).

Example

using Gridap.Arrays

fa(x) = x.*x
fb(x) = sqrt.(x)

x = collect(0:5)

fab = Operation(fa)(fb)
c = evaluate(fab,x)

println(c)

# output
[0.0, 1.0, 2.0, 3.0, 4.0, 5.0]
source
Gridap.Arrays.OperationMapType
OperationMap(f,args)

Returns a mapping that represents the result of applying the function f to the arguments in the tuple args. That is, OperationMap(f,args)(x...) is formally defined as f(map(a->a(x...),args)...)

source
Gridap.Arrays.PosNegPartitionType

struct representing a binary partition of a range of indices

Using this allows one to do a number of important optimizations when working with PosNegReindex

source
Gridap.Arrays.SubVectorType
struct SubVector{T,A<:AbstractVector{T}} <: AbstractVector{T}
  vector::A
  pini::Int
  pend::Int
end

SubVector is deprecated, use view instead.

source
Gridap.Arrays.TableType
 struct Table{T,Vd<:AbstractVector{T},Vp<:AbstractVector} <: AbstractVector{Vector{T}}
    data::Vd
    ptrs::Vp
 end

Type representing a list of lists (i.e., a table) in compressed format.

source
Gridap.Arrays.TableMethod
Table(a::AbstractArray{<:AbstractArray})

Build a table from a vector of vectors. If the inputs are multidimensional arrays instead of vectors, they are flattened.

source
Gridap.Arrays.array_cacheMethod
array_cache(a::AbstractArray)

Returns a cache object to be used in the getindex! function. It defaults to

array_cache(a::T) where T = nothing

for types T such that uses_hash(T) == Val(false), and

function array_cache(a::T) where T
  hash = Dict{UInt,Any}()
  array_cache(hash,a)
end

for types T such that uses_hash(T) == Val(true), see the uses_hash function. In the later case, the type T should implement the following signature:

array_cache(hash::Dict,a::AbstractArray)

where we pass a dictionary (i.e., a hash table) in the first argument. This hash table can be used to test if the object a has already built a cache and re-use it as follows

id = objectid(a)
if haskey(hash,id)
  cache = hash[id] # Reuse cache
else
  cache = ... # Build a new cache depending on your needs
  hash[id] = cache # Register the cache in the hash table
end

This mechanism is needed, e.g., to re-use intermediate results in complex lazy operation trees. In multi-threading computations, a different hash table per thread has to be used in order to avoid race conditions.

source
Gridap.Arrays.block_identity_arrayMethod
block_identity_array(ptrs;T=Int)

Given a vector of pointers of length n+1, returns a vector of length ptrs[end]-1 where the entries are the index of the block to which each entry belongs.

Example


julia> block_identity_array([1,3,7])

6-element Vector{Int64}:
 1
 1
 2
 2
 2
 2
source
Gridap.Arrays.evaluate!Method
evaluate!(cache,f,x...)

Applies the mapping f at the arguments x... using the scratch data provided in the given cache object. The cache object is built with the return_cache function using arguments of the same type as in x. In general, the returned value y can share some part of its state with the cache object. If the result of two or more calls to this function need to be accessed simultaneously (e.g., in multi-threading), create and use several cache objects (e.g., one cache per thread).

source
Gridap.Arrays.evaluateMethod
evaluate(f,x...)

evaluates the mapping f at the arguments in x by creating a temporary cache internally. This functions is equivalent to

cache = return_cache(f,x...)
evaluate!(cache,f,x...)
source
Gridap.Arrays.find_inverse_index_mapFunction
find_inverse_index_map(a_to_b[, nb=maximum(a_to_b)])
find_inverse_index_map!(b_to_a, a_to_b)

Given a vector of indices a_to_b, returns the inverse index map b_to_a.

source
Gridap.Arrays.getindex!Method
getindex!(cache,a::AbstractArray,i...)

Returns the item of the array a associated with index i by (possibly) using the scratch data passed in the cache object.

It defaults to

getindex!(cache,a::AbstractArray,i...) = a[i...]

As for standard Julia arrays, the user needs to implement only one of the following signatures depending on the IndexStyle of the array.

getindex!(cache,a::AbstractArray,i::Integer)
getindex!(cache,a::AbstractArray{T,N},i::Vararg{Integer,N}) where {T,N}

Examples

Iterating over an array using the getindex! function

using Gridap.Arrays

a = collect(10:15)

cache = array_cache(a)
for i in eachindex(a)
  ai = getindex!(cache,a,i)
  println("$i -> $ai")
end

# output
1 -> 10
2 -> 11
3 -> 12
4 -> 13
5 -> 14
6 -> 15
source
Gridap.Arrays.lazy_mapMethod
lazy_map(f,a::AbstractArray...) -> AbstractArray

Applies the Map (or Function) f to the entries of the arrays in a (see the definition of Map).

The resulting array r is such that r[i] equals to evaluate(f,ai...) where ai is the tuple containing the i-th entry of the arrays in a (see function evaluate for more details). In other words, the resulting array is numerically equivalent to:

map( (x...)->evaluate(f,x...), a...)

Examples

Using a function as mapping

using Gridap.Arrays

a = collect(0:5)
b = collect(10:15)

c = lazy_map(+,a,b)

println(c)

# output
[10, 12, 14, 16, 18, 20]

Using a user-defined mapping

using Gridap.Arrays
import Gridap.Arrays: evaluate!

a = collect(0:5)
b = collect(10:15)

struct MySum <: Map end

evaluate!(cache,::MySum,x,y) = x + y

k = MySum()

c = lazy_map(k,a,b)

println(c)

# output
[10, 12, 14, 16, 18, 20]
source
Gridap.Arrays.local_identity_arrayMethod
local_identity_array(ptrs;T=Int)

Given a vector of pointers of length n+1, returns a vector of length ptrs[end]-1 where the entries are the local index of the entry within the block it belongs to.

Example


julia> local_identity_array([1,3,7])

6-element Vector{Int64}:
 1
 2
 1
 2
 3
 4
source
Gridap.Arrays.merge_entriesMethod
merge_entries(a_to_lb_to_b, c_to_la_to_a) -> c_to_lb_to_b

Merge the entries of a_to_lb_to_b, grouping them by c_to_la_to_a. Returns the merged table c_to_lb_to_b.

Accepts the following keyword arguments:

  • acc: Accumulator for the entries of a_to_lb_to_b. Default to a Set, ensuring that the resulting entries are unique.
  • post: Postprocessing function to apply to the accumulator before storing the resulting entries. Defaults to the identity, but can be used to perform local sorts or filters, for example.
source
Gridap.Arrays.return_cacheMethod
return_cache(f,x...)

Returns the cache needed to lazy_map mapping f with arguments of the same type as the objects in x. This function returns nothing by default, i.e., no cache.

source
Gridap.Arrays.return_typeMethod
return_type(f,x...)

Returns the type of the result of calling mapping f with arguments of the types of the objects x.

source
Gridap.Arrays.setsize!Method
setsize!(a, s)

Changes the size of the CachedArray a to the size described the the tuple s. After calling setsize!, the array can store uninitialized values.

source
Gridap.Arrays.setsize_op!Method
setsize_op!(op, c, args)

Sets the size of the CachedArray to accomodate the result of the operation given by op(args...).

source
Gridap.Arrays.test_arrayMethod
test_array(
  a::AbstractArray{T,N}, b::AbstractArray{S,N},cmp=(==)) where {T,S,N}

Checks if the entries in a and b are equal using the comparison function cmp. It also stresses the new methods added to the AbstractArray interface.

source
Gridap.Arrays.test_mapMethod
test_map(y,f,x...;cmp=(==))

Function used to test if the mapping f has been implemented correctly. f is a Map sub-type, x is a tuple in the domain of the mapping and y is the expected result. Function cmp is used to compare the computed result with the expected one. The checks are done with the @test macro.

source
Gridap.Arrays.testargsMethod
testargs(f,x...)

The default implementation of this function is testargs(f,x...) = x. One can overload it in order to use lazy_map with 0-length array and maps with non-trivial domains.

source
Gridap.Arrays.testitemMethod
testitem(a::AbstractArray{T}) -> Any

Returns an arbitrary instance of eltype(a). The default returned value is the first entry in the array if length(a)>0 and testvalue(eltype(a)) if length(a)==0 See the testvalue function.

Examples

using Gridap.Arrays

a = collect(3:10)
ai = testitem(a)

b = Int[]
bi = testitem(b)

(ai, bi)

# output
(3, 0)
source
Gridap.Arrays.testvalueFunction
testvalue(::Type{T}) where T

Returns an arbitrary instance of type T. It defaults to zero(T) for non-array types and to an empty array for array types. It can be overloaded for new types T if zero(T) does not makes sense. This function is used to compute testitem for 0-length arrays.

source
Gridap.Arrays.uses_hashMethod
uses_hash(::Type{<:AbstractArray})

This function is used to specify if the type T uses the hash-based mechanism to reuse caches. It should return either Val(true) or Val(false). It defaults to

uses_hash(::Type{<:AbstractArray}) = Val(false)

Once this function is defined for the type T it can also be called on instances of T.

source