Skip to content

Operators

All discrete operators are built using KernelAbstractions.jl and Cartesian indices, similar to WaterLily.jl. This allows for dimension- and backend-agnostic code. See this blog post for how to write kernels. IncompressibleNavierStokes previously relied on assembling sparse operators to perform the same operations. While being very efficient and also compatible with CUDA (CUSPARSE), storing these matrices in memory is expensive for large 3D problems.

IncompressibleNavierStokes.Offset Type
julia
struct Offset{D}

Cartesian index unit vector in D = 2 or D = 3 dimensions. Calling Offset(D)(i) returns a Cartesian index with 1 in the dimension i and zeros elsewhere.

See https://b-fg.github.io/research/2023-07-05-waterlily-on-gpu.html for writing kernel loops using Cartesian indices.

Fields

source
IncompressibleNavierStokes.apply! Method
julia
apply!(kernel, setup, args...; ndrange, offset)

Apply kernel to args with offset offset. By default, it is applied everywhere except for at the outermost boundary.

source
IncompressibleNavierStokes.applygravity! Method
julia
applygravity!(f, temp, setup, gdir, gravity) -> Any

Compute gravity term (in-place version). add the result to F.

source
IncompressibleNavierStokes.applygravity Method
julia
applygravity(temp, setup, gdir, gravity) -> Any

Compute gravity term (differentiable version).

source
IncompressibleNavierStokes.applypressure! Method
julia
applypressure!(u, p, setup) -> Any

Subtract pressure gradient (in-place version).

source
IncompressibleNavierStokes.avg Method
julia
avg(ϕ, Δ, I, i) -> Any

Average scalar field ϕ in the i-direction.

source
IncompressibleNavierStokes.convection! Method
julia
convection!(f, u, setup) -> Any

Compute convective term (in-place version). Add the result to F.

source
IncompressibleNavierStokes.convection Method
julia
convection(u, setup) -> Any

Compute convective term (differentiable version).

source
IncompressibleNavierStokes.convection_diffusion_temp! Method
julia
convection_diffusion_temp!(
    c,
    u,
    temp,
    setup,
    conductivity
) -> Any

Compute convection-diffusion term for the temperature equation. (in-place version). Add result to c.

source
IncompressibleNavierStokes.convection_diffusion_temp Method
julia
convection_diffusion_temp(
    u,
    temp,
    setup,
    conductivity
) -> Any

Compute convection-diffusion term for the temperature equation. (differentiable version).

source
IncompressibleNavierStokes.convectiondiffusion! Method
julia
convectiondiffusion!(f, u, setup, viscosity) -> Any

Compute diffusive term (in-place version). Add the result to F.

source
IncompressibleNavierStokes.diffusion! Method
julia
diffusion!(f, u, setup, viscosity) -> Any

Compute diffusive term (in-place version). Add the result to F.

source
IncompressibleNavierStokes.diffusion Method
julia
diffusion(u, setup, viscosity) -> Any

Compute diffusive term (differentiable version).

source
IncompressibleNavierStokes.dissipation! Method
julia
dissipation!(diss, u, setup, coeff) -> Any

Compute dissipation term for the temperature equation (in-place version). Add result to diss.

source
IncompressibleNavierStokes.dissipation Method
julia
dissipation(u, setup, coeff) -> Any

Compute dissipation term for the temperature equation (differentiable version).

source
IncompressibleNavierStokes.divergence! Method
julia
divergence!(div, u, setup) -> Any

Compute divergence of velocity field (in-place version).

source
IncompressibleNavierStokes.divergence Method
julia
divergence(u, setup) -> Any

Compute divergence of velocity field (differentiable version).

source
IncompressibleNavierStokes.get_scale_numbers Method
julia
get_scale_numbers(
    u,
    setup,
    viscosity
) -> NamedTuple{(:uavg, , , , :Reλ, :L, , :Re_int), <:NTuple{8, Any}}

Get the following dimensional scale numbers [3]:

  • Velocity uavg=uiui1/2

  • Dissipation rate ϵ=2νSijSij

  • Kolmolgorov length scale η=(ν3ϵ)1/4

  • Taylor length scale λ=(5νϵ)1/2uavg

  • Taylor-scale Reynolds number Reλ=λuavg3ν

  • Integral length scale L=3π2uavg20E(k)kdk

  • Large-eddy turnover time τ=Luavg

