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

Allow users to specify the domain when creating a `RegularCartes… #464

Merged
merged 15 commits into from
Oct 23, 2019
Merged
Show file tree
Hide file tree
Changes from 4 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
51 changes: 31 additions & 20 deletions src/grids.jl
Original file line number Diff line number Diff line change
Expand Up @@ -71,18 +71,29 @@ RegularCartesianGrid{Float32}
domain (Lx, Ly, Lz) = (8.0f0, 8.0f0, 2.0f0)
grid spacing (Δx, Δy, Δz) = (0.25f0, 0.25f0, 0.125f0)
"""
function RegularCartesianGrid(T, N, L)
length(N) == 3 || throw(ArgumentError("N=$N must be a tuple of length 3."))
length(L) == 3 || throw(ArgumentError("L=$L must be a tuple of length 3."))

function RegularCartesianGrid(T; N, x, y, z)
length(N) == 3 || throw(ArgumentError("N=$N must be a tuple of length 3."))
all(isa.(N, Integer)) || throw(ArgumentError("N=$N should contain integers."))
all(isa.(L, Number)) || throw(ArgumentError("L=$L should contain numbers."))

all(N .>= 1) || throw(ArgumentError("N=$N must be nonzero and positive!"))
all(L .> 0) || throw(ArgumentError("L=$L must be nonzero and positive!"))

all(N .>= 1) || throw(ArgumentError("N=$N must be nonzero and positive!"))

function coord2xyz(c)
c == 1 && return "x"
c == 2 && return "y"
c == 3 && return "z"
end

for (i, c) in enumerate((x, y, z))
name = coord2xyz(i)
length(c) == 2 || throw(ArgumentError("$name=$c must be a tuple of length 2."))
ali-ramadhan marked this conversation as resolved.
Show resolved Hide resolved
all(isa.(c, Number)) || throw(ArgumentError("$name=$c should contain numbers."))
c[2] >= c[1] || throw(ArgumentError("$name=$c should be an increasing interval."))
end

x₁, x₂ = x[1], x[2]
y₁, y₂ = y[1], y[2]
z₁, z₂ = z[1], z[2]
Nx, Ny, Nz = N
Lx, Ly, Lz = L
Lx, Ly, Lz = x₂-x₁, y₂-y₁, z₂-z₁

# Right now we only support periodic horizontal boundary conditions and
# usually use second-order advection schemes so halos of size Hx, Hy = 1 are
Expand All @@ -103,29 +114,29 @@ function RegularCartesianGrid(T, N, L)

V = Δx*Δy*Δz

xC = Δx/2:Δx:Lx
yC = Δy/2:Δy:Ly
zC = -Δz/2:-Δz:-Lz
xC = (x₁ + Δx/2) : Δx : x₂
yC = (y₁ + Δy/2) : Δy : y₂
zC = (z₂ - Δz/2) : -Δz : z₁

xF = 0:Δx:Lx
yF = 0:Δy:Ly
zF = 0:-Δz:-Lz
xF = x₁ : Δx : x₂
yF = y₁ : Δy : y₂
zF = z₂ : -Δz : z₁

RegularCartesianGrid{T, typeof(xC)}(Nx, Ny, Nz, Hx, Hy, Hz, Tx, Ty, Tz,
Lx, Ly, Lz, Δx, Δy, Δz, Ax, Ay, Az, V,
xC, yC, zC, xF, yF, zF)
end

RegularCartesianGrid(N, L) = RegularCartesianGrid(Float64, N, L)
RegularCartesianGrid(T, N, L) = RegularCartesianGrid(T; N=N, x=(0, L[1]), y=(0, L[2]), z=(-L[3], 0))

RegularCartesianGrid(T=Float64; N, L) = RegularCartesianGrid(T, N, L)
RegularCartesianGrid(N, L) = RegularCartesianGrid(Float64, N, L)

size(g::RegularCartesianGrid) = (g.Nx, g.Ny, g.Nz)
eltype(g::RegularCartesianGrid{T}) where T = T

show(io::IO, g::RegularCartesianGrid) =
print(io, "RegularCartesianGrid{$(eltype(g))}\n",
print(io, "RegularCartesianGrid{$(eltype(g))}\n",
"domain: x ∈ [$(g.xF[1]), $(g.xF[end])], y ∈ [$(g.yF[1]), $(g.yF[end])], z ∈ [$(g.zF[end]), $(g.zF[1])]", '\n',
" resolution (Nx, Ny, Nz) = ", (g.Nx, g.Ny, g.Nz), '\n',
" halo size (Hx, Hy, Hz) = ", (g.Hx, g.Hy, g.Hz), '\n',
" domain (Lx, Ly, Lz) = ", (g.Lx, g.Ly, g.Lz), '\n',
"grid spacing (Δx, Δy, Δz) = ", (g.Δx, g.Δy, g.Δz))
8 changes: 4 additions & 4 deletions test/test_coriolis_buoyancy.jl
Original file line number Diff line number Diff line change
Expand Up @@ -20,28 +20,28 @@ function instantiate_seawater_buoyancy(T, EquationOfState)
end

function density_perturbation_works(arch, T, eos)
grid = RegularCartesianGrid(T, N=(3, 3, 3), L=(1, 1, 1))
grid = RegularCartesianGrid(T, (3, 3, 3), (1, 1, 1))
C = datatuple(TracerFields(arch, grid))
density_anomaly = ρ′(2, 2, 2, grid, eos, C)
return true
end

function buoyancy_frequency_squared_works(arch, T, buoyancy)
grid = RegularCartesianGrid(N=(3, 3, 3), L=(1, 1, 1))
grid = RegularCartesianGrid((3, 3, 3), (1, 1, 1))
C = datatuple(TracerFields(arch, grid))
N² = buoyancy_frequency_squared(2, 2, 2, grid, buoyancy, C)
return true
end

function thermal_expansion_works(arch, T, eos)
grid = RegularCartesianGrid(T, N=(3, 3, 3), L=(1, 1, 1))
grid = RegularCartesianGrid(T, (3, 3, 3), (1, 1, 1))
C = datatuple(TracerFields(arch, grid))
α = thermal_expansion(2, 2, 2, grid, eos, C)
return true
end

function haline_contraction_works(arch, T, eos)
grid = RegularCartesianGrid(N=(3, 3, 3), L=(1, 1, 1))
grid = RegularCartesianGrid((3, 3, 3), (1, 1, 1))
C = datatuple(TracerFields(arch, grid))
β = haline_contraction(2, 2, 2, grid, eos, C)
return true
Expand Down
2 changes: 1 addition & 1 deletion test/test_dynamics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -169,7 +169,7 @@ function pearson_vortex_test(arch; FT=Float64, N=64, Nt=10)
bottom = BoundaryCondition(Gradient, 0))

model = Model( architecture = arch,
grid = RegularCartesianGrid(FT; N=(Nx, Ny, Nz), L=(Lx, Ly, Lz)),
grid = RegularCartesianGrid(FT, (Nx, Ny, Nz), (Lx, Ly, Lz)),
closure = ConstantIsotropicDiffusivity(FT; ν=1, κ=0), # Turn off diffusivity.
buoyancy = nothing, # turn off buoyancy
boundary_conditions = BoundaryConditions(u=ubcs, v=vbcs))
Expand Down
4 changes: 0 additions & 4 deletions test/test_grids.jl
Original file line number Diff line number Diff line change
Expand Up @@ -50,10 +50,6 @@ end
@test_throws ArgumentError RegularCartesianGrid(FT, (32,), L)
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 64), L)
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32, 16), L)
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32), (100,))
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32), (100, 100))
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32), (100, 100, 1, 1))
@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32), (100, 100, -100))

@test_throws ArgumentError RegularCartesianGrid(FT, (32, 32, 32.0), (1, 1, 1))
@test_throws ArgumentError RegularCartesianGrid(FT, (20.1, 32, 32), (1, 1, 1))
Expand Down
2 changes: 1 addition & 1 deletion test/test_regression.jl
Original file line number Diff line number Diff line change
Expand Up @@ -96,7 +96,7 @@ function run_rayleigh_benard_regression_test(arch)

model = Model(
architecture = arch,
grid = RegularCartesianGrid(N=(Nx, Ny, Nz), L=(Lx, Ly, Lz)),
grid = RegularCartesianGrid((Nx, Ny, Nz), (Lx, Ly, Lz)),
closure = ConstantIsotropicDiffusivity(ν=ν, κ=κ),
buoyancy = BuoyancyTracer(),
boundary_conditions = BoundaryConditions(T=HorizontallyPeriodicBCs(
Expand Down
2 changes: 1 addition & 1 deletion test/test_time_stepping.jl
Original file line number Diff line number Diff line change
Expand Up @@ -114,7 +114,7 @@ function tracer_conserved_in_channel(arch, FT, Nt)
νv, κv = α*νh, α*κh

model = ChannelModel(architecture = arch, float_type = FT,
grid = RegularCartesianGrid(N = (Nx, Ny, Nz), L = (Lx, Ly, Lz)),
grid = RegularCartesianGrid((Nx, Ny, Nz), (Lx, Ly, Lz)),
closure = ConstantAnisotropicDiffusivity(νh=νh, νv=νv, κh=κh, κv=κv))

Ty = 1e-4 # Meridional temperature gradient [K/m].
Expand Down