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

Use fast measure(body) in measure!(flow) #155

Merged
merged 1 commit into from
Jul 31, 2024
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
19 changes: 9 additions & 10 deletions src/AutoBody.jl
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ Base.:-(x::AutoBody, y::AutoBody) = x ∩ -y
"""
d = sdf(body::AutoBody,x,t) = body.sdf(x,t)
"""
sdf(body::AutoBody,x,t) = body.sdf(x,t)
sdf(body::AutoBody,x,t;kwargs...) = body.sdf(x,t)

"""
Bodies(bodies, ops::AbstractVector)
Expand Down Expand Up @@ -96,27 +96,26 @@ end

Computes distance for `Bodies` type.
"""
sdf(a::Bodies,x,t) = sdf_map_d(a.bodies,a.ops,x,t)[end]
sdf(a::Bodies,x,t;kwargs...) = sdf_map_d(a.bodies,a.ops,x,t)[end]

using ForwardDiff
"""
d,n,V = measure(body::AutoBody,x,t;fast=false)
d,n,V = measure(body::Bodies,x,t;fast=false)
d,n,V = measure(body::AutoBody||Bodies,x,t;fastd²=Inf)

Determine the implicit geometric properties from the `sdf` and `map`.
The gradient of `d=sdf(map(x,t))` is used to improve `d` for pseudo-sdfs.
The velocity is determined _solely_ from the optional `map` function.
Using `fast=true` skips the `n,V` calculation when |d|>1.
Skips the `n,V` calculation when `d²>fastd²`.
"""
measure(body::AutoBody,x,t;fast=false) = measure(body.sdf,body.map,x,t,fast)
function measure(a::Bodies,x,t;fast=false)
measure(body::AutoBody,x,t;kwargs...) = measure(body.sdf,body.map,x,t;kwargs...)
function measure(a::Bodies,x,t;kwargs...)
sdf, map, _ = sdf_map_d(a.bodies,a.ops,x,t)
measure(sdf,map,x,t,fast)
measure(sdf,map,x,t;kwargs...)
end
function measure(sdf,map,x,t,fast)
function measure(sdf,map,x,t;fastd²=Inf)
# eval d=f(x,t), and n̂ = ∇f
d = sdf(x,t)
fast && abs(d)>1 && return (d,zero(x),zero(x)) # skip n,V
d^2>fastd² && return (d,zero(x),zero(x)) # skip n,V
n = ForwardDiff.gradient(x->sdf(x,t), x)
any(isnan.(n)) && return (d,zero(x),zero(x))

Expand Down
13 changes: 7 additions & 6 deletions src/Body.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,10 +8,11 @@ Immersed body Abstract Type. Any `AbstractBody` subtype must implement

and

d,n,V = measure(body::AbstractBody, x, t=0)
d,n,V = measure(body::AbstractBody, x, t=0, fastd²=Inf)

where `d` is the signed distance from `x` to the body at time `t`,
and `n` & `V` are the normal and velocity vectors implied at `x`.
A fast-approximate method can return `≈d,zero(x),zero(x)` if `d^2>fastd²`.
"""
abstract type AbstractBody end
"""
Expand All @@ -28,12 +29,12 @@ at time `t` using an immersion kernel of size `ϵ`.
See Maertens & Weymouth, doi:[10.1016/j.cma.2014.09.007](https://doi.org/10.1016/j.cma.2014.09.007).
"""
function measure!(a::Flow{N,T},body::AbstractBody;t=zero(T),ϵ=1) where {N,T}
a.V .= zero(T); a.μ₀ .= one(T); a.μ₁ .= zero(T)
a.V .= zero(T); a.μ₀ .= one(T); a.μ₁ .= zero(T); d²=(2+ϵ)^2
@fastmath @inline function fill!(μ₀,μ₁,V,d,I)
d[I] = sdf(body,loc(0,I,T),t)
if abs(d[I])<2+ϵ
d[I] = sdf(body,loc(0,I,T),t,fastd²=d²)
Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This is the important line for ParametricBodies.

if d[I]^2<d²
for i ∈ 1:N
dᵢ,nᵢ,Vᵢ = measure(body,loc(i,I,T),t)
dᵢ,nᵢ,Vᵢ = measure(body,loc(i,I,T),t,fastd²=d²)
Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This line won't do too much, since the if statement above is already catching most of them. I tried to reduce to d=epsilon, but this broke the code.

V[I,i] = Vᵢ[i]
μ₀[I,i] = WaterLily.μ₀(dᵢ,ϵ)
for j ∈ 1:N
Expand Down Expand Up @@ -64,7 +65,7 @@ end

Uses `sdf(body,x,t)` to fill `a`.
"""
measure_sdf!(a::AbstractArray,body::AbstractBody,t=0) = @inside a[I] = sdf(body,loc(0,I,eltype(a)),t)
measure_sdf!(a::AbstractArray,body::AbstractBody,t=0;kwargs...) = @inside a[I] = sdf(body,loc(0,I,eltype(a)),t;kwargs...)

"""
NoBody
Expand Down
2 changes: 1 addition & 1 deletion src/Metrics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -82,7 +82,7 @@ end
BDIM-masked surface normal.
"""
@inline function nds(body,x,t)
d,n,_ = measure(body,x,t,fast=true)
d,n,_ = measure(body,x,t,fastd²=1)
n*WaterLily.kern(clamp(d,-1,1))
end

Expand Down
4 changes: 4 additions & 0 deletions test/maintests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -223,6 +223,10 @@ end
I = CartesianIndex(2,3)
@test GPUArrays.@allowscalar p[I]≈body1.sdf(loc(0,I,eltype(p)),0.0)
end

# check fast version
@test all(measure(body1,[3.,4.],0.,fastd²=9) .≈ measure(body1,[3.,4.],0.))
@test all(measure(body1,[3.,4.],0.,fastd²=8) .≈ (sdf(body1,[3.,4.],0.,fastd²=9),zeros(2),zeros(2)))
end

function TGVsim(mem;perdir=(1,2),Re=1e8,T=typeof(Re))
Expand Down
Loading