Skip to content
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

fix complex adaptive norms on GPU #347

Merged
merged 2 commits into from
Oct 7, 2019
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
10 changes: 5 additions & 5 deletions src/common_defaults.jl
Original file line number Diff line number Diff line change
@@ -1,12 +1,12 @@
@inline UNITLESS_ABS2(x) = Real(abs2(x)/(typeof(x)(one(x))*typeof(x)(one(x))))
@inline UNITLESS_ABS2(x) = real(abs2(x)/(typeof(x)(one(x))*typeof(x)(one(x))))

@inline ODE_DEFAULT_NORM(u::Union{AbstractFloat,Complex},t) = abs(u)
@inline ODE_DEFAULT_NORM(u::Array{T},t) where T<:Union{AbstractFloat,Complex} =
@fastmath sqrt(sum(abs2,u) / length(u))
@fastmath sqrt(real(sum(abs2,u)) / length(u))
@inline ODE_DEFAULT_NORM(u::RecursiveArrayTools.AbstractVectorOfArray,t) =
sum(sqrt(sum(UNITLESS_ABS2,_u) / length(_u)) for _u in u.u)
@inline ODE_DEFAULT_NORM(u::AbstractArray,t) = sqrt(sum(UNITLESS_ABS2,u) / length(u))
@inline ODE_DEFAULT_NORM(u::AbstractArray{T,N},t) where {T<:AbstractArray,N} = sqrt(sum(x->ODE_DEFAULT_NORM(x[1],x[2]),zip(u,t)) / length(u))
sum(sqrt(real(sum(UNITLESS_ABS2,_u)) / length(_u)) for _u in u.u)
@inline ODE_DEFAULT_NORM(u::AbstractArray,t) = sqrt(real(sum(UNITLESS_ABS2,u)) / length(u))
@inline ODE_DEFAULT_NORM(u::AbstractArray{T,N},t) where {T<:AbstractArray,N} = sqrt(real(sum(x->ODE_DEFAULT_NORM(x[1],x[2]),zip(u,t))) / length(u))
@inline ODE_DEFAULT_NORM(u,t) = norm(u)

@inline ODE_DEFAULT_ISOUTOFDOMAIN(u,p,t) = false
Expand Down
6 changes: 6 additions & 0 deletions test/gpu/simple_gpu.jl
Original file line number Diff line number Diff line change
Expand Up @@ -55,3 +55,9 @@ u0 = [1.0;0.0;0.0]
tspan = (0.0,100.0)
prob_num = ODEProblem(ff2,u0,tspan)
sol = solve(prob_num,Rosenbrock23(linsolve=LinSolveGPUFactorize()))

# Complex Numbers Adaptivity DifferentialEquations.jl#460
f_complex(u,nothing,t) = 1/2 .*u
u0 = cu(rand(32,32).+ 1im*rand(32,32));
prob = ODEProblem(f_complex,u0,(0.0,1.0))
sol = solve(prob,Tsit5())