WaterLily
WaterLily.AbstractBody
— TypeAbstractBody
Abstract body geometry data structure. This solver will call
`measure(body::AbstractBody,x::Vector,t::Real)`
to query the body geometry for these properties at location x
and time t
:
`d :: Real`, Signed distance to surface
`n̂ :: Vector`, Outward facing unit normal
`V :: Vector`, Body velocity
`H,K :: Real`, Mean and Gaussian curvature
WaterLily.AbstractPoisson
— TypePoisson{N,M}
Composite type for conservative variable coefficient Poisson equations:
∮ds β ∂x/∂n = σ
The resulting linear system is
Ax = [L+D+L']x = b
where A is symmetric, block-tridiagonal and extremely sparse. Implemented on a structured grid of dimension N, then L has dimension M=N+1 and size(L,M)=N. Moreover, D[I]=-∑ᵢ(L[I,i]+L'[I,i]). This means matrix storage, multiplication, ect can be easily implemented and optimized without external libraries.
To help iteratively solve the system above, the Poisson structure holds helper arrays for inv(D), the error ϵ, and residual r=b-Ax. An iterative solution method then estimates the error ϵ=̃A⁻¹r and increments x+=ϵ, r-=Aϵ.
WaterLily.AutoBody
— TypeAutoBody(sdf,map=(x,t)->x) <: AbstractBody
- sdf(x::AbstractVector,t::Real)::Real: signed distance function
- map(x::AbstractVector,t::Real)::AbstractVector: coordinate mapping function
Define a geometry by its sdf
and optional coordinate map
. All other properties are determined using Automatic Differentiation. Note: the map
is composed automatically if provided, ie sdf(x,t) = sdf(map(x,t),t)
.
WaterLily.NoBody
— TypeNoBody
Use for a simulation without a body
WaterLily.Simulation
— TypeSimulation(dims::Tuple, u_BC::Vector, L::Number;
U=norm2(u_BC), Δt=0.25, ν=0., ϵ = 1,
uλ::Function=(i,x)->u_BC[i],
body::AbstractBody=NoBody())
Constructor for a WaterLily.jl simulation: - dims
: Simulation domain dimensions. - u_BC
: Simulation domain velocity boundary conditions, u_BC[i]=uᵢ, i=1,2...
. - L
: Simulation length scale. - U
: Simulation velocity scale. - ϵ
: BDIM kernel width. - Δt
: Initial time step. - ν
: Scaled viscosity (Re=UL/ν
) - uλ
: Function to generate the initial velocity field. - body
: Embedded geometry
See files in examples
folder for examples.
WaterLily.GS!
— MethodGS!(p::Poisson;it=0)
Gauss-Sidel smoother. When it=0, the function serves as a Jacobi preconditioner.
WaterLily.SOR!
— MethodSOR!(p::Poisson;ω=1.5)
Successive Over Relaxation preconditioner. The routine uses backsubstitution to compute ϵ=̃A⁻¹r, where ̃A=[D/ω+L], and then increments x,r.
WaterLily.apply!
— Methodapply!(f, c)
Apply a vector function f(i,x)
to the faces of a uniform staggered array c
.
WaterLily.curl
— Methodcurl(i,I,u)
Compute component i
of ∇×u at the edge of cell I
. For example curl(3,CartesianIndex(2,2,2),u)
will compute ω₃(x=1.5,y=1.5,z=2)
as this edge produces the highest accuracy for this mix of cross derivatives on a staggered grid.
WaterLily.curvature
— Methodcurvature(A::AbstractMatrix)
Return H,K
the mean and Gaussian curvature from A=hessian(sdf)
. K=tr(minor(A)) in 3D and K=0 in 2D.
WaterLily.ke
— Methodke(I::CartesianIndex,u,U=0)
Compute ½|u-U|² at center of cell I
where U
can be used to subtract a background flow.
WaterLily.measure!
— Functionmeasure!(sim::Simulation,t=time(sim))
Measure a dynamic body
to update the flow
and pois
coefficients.
WaterLily.measure!
— Methodmeasure(a::Flow,body::AbstractBody;t=0,ϵ=1)
Measure the body
properties on flow
using a kernel size ϵ
. Weymouth & Yue, JCP, 2011
WaterLily.measure
— Methodmeasure(body::AutoBody,x::Vector,t::Real)
ForwardDiff is used to determine the geometric properties from the sdf
. Note: The velocity is determined soley from the optional map
function.
WaterLily.sim_step!
— Methodsim_step!(sim::Simulation,t_end;verbose=false)
Integrate the simulation sim
up to dimensionless time t_end
. If verbose=true
the time tU/L
and adaptive time step Δt
are printed every time step.
WaterLily.sim_time
— Methodsim_time(sim::Simulation)
Return the current dimensionless time of the simulation tU/L
where t=sum(Δt)
, and U
,L
are the simulation velocity and length scales.
WaterLily.slice
— Methodslice(N,s,dims) -> R
Return CartesianIndices
slicing through an array of size N
.
WaterLily.λ₂
— Methodλ₂(I::CartesianIndex{3},u)
λ₂ is a deformation tensor metric to identify vortex cores. See https://en.wikipedia.org/wiki/Lambda2_method and Jeong, J., & Hussain, F. doi:10.1017/S0022112095000462
WaterLily.ω
— Methodω(I::CartesianIndex{3},u)
Compute 3-vector ω=∇×u at the center of cell I
.
WaterLily.ω_mag
— Methodω_mag(I::CartesianIndex{3},u)
Compute |ω| at the center of cell I
.
WaterLily.ω_θ
— Methodω_θ(I::CartesianIndex{3},z,center,u)
Compute ω⋅θ at the center of cell I
where θ is the azimuth direction around vector z
passing through center
.
WaterLily.∂
— Method∂(i,j,I,u)
Compute ∂uᵢ/∂xⱼ at center of cell I
. Cross terms are computed less accurately than inline terms because of the staggered grid.
WaterLily.AbstractBody
WaterLily.AbstractPoisson
WaterLily.AutoBody
WaterLily.NoBody
WaterLily.Simulation
WaterLily.GS!
WaterLily.SOR!
WaterLily.apply!
WaterLily.curl
WaterLily.curvature
WaterLily.ke
WaterLily.measure
WaterLily.measure!
WaterLily.measure!
WaterLily.sim_step!
WaterLily.sim_time
WaterLily.slice
WaterLily.λ₂
WaterLily.ω
WaterLily.ω_mag
WaterLily.ω_θ
WaterLily.∂