/
Oceananigans.jl
300 lines (243 loc) 路 7.89 KB
/
Oceananigans.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
"""
Main module for `Oceananigans.jl` -- a Julia software for fast, friendly, flexible,
data-driven, ocean-flavored fluid dynamics on CPUs and GPUs.
"""
module Oceananigans
if VERSION < v"1.8"
@warn "Oceananigans is tested on Julia v1.8 and therefore it is strongly recommended you run Oceananigans on Julia v1.8 or newer."
end
export
# Architectures
CPU, GPU,
# Logging
OceananigansLogger,
# Grids
Center, Face,
Periodic, Bounded, Flat,
FullyConnected, LeftConnected, RightConnected,
RectilinearGrid,
LatitudeLongitudeGrid,
ConformalCubedSphereFaceGrid,
xnodes, ynodes, znodes, nodes,
# Immersed boundaries
ImmersedBoundaryGrid, GridFittedBoundary, GridFittedBottom, ImmersedBoundaryCondition,
# Advection schemes
Centered, CenteredSecondOrder, CenteredFourthOrder,
UpwindBiased, UpwindBiasedFirstOrder, UpwindBiasedThirdOrder, UpwindBiasedFifthOrder,
WENO, WENOThirdOrder, WENOFifthOrder,
VectorInvariant, EnergyConservingScheme, EnstrophyConservingScheme,
# Boundary conditions
BoundaryCondition,
FluxBoundaryCondition, ValueBoundaryCondition, GradientBoundaryCondition, OpenBoundaryCondition,
FieldBoundaryConditions,
# Fields and field manipulation
Field, CenterField, XFaceField, YFaceField, ZFaceField,
Average, Integral, Reduction, BackgroundField,
interior, set!, compute!, regrid!, location,
# Forcing functions
Forcing, Relaxation, LinearTarget, GaussianMask, AdvectiveForcing,
# Coriolis forces
FPlane, ConstantCartesianCoriolis, BetaPlane, NonTraditionalBetaPlane,
# BuoyancyModels and equations of state
Buoyancy, BuoyancyTracer, SeawaterBuoyancy,
LinearEquationOfState, TEOS10,
BuoyancyField,
# Surface wave Stokes drift via Craik-Leibovich equations
UniformStokesDrift,
# Turbulence closures
VerticalScalarDiffusivity,
HorizontalScalarDiffusivity,
ScalarDiffusivity,
VerticalScalarBiharmonicDiffusivity,
HorizontalScalarBiharmonicDiffusivity,
ScalarBiharmonicDiffusivity,
SmagorinskyLilly,
AnisotropicMinimumDissipation,
ConvectiveAdjustmentVerticalDiffusivity,
RiBasedVerticalDiffusivity,
IsopycnalSkewSymmetricDiffusivity,
FluxTapering,
VerticallyImplicitTimeDiscretization,
# Lagrangian particle tracking
LagrangianParticles,
# Models
NonhydrostaticModel,
HydrostaticFreeSurfaceModel,
ShallowWaterModel, ConservativeFormulation, VectorInvariantFormulation,
PressureField,
fields,
# Hydrostatic free surface model stuff
VectorInvariant, ExplicitFreeSurface, ImplicitFreeSurface, SplitExplicitFreeSurface,
HydrostaticSphericalCoriolis,
PrescribedVelocityFields,
# Time stepping
Clock, TimeStepWizard, time_step!,
# Simulations
Simulation, run!, Callback, iteration, stopwatch,
iteration_limit_exceeded, stop_time_exceeded, wall_time_limit_exceeded,
erroring_NaNChecker!,
TimeStepCallsite, TendencyCallsite, UpdateStateCallsite,
# Diagnostics
StateChecker, CFL, AdvectiveCFL, DiffusiveCFL,
# Output writers
NetCDFOutputWriter, JLD2OutputWriter, Checkpointer,
TimeInterval, IterationInterval, AveragedTimeInterval, SpecifiedTimes,
AndSchedule, OrSchedule,
# Output readers
FieldTimeSeries, FieldDataset, InMemory, OnDisk,
# Abstract operations
鈭倄, 鈭倅, 鈭倆, @at, KernelFunctionOperation,
# MultiRegion and Cubed sphere
MultiRegionGrid, XPartition,
ConformalCubedSphereGrid,
# Utils
prettytime, apply_regionally!, construct_regionally, @apply_regionally, MultiRegionObject
using Printf
using Logging
using Statistics
using LinearAlgebra
using CUDA
using Adapt
using DocStringExtensions
using OffsetArrays
using FFTW
using JLD2
using Base: @propagate_inbounds
using Statistics: mean
import Base:
+, -, *, /,
size, length, eltype,
iterate, similar, show,
getindex, lastindex, setindex!,
push!
# TODO: find a way to check whether the libraries for AMGX and NETCDF
# (libamgxsh and libnetcdf, respectively) are installed on the machine
"Boolean denoting whether AMGX.jl can be loaded on machine."
const hasamgx = @static (Sys.islinux() && Sys.ARCH == :x86_64) ? true : false
"Boolean denoting whether NCDatasets.jl can be loaded on machine."
const hasnetcdf = @static (Sys.islinux() && Sys.ARCH == :x86_64) ? true : false
"""
@ifhasamgx expr
Evaluate `expr` only if `hasamgx == true`.
"""
macro ifhasamgx(expr)
hasamgx ? :($(esc(expr))) : :(nothing)
end
"""
@ifnetcdf expr
Evaluate `expr` only if `hasnetcdf == true`.
"""
macro ifhasnetcdf(expr)
hasnetcdf ? :($(esc(expr))) : :(nothing)
end
#####
##### Abstract types
#####
"""
AbstractModel
Abstract supertype for models.
"""
abstract type AbstractModel{TS} end
"""
AbstractDiagnostic
Abstract supertype for diagnostics that compute information from the current
model state.
"""
abstract type AbstractDiagnostic end
"""
AbstractOutputWriter
Abstract supertype for output writers that write data to disk.
"""
abstract type AbstractOutputWriter end
struct TimeStepCallsite end
struct TendencyCallsite end
struct UpdateStateCallsite end
#####
##### Place-holder functions
#####
function run_diagnostic! end
function write_output! end
function location end
function instantiated_location end
function tupleit end
function fields end
function prognostic_fields end
function tracer_tendency_kernel_function end
#####
##### Include all the submodules
#####
# Basics
include("Architectures.jl")
include("Units.jl")
include("Grids/Grids.jl")
include("Utils/Utils.jl")
include("Logger.jl")
include("Operators/Operators.jl")
include("BoundaryConditions/BoundaryConditions.jl")
include("Fields/Fields.jl")
include("AbstractOperations/AbstractOperations.jl")
include("Advection/Advection.jl")
include("Solvers/Solvers.jl")
include("Distributed/Distributed.jl")
# Physics, time-stepping, and models
include("Coriolis/Coriolis.jl")
include("BuoyancyModels/BuoyancyModels.jl")
include("StokesDrift.jl")
include("TurbulenceClosures/TurbulenceClosures.jl")
include("Forcings/Forcings.jl")
include("ImmersedBoundaries/ImmersedBoundaries.jl")
include("LagrangianParticleTracking/LagrangianParticleTracking.jl")
include("TimeSteppers/TimeSteppers.jl")
include("Models/Models.jl")
# Output and Physics, time-stepping, and models
include("Diagnostics/Diagnostics.jl")
include("OutputWriters/OutputWriters.jl")
include("OutputReaders/OutputReaders.jl")
include("Simulations/Simulations.jl")
# Abstractions for distributed and multi-region models
include("MultiRegion/MultiRegion.jl")
include("CubedSpheres/CubedSpheres.jl")
#####
##### Needed so we can export names from sub-modules at the top-level
#####
using .Logger
using .Architectures
using .Utils
using .Advection
using .Grids
using .BoundaryConditions
using .Fields
using .Coriolis
using .BuoyancyModels
using .StokesDrift
using .TurbulenceClosures
using .LagrangianParticleTracking
using .Solvers
using .Forcings
using .ImmersedBoundaries
using .Distributed
using .Models
using .TimeSteppers
using .Diagnostics
using .OutputWriters
using .OutputReaders
using .Simulations
using .AbstractOperations
using .MultiRegion
using .CubedSpheres
function __init__()
threads = Threads.nthreads()
if threads > 1
@info "Oceananigans will use $threads threads"
# See: https://github.com/CliMA/Oceananigans.jl/issues/1113
FFTW.set_num_threads(4*threads)
end
if CUDA.has_cuda()
@debug "CUDA-enabled GPU(s) detected:"
for (gpu, dev) in enumerate(CUDA.devices())
@debug "$dev: $(CUDA.name(dev))"
end
CUDA.allowscalar(false)
end
end
end # module