Model setup
This section describes all the options and features that can be used to set up a model. For more detailed information consult the API documentation.
Each structure covered in this section can be constructed and passed to the models' constructors. For examples of model construction, see the examples. The validation experiments provide more advanced examples.
For reference, here are all the option or keyword arguments that can be passed to the currently implemented models. See the different sections on the sidebar for more details and examples for each keyword argument.
NonhydrostaticModel
Oceananigans.Models.NonhydrostaticModels.NonhydrostaticModel
— TypeNonhydrostaticModel(; grid,
clock = Clock{eltype(grid)}(time = 0),
advection = Centered(),
buoyancy = nothing,
coriolis = nothing,
stokes_drift = nothing,
forcing::NamedTuple = NamedTuple(),
closure = nothing,
boundary_conditions::NamedTuple = NamedTuple(),
tracers = (),
timestepper = :RungeKutta3,
background_fields::NamedTuple = NamedTuple(),
particles::ParticlesOrNothing = nothing,
biogeochemistry::AbstractBGCOrNothing = nothing,
velocities = nothing,
nonhydrostatic_pressure = CenterField(grid),
hydrostatic_pressure_anomaly = DefaultHydrostaticPressureAnomaly(),
diffusivity_fields = nothing,
pressure_solver = nothing,
auxiliary_fields = NamedTuple())
Construct a model for a non-hydrostatic, incompressible fluid on grid
, using the Boussinesq approximation when buoyancy != nothing
. By default, all Bounded directions are rigid and impenetrable.
Keyword arguments
grid
: (required) The resolution and discrete geometry on which themodel
is solved. The architecture (CPU/GPU) that the model is solved on is inferred from the architecture of thegrid
. Note that the grid needs to be regularly spaced in the horizontal dimensions, $x$ and $y$.advection
: The scheme that advects velocities and tracers. SeeOceananigans.Advection
.buoyancy
: The buoyancy model. SeeOceananigans.BuoyancyFormulations
.coriolis
: Parameters for the background rotation rate of the model.stokes_drift
: Parameters for Stokes drift fields associated with surface waves. Default:nothing
.forcing
:NamedTuple
of user-defined forcing functions that contribute to solution tendencies.closure
: The turbulence closure formodel
. SeeOceananigans.TurbulenceClosures
.boundary_conditions
:NamedTuple
containing field boundary conditions.tracers
: A tuple of symbols defining the names of the modeled tracers, or aNamedTuple
of preallocatedCenterField
s.timestepper
: A symbol that specifies the time-stepping method. Either:QuasiAdamsBashforth2
or:RungeKutta3
(default).background_fields
:NamedTuple
with background fields (e.g., background flow). Default:nothing
.particles
: Lagrangian particles to be advected with the flow. Default:nothing
.biogeochemistry
: Biogeochemical model fortracers
.velocities
: The model velocities. Default:nothing
.nonhydrostatic_pressure
: The nonhydrostatic pressure field. Default:CenterField(grid)
.hydrostatic_pressure_anomaly
: An optional field that stores the part of the nonhydrostatic pressure in hydrostatic balance with the buoyancy field. IfCenterField(grid)
(default), the anomaly is precomputed by vertically integrating the buoyancy field. In this case, thenonhydrostatic_pressure
represents only the part of pressure that deviates from the hydrostatic anomaly. Ifnothing
, the anomaly is not computed.diffusivity_fields
: Diffusivity fields. Default:nothing
.pressure_solver
: Pressure solver to be used in the model. Ifnothing
(default), the model constructor chooses the default based on thegrid
provide.auxiliary_fields
:NamedTuple
of auxiliary fields. Default:nothing
HydrostaticFreeSurfaceModel
Oceananigans.Models.HydrostaticFreeSurfaceModels.HydrostaticFreeSurfaceModel
— TypeHydrostaticFreeSurfaceModel(; grid,
clock = Clock{eltype(grid)}(time = 0),
momentum_advection = VectorInvariant(),
tracer_advection = Centered(),
buoyancy = SeawaterBuoyancy(eltype(grid)),
coriolis = nothing,
free_surface = default_free_surface(grid, gravitational_acceleration=g_Earth),
forcing::NamedTuple = NamedTuple(),
closure = nothing,
timestepper = :QuasiAdamsBashforth2,
boundary_conditions::NamedTuple = NamedTuple(),
tracers = (:T, :S),
particles::ParticlesOrNothing = nothing,
biogeochemistry::AbstractBGCOrNothing = nothing,
velocities = nothing,
pressure = nothing,
diffusivity_fields = nothing,
auxiliary_fields = NamedTuple())
Construct a hydrostatic model with a free surface on grid
.
Keyword arguments
grid
: (required) The resolution and discrete geometry on whichmodel
is solved. The architecture (CPU/GPU) that the model is solved is inferred from the architecture of thegrid
.momentum_advection
: The scheme that advects velocities. SeeOceananigans.Advection
.tracer_advection
: The scheme that advects tracers. SeeOceananigans.Advection
.buoyancy
: The buoyancy model. SeeOceananigans.BuoyancyFormulations
.coriolis
: Parameters for the background rotation rate of the model.free_surface
: The free surface model. The default free-surface solver depends on the geometry of thegrid
. If thegrid
is aRectilinearGrid
that is regularly spaced in the horizontal the default is anImplicitFreeSurface
solver withsolver_method = :FFTBasedPoissonSolver
. In all other cases, the default is aSplitExplicitFreeSurface
.tracers
: A tuple of symbols defining the names of the modeled tracers, or aNamedTuple
of preallocatedCenterField
s.forcing
:NamedTuple
of user-defined forcing functions that contribute to solution tendencies.closure
: The turbulence closure formodel
. SeeOceananigans.TurbulenceClosures
.timestepper
: A symbol that specifies the time-stepping method. Either:QuasiAdamsBashforth2
(default) or:SplitRungeKutta3
.boundary_conditions
:NamedTuple
containing field boundary conditions.particles
: Lagrangian particles to be advected with the flow. Default:nothing
.biogeochemistry
: Biogeochemical model fortracers
.velocities
: The model velocities. Default:nothing
.pressure
: Hydrostatic pressure field. Default:nothing
.diffusivity_fields
: Diffusivity fields. Default:nothing
.auxiliary_fields
:NamedTuple
of auxiliary fields. Default:nothing
.
ShallowWaterModel
Oceananigans.Models.ShallowWaterModels.ShallowWaterModel
— TypeShallowWaterModel(; grid,
gravitational_acceleration,
clock = Clock{eltype(grid)}(time = 0),
momentum_advection = UpwindBiased(order=5),
tracer_advection = WENO(),
mass_advection = WENO(),
coriolis = nothing,
forcing::NamedTuple = NamedTuple(),
closure = nothing,
bathymetry = nothing,
tracers = (),
diffusivity_fields = nothing,
boundary_conditions::NamedTuple = NamedTuple(),
timestepper::Symbol = :RungeKutta3,
formulation = ConservativeFormulation())
Construct a shallow water model on grid
with gravitational_acceleration
constant.
Keyword arguments
grid
: (required) The resolution and discrete geometry on whichmodel
is solved. The architecture (CPU/GPU) that the model is solve is inferred from the architecture of the grid.gravitational_acceleration
: (required) The gravitational acceleration constant.clock
: Theclock
for the model.momentum_advection
: The scheme that advects velocities. SeeOceananigans.Advection
. Default:UpwindBiased(order=5)
.tracer_advection
: The scheme that advects tracers. SeeOceananigans.Advection
. Default:WENO()
.mass_advection
: The scheme that advects the mass equation. SeeOceananigans.Advection
. Default:WENO()
.coriolis
: Parameters for the background rotation rate of the model.forcing
:NamedTuple
of user-defined forcing functions that contribute to solution tendencies.closure
: The turbulence closure formodel
. SeeOceananigans.TurbulenceClosures
.bathymetry
: The bottom bathymetry.tracers
: A tuple of symbols defining the names of the modeled tracers, or aNamedTuple
of preallocatedCenterField
s.diffusivity_fields
: Stores diffusivity fields when the closures require a diffusivity to be calculated at each timestep.boundary_conditions
:NamedTuple
containing field boundary conditions.timestepper
: A symbol that specifies the time-stepping method. Either:QuasiAdamsBashforth2
or:RungeKutta3
(default).formulation
: Whether the dynamics are expressed in conservative form (ConservativeFormulation()
; default) or in non-conservative form with a vector-invariant formulation for the non-linear terms (VectorInvariantFormulation()
).
The ConservativeFormulation()
requires RectilinearGrid
. Use VectorInvariantFormulation()
with LatitudeLongitudeGrid
.