Skip to content

Fix typos in cell_diffusion_timescale methods. #557

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

Merged
merged 3 commits into from
Dec 14, 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
13 changes: 9 additions & 4 deletions src/TurbulenceClosures/turbulence_closure_diagnostics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,7 @@ function cell_diffusion_timescale(closure::ConstantIsotropicDiffusivity, diffusi
return min(Δ^2 / closure.ν, Δ^2 / max_κ)
end

function cell_diffusion_timescale(closure::ConstantAnisotropicDiffusivity, diffusivies, grid)
function cell_diffusion_timescale(closure::ConstantAnisotropicDiffusivity, diffusivities, grid)
Δh = min_Δxy(grid)
Δz = min_Δz(grid)
max_κh = maximum(closure.κh)
Expand All @@ -21,7 +21,7 @@ function cell_diffusion_timescale(closure::ConstantAnisotropicDiffusivity, diffu
Δz^2 / max_κv, Δh^2 / max_κh)
end

function cell_diffusion_timescale(closure::AnisotropicBiharmonicDiffusivity, diffusivies, grid)
function cell_diffusion_timescale(closure::AnisotropicBiharmonicDiffusivity, diffusivities, grid)
Δh = min_Δxy(grid)
Δz = min_Δz(grid)
max_κh = maximum(closure.κh)
Expand All @@ -30,18 +30,23 @@ function cell_diffusion_timescale(closure::AnisotropicBiharmonicDiffusivity, dif
Δz^4 / max_κv, Δh^4 / max_κh)
end

function cell_diffusion_timescale(closure::Union{AbstractSmagorinsky, AbstractLeith}, diffusivies, grid)
function cell_diffusion_timescale(closure::AbstractSmagorinsky, diffusivities, grid)
Δ = min_Δxyz(grid)
min_Pr = minimum(closure.Pr)
max_κ = maximum(closure.κ)
max_νκ = maximum(diffusivities.νₑ.data.parent) * max(1, 1/min_Pr)
return min(Δ^2 / max_νκ, Δ^2 / max_κ)
end

function cell_diffusion_timescale(closure::AbstractAnisotropicMinimumDissipation, diffusivies, grid)
function cell_diffusion_timescale(closure::AbstractAnisotropicMinimumDissipation, diffusivities, grid)
Δ = min_Δxyz(grid)
max_ν = maximum(diffusivities.νₑ.data.parent)
max_κ = max(Tuple(maximum(κₑ.data.parent) for κₑ in diffusivities.κₑ)...)
return min(Δ^2 / max_ν, Δ^2 / max_κ)
end

function cell_diffusion_timescale(closure::AbstractLeith, diffusivities, grid)
Δ = min_Δxyz(grid)
max_ν = maximum(diffusivities.νₑ.data.parent)
return Δ^2 / max_ν
end
Original file line number Diff line number Diff line change
Expand Up @@ -38,12 +38,12 @@ and `C` is a model constant.
Keyword arguments
=================
- `C` : Model constant
- `C_Redi` : Coefficient for down-gradient tracer diffusivity for each tracer. i
Either a constant applied to every
tracer, or a `NamedTuple` with fields for each tracer individually.
- `C_GM` : Coefficient for down-gradient tracer diffusivity for each tracer. i
Either a constant applied to every
tracer, or a `NamedTuple` with fields for each tracer individually.
- `C_Redi` : Coefficient for down-gradient tracer diffusivity for each tracer.
Either a constant applied to every tracer, or a `NamedTuple` with fields
for each tracer individually.
- `C_GM` : Coefficient for down-gradient tracer diffusivity for each tracer.
Either a constant applied to every tracer, or a `NamedTuple` with fields
for each tracer individually.

References
==========
Expand Down
18 changes: 9 additions & 9 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -50,15 +50,15 @@ archs = (CPU(),)
@hascuda archs = (CPU(), GPU())

closures = (
:ConstantIsotropicDiffusivity,
:ConstantAnisotropicDiffusivity,
:AnisotropicBiharmonicDiffusivity,
:TwoDimensionalLeith,
:SmagorinskyLilly,
:BlasiusSmagorinsky,
:RozemaAnisotropicMinimumDissipation,
:VerstappenAnisotropicMinimumDissipation
)
:ConstantIsotropicDiffusivity,
:ConstantAnisotropicDiffusivity,
:AnisotropicBiharmonicDiffusivity,
:TwoDimensionalLeith,
:SmagorinskyLilly,
:BlasiusSmagorinsky,
:RozemaAnisotropicMinimumDissipation,
:VerstappenAnisotropicMinimumDissipation
)

