-
Notifications
You must be signed in to change notification settings - Fork 14
/
Copy pathacc_regional_simulation.jl
219 lines (172 loc) · 7.12 KB
/
acc_regional_simulation.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
using Oceananigans
using Oceananigans.Units
using Oceananigans.Grids: φnode
using ClimaOcean
using ClimaOcean.ECCO
using Printf
using CairoMakie
using CFTime
using Dates
arch = GPU()
Nx = 1440
Ny = 600
Nz = 40
z_faces = exponential_z_faces(Nz=Nz+1, depth=6000)
grid = LatitudeLongitudeGrid(arch;
size = (Nx, Ny, Nz),
halo = (7, 7, 7),
z = z_faces,
latitude = (-80, -20),
longitude = (0, 360))
bottom_height = regrid_bathymetry(grid;
minimum_depth = 10,
interpolation_passes = 5,
major_basins = 1)
grid = ImmersedBoundaryGrid(grid, GridFittedBottom(bottom_height), active_cells_map=true)
dates = DateTimeProlepticGregorian(1993, 1, 1) : Month(1) : DateTimeProlepticGregorian(1993, 12, 1)
#
# Restoring force
# φS φN -20
# -------------- | ------------------ | ------------ |
# no restoring 0 linear mask 1 mask = 1 1
#
const φN₁ = -23
const φN₂ = -25
const φS₁ = -78
const φS₂ = -75
@inline northern_mask(φ) = min(max((φ - φN₂) / (φN₁ - φN₂), zero(φ)), one(φ))
@inline southern_mask(φ, z) = ifelse(z > -20,
min(max((φ - φS₂) / (φS₁ - φS₂), zero(φ)), one(φ)),
zero(φ))
@inline function tracer_mask(λ, φ, z, t)
n = northern_mask(φ)
s = southern_mask(φ, z)
return max(s, n)
end
@inline function u_restoring(i, j, k, grid, clock, fields, p)
φ = φnode(i, j, k, grid, Face(), Center(), Center())
return - p.rate * fields.u[i, j, k] * northern_mask(φ)
end
@inline function v_restoring(i, j, k, grid, clock, fields, p)
φ = φnode(i, j, k, grid, Center(), Face(), Center())
return - p.rate * fields.v[i, j, k] * northern_mask(φ)
end
forcing = (T=ECCORestoring(:temperature, grid; dates, rate=1/5days, mask=tracer_mask),
S=ECCORestoring(:salinity, grid; dates, rate=1/5days, mask=tracer_mask),
u=Forcing(u_restoring; discrete_form=true, parameters=(; rate=1/5days)),
v=Forcing(v_restoring; discrete_form=true, parameters=(; rate=1/5days)))
ocean = ocean_simulation(grid; forcing)
model = ocean.model
set!(model,
T = ECCOMetadata(:temperature; dates = dates[1]),
S = ECCOMetadata(:salinity; dates = dates[1]))
backend = JRA55NetCDFBackend(41)
atmosphere = JRA55PrescribedAtmosphere(arch; backend)
radiation = Radiation()
coupled_model = OceanSeaIceModel(ocean; atmosphere, radiation)
coupled_simulation = Simulation(coupled_model; Δt=10minutes, stop_time = 10days)
wall_time = [time_ns()]
function progress(sim)
ocean = sim.model.ocean
u, v, w = ocean.model.velocities
T = ocean.model.tracers.T
Tmax, Tmin = maximum(interior(T)), minimum(interior(T))
umax = maximum(abs, interior(u)), maximum(abs, interior(v)), maximum(abs, interior(w))
step_time = 1e-9 * (time_ns() - wall_time[1])
@info @sprintf("Time: %s, Iteration %d, Δt %s, max(vel): (%.2e, %.2e, %.2e), max(T): %.2f, min(T): %.2f, wtime: %s \n",
prettytime(ocean.model.clock.time),
ocean.model.clock.iteration,
prettytime(ocean.Δt),
umax..., Tmax, Tmin, prettytime(step_time))
wall_time[1] = time_ns()
end
coupled_simulation.callbacks[:progress] = Callback(progress, TimeInterval(4hours))
ocean.output_writers[:surface] = JLD2OutputWriter(model, merge(model.tracers, model.velocities);
schedule = TimeInterval(5days),
filename = "surface",
indices = (:, :, grid.Nz),
overwrite_existing = true,
array_type = Array{Float32})
nothing #hide
# ### Spinning up the simulation
#
# As an initial condition, we have interpolated ECCO tracer fields onto our custom grid.
# The bathymetry of the original ECCO data may differ from our grid, so the initialization of the velocity
# field might cause shocks if a large time step is used.
#
# Therefore, we spin up the simulation with a small time step to ensure that the interpolated initial
# conditions adapt to the model numerics and parameterization without causing instability. A 10-day
# integration with a maximum time step of 1 minute should be sufficient to dissipate spurious
# initialization shocks.
coupled_simulation.stop_time = 10days
coupled_simulation.Δt = 2minutes
run!(coupled_simulation)
nothing #hide
# ### Running the simulation
#
# Now that the simulation has spun up, we can run it for the full 2 years.
# We increase the maximum time step size to 10 minutes and let the simulation run for 2 years.
coupled_simulation.stop_time = 2*365days
coupled_simulation.Δt = 10minutes
run!(coupled_simulation)
nothing #hide
# ## Visualizing the results
#
# The simulation has finished, let's visualize the results.
# In this section we pull up the saved data and create visualizations using the CairoMakie.jl package.
# In particular, we generate an animation of the evolution of surface fields:
# surface speed (s), surface temperature (T), and turbulent kinetic energy (e).
u = FieldTimeSeries("surface.jld2", "u"; backend = OnDisk())
v = FieldTimeSeries("surface.jld2", "v"; backend = OnDisk())
T = FieldTimeSeries("surface.jld2", "T"; backend = OnDisk())
e = FieldTimeSeries("surface.jld2", "e"; backend = OnDisk())
times = u.times
Nt = length(times)
iter = Observable(Nt)
Ti = @lift begin
Ti = interior(T[$iter], :, :, 1)
Ti[Ti .== 0] .= NaN
Ti
end
ei = @lift begin
ei = interior(e[$iter], :, :, 1)
ei[ei .== 0] .= NaN
ei
end
si = @lift begin
s = Field(sqrt(u[$iter]^2 + v[$iter]^2))
compute!(s)
s = interior(s, :, :, 1)
s[s .== 0] .= NaN
s
end
fig = Figure(size = (800, 400))
ax = Axis(fig[1, 1])
hm = heatmap!(ax, si, colorrange = (0, 0.5), colormap = :deep)
cb = Colorbar(fig[0, 1], hm, vertical = false, label = "Surface speed (ms⁻¹)")
hidedecorations!(ax)
CairoMakie.record(fig, "near_global_ocean_surface_s.mp4", 1:Nt, framerate = 8) do i
iter[] = i
end
nothing #hide
# ![](near_global_ocean_surface_s.mp4)
fig = Figure(size = (800, 400))
ax = Axis(fig[1, 1])
hm = heatmap!(ax, Ti, colorrange = (-1, 30), colormap = :magma)
cb = Colorbar(fig[0, 1], hm, vertical = false, label = "Surface Temperature (Cᵒ)")
hidedecorations!(ax)
CairoMakie.record(fig, "near_global_ocean_surface_T.mp4", 1:Nt, framerate = 8) do i
iter[] = i
end
nothing #hide
# ![](near_global_ocean_surface_T.mp4)
fig = Figure(size = (800, 400))
ax = Axis(fig[1, 1])
hm = heatmap!(ax, ei, colorrange = (0, 1e-3), colormap = :solar)
cb = Colorbar(fig[0, 1], hm, vertical = false, label = "Turbulent Kinetic Energy (m²s⁻²)")
hidedecorations!(ax)
CairoMakie.record(fig, "near_global_ocean_surface_e.mp4", 1:Nt, framerate = 8) do i
iter[] = i
end
nothing #hide
# ![](near_global_ocean_surface_e.mp4)