source
IncompressibleNavierStokes.gridsize Method
julia
gridsize(setup, I::CartesianIndex{2}) -> Any

Gridsize based on the length of the diagonal of the cell.

source
IncompressibleNavierStokes.gridsize_vol Method
julia
gridsize_vol(setup, I::CartesianIndex{2}) -> Any

Grid size based on the volume of the cell.

source
IncompressibleNavierStokes.interpolate_u_p! Method
julia
interpolate_u_p!(up, u, setup) -> Any

Interpolate velocity to pressure points (in-place version).

source
IncompressibleNavierStokes.interpolate_u_p Method
julia
interpolate_u_p(u, setup) -> Any

Interpolate velocity to pressure points (differentiable version).

source
IncompressibleNavierStokes.interpolate_ω_p! Method
julia
interpolate_ω_p!(ωp, ω, setup) -> Any

Interpolate vorticity to pressure points (in-place version).

source
IncompressibleNavierStokes.interpolate_ω_p Method
julia
interpolate_ω_p(ω, setup) -> Any

Interpolate vorticity to pressure points (differentiable version).

source
IncompressibleNavierStokes.kinetic_energy! Method
julia
kinetic_energy!(ke, u, setup; interpolate_first) -> Any

Compute kinetic energy field (in-place version).

source
IncompressibleNavierStokes.kinetic_energy Method
julia
kinetic_energy(u, setup; kwargs...) -> Any

Compute kinetic energy field k (in-place version). If interpolate_first is true, it is given by

kI=18α(uI+hαα+uIhαα)2.

Otherwise, it is given by

kI=14α((uI+hαα)2+(uIhαα)2),

as in [4].

Differentiable version.

source
IncompressibleNavierStokes.laplacian! Method
julia
laplacian!(L, p, setup) -> Any

Compute Laplacian of pressure field (in-place version).

source
IncompressibleNavierStokes.laplacian Method
julia
laplacian(p, setup) -> Any

Compute Laplacian of pressure field (differentiable version).

source
IncompressibleNavierStokes.left Method
julia
left(I::CartesianIndex{D}, i) -> Any
left(I::CartesianIndex{D}, i, n) -> Any

Left index n times away in direction i.

source
IncompressibleNavierStokes.pressuregradient! Method
julia
pressuregradient!(G, p, setup) -> Any

Compute pressure gradient (in-place version).

source
IncompressibleNavierStokes.pressuregradient Method
julia
pressuregradient(p, setup) -> Any

Compute pressure gradient (differentiable version).

source
IncompressibleNavierStokes.qcrit Method
julia
qcrit(u, setup) -> Any

Compute Q, the second invariant of the velocity gradient tensor.

source
IncompressibleNavierStokes.right Method
julia
right(I::CartesianIndex{D}, i) -> Any
right(I::CartesianIndex{D}, i, n) -> Any

Right index n times away in direction i.

source
IncompressibleNavierStokes.scalewithvolume! Method
julia
scalewithvolume!(p, setup) -> Any

Scale scalar field with volume sizes (in-place version).

source
IncompressibleNavierStokes.scalewithvolume Method
julia
scalewithvolume(p, setup) -> Any

Scale scalar field p with volume sizes (differentiable version).

source
IncompressibleNavierStokes.total_kinetic_energy Method
julia
total_kinetic_energy(u, setup; kwargs...) -> Any

Compute total kinetic energy. The velocity components are interpolated to the volume centers and squared.

source
IncompressibleNavierStokes.unit_cartesian_indices Method
julia
unit_cartesian_indices(D) -> Any

Get tuple of all unit vectors as Cartesian indices.

source
IncompressibleNavierStokes.vorticity! Method
julia
vorticity!(ω, u, setup) -> Any

Compute vorticity field (in-place version).

source
IncompressibleNavierStokes.vorticity Method
julia
vorticity(u, setup) -> Any

Compute vorticity field (differentiable version).

source
IncompressibleNavierStokes.δ Method
julia
δ(setup, p, i, I) -> Any

Differentiate scalar p in direction ei.

source
IncompressibleNavierStokes.δ Method
julia
δ(setup, u, i, j, I) -> Any

Differentiate vector ui in direction ej. Make sure that i and j are known at compile-time to remove the if-statement.

source