@testset "Oceananigans" begin
include("test_grids.jl")
Expand Down
76 changes: 47 additions & 29 deletions test/test_turbulence_closures.jl
Original file line number Diff line number Diff line change
@@ -1,3 +1,5 @@
using Oceananigans.Diagnostics

for closure in closures
@eval begin
using Oceananigans.TurbulenceClosures: $closure
Expand All @@ -11,23 +13,6 @@ function test_closure_instantiation(FT, closurename)
return eltype(closure) == FT
end

function test_calculate_diffusivities(arch, closurename, FT=Float64; kwargs...)
tracernames = (:b,)
closure = getproperty(TurbulenceClosures, closurename)(FT; kwargs...)
closure = with_tracers(tracernames, closure)
grid = RegularCartesianGrid(FT; size=(3, 3, 3), length=(3, 3, 3))
diffusivities = TurbulentDiffusivities(arch, grid, tracernames, closure)
buoyancy = BuoyancyTracer()
velocities = Oceananigans.VelocityFields(arch, grid)
tracers = Oceananigans.TracerFields(arch, grid, tracernames)

U, C, K = datatuples(velocities, tracers, diffusivities)

calculate_diffusivities!(K, arch, grid, closure, buoyancy, U, C)

return true
end

function test_constant_isotropic_diffusivity_basic(T=Float64; ν=T(0.3), κ=T(0.7))
closure = ConstantIsotropicDiffusivity(T; κ=(T=κ, S=κ), ν=ν)
return closure.ν == ν && closure.κ.T == κ
Expand Down Expand Up @@ -96,6 +81,23 @@ function test_anisotropic_diffusivity_fluxdiv(FT=Float64; νh=FT(0.3), κh=FT(0.
∂ⱼ_2ν_Σ₃ⱼ(2, 1, 3, grid, closure, U) == 6νh + 8νv)
end

function test_calculate_diffusivities(arch, closurename, FT=Float64; kwargs...)
tracernames = (:b,)
closure = getproperty(TurbulenceClosures, closurename)(FT; kwargs...)
closure = with_tracers(tracernames, closure)
grid = RegularCartesianGrid(FT; size=(3, 3, 3), length=(3, 3, 3))
diffusivities = TurbulentDiffusivities(arch, grid, tracernames, closure)
buoyancy = BuoyancyTracer()
velocities = Oceananigans.VelocityFields(arch, grid)
tracers = Oceananigans.TracerFields(arch, grid, tracernames)

U, C, K = datatuples(velocities, tracers, diffusivities)

calculate_diffusivities!(K, arch, grid, closure, buoyancy, U, C)

return true
end

function time_step_with_tupled_closure(FT, arch)
closure_tuple = (AnisotropicMinimumDissipation(FT), ConstantAnisotropicDiffusivity(FT))

Expand All @@ -106,6 +108,15 @@ function time_step_with_tupled_closure(FT, arch)
return true
end

function compute_closure_specific_diffusive_cfl(closurename)
grid = RegularCartesianGrid(size=(16, 16, 16), length=(1, 2, 3))
closure = getproperty(TurbulenceClosures, closurename)()
model = Model(grid=grid, closure=closure)
dcfl = DiffusiveCFL(0.1)
dcfl(model)
return true # Just make sure dcfl(model) does not error.
end

@testset "Turbulence closures" begin
println("Testing turbulence closures...")

Expand All @@ -118,18 +129,6 @@ end
end
end

@testset "Calculation of nonlinear diffusivities" begin
println(" Testing calculation of nonlinear diffusivities...")
for T in float_types
for arch in archs
for closure in closures
println(" Calculating diffusivities for $closure [$T, $arch]")
@test test_calculate_diffusivities(arch, closure, T)
end
end
end
end

@testset "Constant isotropic diffusivity" begin
println(" Testing constant isotropic diffusivity...")
for T in float_types
Expand All @@ -146,6 +145,18 @@ end
end
end

@testset "Calculation of nonlinear diffusivities" begin
println(" Testing calculation of nonlinear diffusivities...")
for T in float_types
for arch in archs
for closure in closures
println(" Calculating diffusivities for $closure [$T, $arch]")
@test test_calculate_diffusivities(arch, closure, T)
end
end
end
end

@testset "Closure tuples" begin
println(" Testing time-stepping with a tuple of closures...")
for arch in archs
Expand All @@ -154,4 +165,11 @@ end
end
end
end

@testset "Diagnostics" begin
println(" Testing turbulence closure diagnostics...")
for closure in closures
@test compute_closure_specific_diffusive_cfl(closure)
end
end
end