Skip to content

Make ArrayPartitions of GPU arrays work in DifferentialEquations solvers #416

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 2 commits into from
Jan 29, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
7 changes: 7 additions & 0 deletions src/array_partition.jl
Original file line number Diff line number Diff line change
Expand Up @@ -209,6 +209,13 @@ function Base.copyto!(A::ArrayPartition, src::ArrayPartition)
A
end

function recursivefill!(b::ArrayPartition, a::T2) where {T2 <: Union{Number, Bool}}
unrolled_foreach!(b.x) do x
fill!(x, a)
end
end


## indexing

# Interface for the linear indexing. This is just a view of the underlying nested structure
Expand Down
5 changes: 5 additions & 0 deletions src/utils.jl
Original file line number Diff line number Diff line change
@@ -1,3 +1,7 @@
unrolled_foreach!(f, t::Tuple) = (f(t[1]); unrolled_foreach!(f, Base.tail(t)))
unrolled_foreach!(f, ::Tuple{}) = nothing


"""
```julia
recursivecopy(a::Union{AbstractArray{T, N}, AbstractVectorOfArray{T, N}})
Expand Down Expand Up @@ -127,6 +131,7 @@ function recursivefill!(bs::AbstractVectorOfArray{T, N},
end
end


for type in [AbstractArray, AbstractVectorOfArray]
@eval function recursivefill!(b::$type{T, N}, a::T2) where {T <: Enum, T2 <: Enum, N}
fill!(b, a)
Expand Down
16 changes: 16 additions & 0 deletions test/gpu/arraypartition_gpu.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,16 @@
using RecursiveArrayTools, CUDA, Test
CUDA.allowscalar(false)


# Test indexing with colon
a = (CUDA.zeros(5), CUDA.zeros(5))
pA = ArrayPartition(a)
pA[:, :]

# Indexing with boolean masks does not work yet
mask = pA .> 0
# pA[mask]

# Test recursive filling is done using GPU kernels and not scalar indexing
RecursiveArrayTools.recursivefill!(pA, true)
@test all(pA .== true)
1 change: 1 addition & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -54,5 +54,6 @@ end
if GROUP == "GPU"
activate_gpu_env()
@time @safetestset "VectorOfArray GPU" include("gpu/vectorofarray_gpu.jl")
@time @safetestset "ArrayPartition GPU" include("gpu/arraypartition_gpu.jl")
end
end
Loading