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

Non-uniform boundary conditions #165

Draft
wants to merge 5 commits into
base: master
Choose a base branch
from
Draft
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
38 changes: 16 additions & 22 deletions src/Flow.jl
Original file line number Diff line number Diff line change
Expand Up @@ -65,19 +65,13 @@ using EllipsisNotation

Add a uniform acceleration `gᵢ+dUᵢ/dt` at time `t=sum(dt)` to field `r`.
"""
accelerate!(r,dt,g::Function,::Tuple,t=sum(dt)) = for i ∈ 1:last(size(r))
r[..,i] .+= g(i,t)
accelerate!(r,U::Function,t) = for i ∈ 1:last(size(r))
@loop r[I,i] += U(i,loc(i,I,eltype(r)),t) over I ∈ CartesianIndices(Base.front(size(r)))
end
accelerate!(r,dt,g::Nothing,U::Function) = accelerate!(r,dt,(i,t)->ForwardDiff.derivative(τ->U(i,τ),t),())
accelerate!(r,dt,g::Function,U::Function) = accelerate!(r,dt,(i,t)->g(i,t)+ForwardDiff.derivative(τ->U(i,τ),t),())
accelerate!(r,dt,::Nothing,::Tuple) = nothing
"""
BCTuple(U,dt,N)

