-
Notifications
You must be signed in to change notification settings - Fork 189
/
near_global_quarter_degree.jl
320 lines (235 loc) · 11 KB
/
near_global_quarter_degree.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
using Statistics
using JLD2
using Printf
using Plots
using Oceananigans
using Oceananigans.Units
using Oceananigans.Fields: interpolate, Field
using Oceananigans.Advection: VelocityStencil
using Oceananigans.Architectures: arch_array
using Oceananigans.Coriolis: HydrostaticSphericalCoriolis
using Oceananigans.BoundaryConditions
using Oceananigans.ImmersedBoundaries: ImmersedBoundaryGrid, GridFittedBottom, inactive_node, peripheral_node
using CUDA: @allowscalar, device!
using Oceananigans.Operators
using Oceananigans.Operators: Δzᵃᵃᶜ
using Oceananigans: prognostic_fields
device!(3)
@inline function visualize(field, lev, dims)
(dims == 1) && (idx = (lev, :, :))
(dims == 2) && (idx = (:, lev, :))
(dims == 3) && (idx = (:, :, lev))
r = deepcopy(Array(interior(field)))[idx...]
r[ r.==0 ] .= NaN
return r
end
#####
##### Grid
#####
arch = GPU()
reference_density = 1029
latitude = (-75, 75)
# 0.25 degree resolution
Nx = 1440
Ny = 600
Nz = 48
const Nyears = 1
const Nmonths = 12
const thirty_days = 30days
output_prefix = "near_global_lat_lon_$(Nx)_$(Ny)_$(Nz)_fine"
pickup_file = false
#####
##### Load forcing files and inital conditions from ECCO version 4
##### https://ecco.jpl.nasa.gov/drive/files
##### Bathymetry is interpolated from ETOPO1 https://www.ngdc.noaa.gov/mgg/global/
#####
using DataDeps
path = "https://github.com/CliMA/OceananigansArtifacts.jl/raw/ss/new_hydrostatic_data_after_cleared_bugs/quarter_degree_near_global_input_data/"
datanames = ["z_faces-50-levels",
"bathymetry-1440x600",
"temp-1440x600-latitude-75",
"salt-1440x600-latitude-75",
"tau_x-1440x600-latitude-75",
"tau_y-1440x600-latitude-75",
"initial_conditions"]
dh = DataDep("quarter_degree_near_global_lat_lon",
"Forcing data for global latitude longitude simulation",
[path * data * ".jld2" for data in datanames]
)
DataDeps.register(dh)
datadep"quarter_degree_near_global_lat_lon"
files = [:file_z_faces, :file_bathymetry, :file_temp, :file_salt, :file_tau_x, :file_tau_y, :file_init]
for (data, file) in zip(datanames, files)
datadep_path = @datadep_str "quarter_degree_near_global_lat_lon/" * data * ".jld2"
@eval $file = jldopen($datadep_path)
end
bathymetry = file_bathymetry["bathymetry"]
τˣ = zeros(Nx, Ny, Nmonths)
τʸ = zeros(Nx, Ny, Nmonths)
T★ = zeros(Nx, Ny, Nmonths)
S★ = zeros(Nx, Ny, Nmonths)
# Files contain 1 year (1992) of 12 monthly averages
τˣ = file_tau_x["field"] ./ reference_density
τʸ = file_tau_y["field"] ./ reference_density
T★ = file_temp["field"]
S★ = file_salt["field"]
# Remember the convention!! On the surface a negative flux increases a positive decreases
bathymetry = arch_array(arch, bathymetry)
# Stretched faces taken from ECCO Version 4 (50 levels in the vertical)
z_faces = file_z_faces["z_faces"][3:end]
# A spherical domain
@show underlying_grid = LatitudeLongitudeGrid(arch,
size = (Nx, Ny, Nz),
longitude = (-180, 180),
latitude = latitude,
halo = (5, 5, 5),
z = z_faces,
precompute_metrics = true)
grid = ImmersedBoundaryGrid(underlying_grid, GridFittedBottom(bathymetry))
τˣ = arch_array(arch, - τˣ)
τʸ = arch_array(arch, - τʸ)
target_sea_surface_temperature = T★ = arch_array(arch, T★)
target_sea_surface_salinity = S★ = arch_array(arch, S★)
#####
##### Physics and model setup
#####
νz = 5e-3
κz = 1e-4
using Oceananigans.Operators: Δx, Δy
using Oceananigans.TurbulenceClosures
@inline νhb(i, j, k, grid, lx, ly, lz) = (1 / (1 / Δx(i, j, k, grid, lx, ly, lz)^2 + 1 / Δy(i, j, k, grid, lx, ly, lz)^2 ))^2 / 5days
vertical_diffusivity = VerticalScalarDiffusivity(VerticallyImplicitTimeDiscretization(), ν=νz, κ=κz)
convective_adjustment = ConvectiveAdjustmentVerticalDiffusivity(VerticallyImplicitTimeDiscretization(), convective_κz = 1.0)
biharmonic_viscosity = HorizontalDivergenceScalarBiharmonicDiffusivity(ν=νhb, discrete_form=true)
closures = (vertical_diffusivity, convective_adjustment, biharmonic_viscosity)
#####
##### Boundary conditions / time-dependent fluxes
#####
@inline current_time_index(time, tot_months) = mod(unsafe_trunc(Int32, time / thirty_days), tot_months) + 1
@inline next_time_index(time, tot_months) = mod(unsafe_trunc(Int32, time / thirty_days) + 1, tot_months) + 1
@inline cyclic_interpolate(u₁::Number, u₂, time) = u₁ + mod(time / thirty_days, 1) * (u₂ - u₁)
Δz_top = @allowscalar Δzᵃᵃᶜ(1, 1, grid.Nz, grid.underlying_grid)
Δz_bottom = @allowscalar Δzᵃᵃᶜ(1, 1, 1, grid.underlying_grid)
@inline function surface_wind_stress(i, j, grid, clock, fields, τ)
time = clock.time
n₁ = current_time_index(time, Nmonths)
n₂ = next_time_index(time, Nmonths)
@inbounds begin
τ₁ = τ[i, j, n₁]
τ₂ = τ[i, j, n₂]
end
return cyclic_interpolate(τ₁, τ₂, time)
end
u_wind_stress_bc = FluxBoundaryCondition(surface_wind_stress, discrete_form = true, parameters = τˣ)
v_wind_stress_bc = FluxBoundaryCondition(surface_wind_stress, discrete_form = true, parameters = τʸ)
# Linear bottom drag:
μ = 0.001 # ms⁻¹
@inline u_bottom_drag(i, j, grid, clock, fields, μ) = @inbounds - μ * fields.u[i, j, 1]
@inline v_bottom_drag(i, j, grid, clock, fields, μ) = @inbounds - μ * fields.v[i, j, 1]
@inline u_immersed_bottom_drag(i, j, k, grid, clock, fields, μ) = @inbounds - μ * fields.u[i, j, k]
@inline v_immersed_bottom_drag(i, j, k, grid, clock, fields, μ) = @inbounds - μ * fields.v[i, j, k]
u_immersed_bc = ImmersedBoundaryCondition(bottom = FluxBoundaryCondition(u_immersed_bottom_drag, discrete_form = true, parameters = μ))
v_immersed_bc = ImmersedBoundaryCondition(bottom = FluxBoundaryCondition(v_immersed_bottom_drag, discrete_form = true, parameters = μ))
u_bottom_drag_bc = FluxBoundaryCondition(u_bottom_drag, discrete_form = true, parameters = μ)
v_bottom_drag_bc = FluxBoundaryCondition(v_bottom_drag, discrete_form = true, parameters = μ)
@inline function surface_temperature_relaxation(i, j, grid, clock, fields, p)
time = clock.time
n₁ = current_time_index(time, Nmonths)
n₂ = next_time_index(time, Nmonths)
@inbounds begin
T★₁ = p.T★[i, j, n₁]
T★₂ = p.T★[i, j, n₂]
T_surface = fields.T[i, j, grid.Nz]
end
T★ = cyclic_interpolate(T★₁, T★₂, time)
return p.λ * (T_surface - T★)
end
@inline function surface_salinity_relaxation(i, j, grid, clock, fields, p)
time = clock.time
n₁ = current_time_index(time, Nmonths)
n₂ = next_time_index(time, Nmonths)
@inbounds begin
S★₁ = p.S★[i, j, n₁]
S★₂ = p.S★[i, j, n₂]
S_surface = fields.S[i, j, grid.Nz]
end
S★ = cyclic_interpolate(S★₁, S★₂, time)
return p.λ * (S_surface - S★)
end
T_surface_relaxation_bc = FluxBoundaryCondition(surface_temperature_relaxation,
discrete_form = true,
parameters = (λ = Δz_top/7days, T★ = target_sea_surface_temperature))
S_surface_relaxation_bc = FluxBoundaryCondition(surface_salinity_relaxation,
discrete_form = true,
parameters = (λ = Δz_top/7days, S★ = target_sea_surface_salinity))
u_bcs = FieldBoundaryConditions(top = u_wind_stress_bc, bottom = u_bottom_drag_bc, immersed = u_immersed_bc)
v_bcs = FieldBoundaryConditions(top = v_wind_stress_bc, bottom = v_bottom_drag_bc, immersed = v_immersed_bc)
T_bcs = FieldBoundaryConditions(top = T_surface_relaxation_bc)
S_bcs = FieldBoundaryConditions(top = S_surface_relaxation_bc)
free_surface = ImplicitFreeSurface(solver_method=:HeptadiagonalIterativeSolver)
buoyancy = SeawaterBuoyancy(equation_of_state=LinearEquationOfState())
model = HydrostaticFreeSurfaceModel(grid = grid,
free_surface = free_surface,
momentum_advection = WENO(vector_invariant = VelocityStencil()),
coriolis = HydrostaticSphericalCoriolis(),
buoyancy = buoyancy,
tracers = (:T, :S),
closure = (vertical_diffusivity, convective_adjustment, biharmonic_viscosity),
boundary_conditions = (u=u_bcs, v=v_bcs, T=T_bcs, S=S_bcs),
tracer_advection = WENO(underlying_grid))
#####
##### Initial condition:
#####
u, v, w = model.velocities
η = model.free_surface.η
T = model.tracers.T
S = model.tracers.S
@info "Reading initial conditions"
T_init = file_init["T"]
S_init = file_init["S"]
set!(model, T=T_init, S=S_init)
fill_halo_regions!(T)
fill_halo_regions!(S)
@info "model initialized"
#####
##### Simulation setup
#####
Δt = 6minutes # for initialization, then we can go up to 6 minutes?
simulation = Simulation(model, Δt = Δt, stop_time = Nyears*years)
start_time = [time_ns()]
using Oceananigans.Utils
function progress(sim)
wall_time = (time_ns() - start_time[1]) * 1e-9
u = sim.model.velocities.u
η = sim.model.free_surface.η
@info @sprintf("Time: % 12s, iteration: %d, max(|u|): %.2e ms⁻¹, wall time: %s",
prettytime(sim.model.clock.time),
sim.model.clock.iteration, maximum(abs, u),
prettytime(wall_time))
start_time[1] = time_ns()
return nothing
end
simulation.callbacks[:progress] = Callback(progress, IterationInterval(10))
u, v, w = model.velocities
T = model.tracers.T
S = model.tracers.S
η = model.free_surface.η
output_fields = (; u, v, T, S, η)
save_interval = 5days
simulation.output_writers[:surface_fields] = JLD2OutputWriter(model, (; u, v, T, S, η),
schedule = TimeInterval(save_interval),
filename = output_prefix * "_surface",
indices = (:, :, grid.Nz),
overwrite_existing = true)
simulation.output_writers[:checkpointer] = Checkpointer(model,
schedule = TimeInterval(1year),
prefix = output_prefix * "_checkpoint",
overwrite_existing = true)
# Let's goo!
@info "Running with Δt = $(prettytime(simulation.Δt))"
run!(simulation, pickup = pickup_file)
@info """
Simulation took $(prettytime(simulation.run_wall_time))
Free surface: $(typeof(model.free_surface).name.wrapper)
Time step: $(prettytime(Δt))
"""