Return BC tuple `U(i∈1:N, t=sum(dt))`.
"""
BCTuple(f::Function,dt,N,t=sum(dt))=ntuple(i->f(i,t),N)
BCTuple(f::Tuple,dt,N)=f
accelerate!(r,t,g::Nothing,U::Function) = accelerate!(r,(i,x,t)->ForwardDiff.derivative(τ->U(i,x,τ),t),t)
accelerate!(r,t,g::Function,U::Function) = accelerate!(r,(i,x,t)->g(i,t)+ForwardDiff.derivative(τ->U(i,x,τ),t),t)
accelerate!(r,t,g::Function,::Tuple) = accelerate!(r,(i,x,t)->g(i,t),t)
accelerate!(r,t,::Nothing,::Tuple) = nothing

"""
Flow{D::Int, T::Float, Sf<:AbstractArray{T,D}, Vf<:AbstractArray{T,D+1}, Tf<:AbstractArray{T,D+2}}
Expand Down Expand Up @@ -112,7 +106,7 @@ struct Flow{D, T, Sf<:AbstractArray{T}, Vf<:AbstractArray{T}, Tf<:AbstractArray{
Ng = N .+ 2
Nd = (Ng..., D)
u = Array{T}(undef, Nd...) |> f; apply!(uλ, u);
BC!(u,BCTuple(U,0.,D),exitBC,perdir); exitBC!(u,u,BCTuple(U,0.,D),0.)
BC!(u,U,exitBC,perdir); exitBC!(u,u,0.)
u⁰ = copy(u);
fv, p, σ = zeros(T, Nd) |> f, zeros(T, Ng) |> f, zeros(T, Ng) |> f
V, μ₀, μ₁ = zeros(T, Nd) |> f, ones(T, Nd) |> f, zeros(T, Ng..., D, D) |> f
Expand Down Expand Up @@ -153,18 +147,18 @@ and the `AbstractPoisson` pressure solver to project the velocity onto an incomp
@fastmath function mom_step!(a::Flow{N},b::AbstractPoisson) where N
a.u⁰ .= a.u; scale_u!(a,0)
# predictor u → u'
U = BCTuple(a.U,@view(a.Δt[1:end-1]),N)
t = sum(@view(a.Δt[1:end-1]))
conv_diff!(a.f,a.u⁰,a.σ,ν=a.ν,perdir=a.perdir)
accelerate!(a.f,@view(a.Δt[1:end-1]),a.g,a.U)
BDIM!(a); BC!(a.u,U,a.exitBC,a.perdir)
a.exitBC && exitBC!(a.u,a.u⁰,U,a.Δt[end]) # convective exit
project!(a,b); BC!(a.u,U,a.exitBC,a.perdir)
accelerate!(a.f,t,a.g,a.U)
BDIM!(a); BC!(a.u,a.U,a.exitBC,a.perdir,t)
a.exitBC && exitBC!(a.u,a.u⁰,a.Δt[end]) # convective exit
project!(a,b); BC!(a.u,a.U,a.exitBC,a.perdir,t)
# corrector u → u¹
U = BCTuple(a.U,a.Δt,N)
t = sum(a.Δt)
conv_diff!(a.f,a.u,a.σ,ν=a.ν,perdir=a.perdir)
accelerate!(a.f,a.Δt,a.g,a.U)
BDIM!(a); scale_u!(a,0.5); BC!(a.u,U,a.exitBC,a.perdir)
project!(a,b,0.5); BC!(a.u,U,a.exitBC,a.perdir)
accelerate!(a.f,t,a.g,a.U)
BDIM!(a); scale_u!(a,0.5); BC!(a.u,a.U,a.exitBC,a.perdir,t)
project!(a,b,0.5); BC!(a.u,a.U,a.exitBC,a.perdir,t)
push!(a.Δt,CFL(a))
end
scale_u!(a,scale) = @loop a.u[Ii] *= scale over Ii ∈ inside_u(size(a.p))
Expand Down
6 changes: 3 additions & 3 deletions src/WaterLily.jl
Original file line number Diff line number Diff line change
Expand Up @@ -39,7 +39,7 @@ Constructor for a WaterLily.jl simulation:

- `dims`: Simulation domain dimensions.
- `u_BC`: Simulation domain velocity boundary conditions, either a
tuple `u_BC[i]=uᵢ, i=eachindex(dims)`, or a time-varying function `f(i,t)`
tuple `u_BC[i]=uᵢ, i=eachindex(dims)`, or a time and space-varying function `f(i,x,t)`
- `L`: Simulation length scale.
- `U`: Simulation velocity scale.
- `Δt`: Initial time step.
Expand Down Expand Up @@ -68,8 +68,8 @@ mutable struct Simulation
T=Float32, mem=Array) where N
@assert !(isa(u_BC,Function) && isa(uλ,Function)) "`u_BC` and `uλ` cannot be both specified as Function"
@assert !(isnothing(U) && isa(u_BC,Function)) "`U` must be specified if `u_BC` is a Function"
isa(u_BC,Function) && @assert all(typeof.(ntuple(i->u_BC(i,zero(T)),N)).==T) "`u_BC` is not type stable"
uλ = isnothing(uλ) ? ifelse(isa(u_BC,Function),(i,x)->u_BC(i,0.),(i,x)->u_BC[i]) : uλ
isa(u_BC,Function) && @assert all(typeof.(ntuple(i->u_BC(i,zeros(SVector{N}),zero(T)),N)).==T) "`u_BC` is not type stable"
uλ = isnothing(uλ) ? ifelse(isa(u_BC,Function),(i,x)->u_BC(i,x,0.),(i,x)->u_BC[i]) : uλ
U = isnothing(U) ? √sum(abs2,u_BC) : U # default if not specified
flow = Flow(dims,u_BC;uλ,Δt,ν,g,T,f=mem,perdir,exitBC)
measure!(flow,body;ϵ)
Expand Down
16 changes: 10 additions & 6 deletions src/util.jl
Original file line number Diff line number Diff line change
Expand Up @@ -166,7 +166,7 @@ condition `a[I,i]=A[i]` is applied to the vector component _normal_ to the domai
boundary. For example `aₓ(x)=Aₓ ∀ x ∈ minmax(X)`. A zero Neumann condition
is applied to the tangential components.
"""
function BC!(a,A,saveexit=false,perdir=())
function BC!(a,A,saveexit=false,perdir=(),t=0)
N,n = size_u(a)
for i ∈ 1:n, j ∈ 1:n
if j in perdir
Expand All @@ -175,26 +175,30 @@ function BC!(a,A,saveexit=false,perdir=())
else
if i==j # Normal direction, Dirichlet
for s ∈ (1,2)
@loop a[I,i] = A[i] over I ∈ slice(N,s,j)
@loop a[I,i] = uBC(i,A,I,t) over I ∈ slice(N,s,j)
end
(!saveexit || i>1) && (@loop a[I,i] = A[i] over I ∈ slice(N,N[j],j)) # overwrite exit
(!saveexit || i>1) && (@loop a[I,i] = uBC(i,A,I,t) over I ∈ slice(N,N[j],j)) # overwrite exit
else # Tangential directions, Neumann
@loop a[I,i] = a[I+δ(j,I),i] over I ∈ slice(N,1,j)
@loop a[I,i] = a[I-δ(j,I),i] over I ∈ slice(N,N[j],j)
end
end
end
end
uBC(i,A,I,t) = A[i]
uBC(i,A::Function,I,t) = A(i,loc(i,I),t)

"""
exitBC!(u,u⁰,U,Δt)

Apply a 1D convection scheme to fill the ghost cell on the exit of the domain.
"""
function exitBC!(u,u⁰,U,Δt)
function exitBC!(u,u⁰,Δt)
N,_ = size_u(u)
exitR = slice(N.-1,N[1],1,2) # exit slice excluding ghosts
@loop u[I,1] = u⁰[I,1]-U[1]*Δt*(u⁰[I,1]-u⁰[I-δ(1,I),1]) over I ∈ exitR
∮u = sum(u[exitR,1])/length(exitR)-U[1] # mass flux imbalance
U = sum(@view(u[slice(N.-1,2,1,2),1]))/length(exitR) # inflow mass flux
@loop u[I,1] = u⁰[I,1]-U*Δt*(u⁰[I,1]-u⁰[I-δ(1,I),1]) over I ∈ exitR
∮u = sum(@view(u[exitR,1]))/length(exitR)-U # mass flux imbalance
@loop u[I,1] -= ∮u over I ∈ exitR # correct flux
end
"""
Expand Down
17 changes: 9 additions & 8 deletions test/maintests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -42,9 +42,13 @@ using ReadVTK, WriteVTK
BC!(u,U,true) # save exit values
@test GPUArrays.@allowscalar all(u[end, :, 1] .== 3)

WaterLily.exitBC!(u,u,U,0) # conservative exit check
WaterLily.exitBC!(u,u,0) # conservative exit check
@test GPUArrays.@allowscalar all(u[end,2:end-1, 1] .== U[1])

A = [1.,2.,3]; Af(i,x,t) = i
@test all([WaterLily.uBC(i,A,CartesianIndex(1,1,1),0) for i in 1:3] .== A)
@test all([WaterLily.uBC(i,Af,CartesianIndex(1,1,1),0) for i in 1:3] .== A)

BC!(u,U,true,(2,)) # periodic in y and save exit values
@test GPUArrays.@allowscalar all(u[:, 1:2, 1] .== u[:, end-1:end, 1]) && all(u[:, 1:2, 1] .== u[:,end-1:end,1])
WaterLily.perBC!(σ,(1,2)) # periodic in two directions
Expand Down Expand Up @@ -154,19 +158,16 @@ end
Ip = WaterLily.CIj(1,I,length(f)-2); # make periodic
@test ϕuP(1,Ip,I,f,1)==λ(f[Ip],f[I-δ(1,I)],f[I])

@test all(WaterLily.BCTuple((1,2,3),[0],3).==WaterLily.BCTuple((i,t)->i,0,3))
@test all(WaterLily.BCTuple((i,t)->t,[1.234],3).==ntuple(i->1.234,3))

# check applying acceleration
for f ∈ arrays
N = 4; a = zeros(N,N,2) |> f
WaterLily.accelerate!(a,[1],nothing,())
WaterLily.accelerate!(a,1,nothing,())
@test all(a .== 0)
WaterLily.accelerate!(a,[1],(i,t) -> i==1 ? t : 2*t,())
WaterLily.accelerate!(a,1,(i,t) -> i==1 ? t : 2*t,())
@test all(a[:,:,1] .== 1) && all(a[:,:,2] .== 2)
WaterLily.accelerate!(a,[1],nothing,(i,t) -> i==1 ? -t : -2*t)
WaterLily.accelerate!(a,1,nothing,(i,x,t) -> i==1 ? -t : -2*t)
@test all(a[:,:,1] .== 0) && all(a[:,:,2] .== 0)
WaterLily.accelerate!(a,[1],(i,t) -> i==1 ? t : 2*t,(i,t) -> i==1 ? -t : -2*t)
WaterLily.accelerate!(a,1,(i,t) -> i==1 ? t : 2*t,(i,x,t) -> i==1 ? -t : -2*t)
@test all(a[:,:,1] .== 0) && all(a[:,:,2] .== 0)
end
# Impulsive flow in a box
Expand Down