Library
Documenting the public user interface.
Boundary conditions
Oceananigans.BoundaryCondition
— TypeBoundaryCondition{C<:BCType}(condition)
Construct a boundary condition of type C
with a condition
that may be given by a number, an array, or a function with signature:
condition(i, j, grid, time, iteration, U, Φ, parameters) = # function definition
that returns a number and where i
and j
are indices along the boundary.
Boundary condition types include Periodic
, Flux
, Value
, Gradient
, and NoPenetration
.
Oceananigans.BoundaryFunction
— TypeBoundaryFunction{B, X1, X2, F}
A wrapper for user-defined boundary condition functions.
Oceananigans.CoordinateBoundaryConditions
— TypeCoordinateBoundaryConditions(left, right)
A set of two BoundaryCondition
s to be applied along a coordinate x, y, or z.
The left
boundary condition is applied on the negative or lower side of the coordinate while the right
boundary condition is applied on the positive or higher side.
Oceananigans.Dirchlet
— TypeDirchlet
An alias for the Value
boundary condition type.
Oceananigans.FieldBoundaryConditions
— TypeFieldBoundaryConditions
An alias for NamedTuple{(:x, :y, :z)}
that represents a set of three CoordinateBoundaryCondition
s applied to a field along x, y, and z.
Oceananigans.FieldBoundaryConditions
— MethodFieldBoundaryConditions(x, y, z)
Construct a FieldBoundaryConditions
using a CoordinateBoundaryCondition
for each of the x
, y
, and z
coordinates.
Oceananigans.Flux
— TypeFlux
A type specifying a boundary condition on the flux of a field.
Oceananigans.Gradient
— TypeGradient
A type specifying a boundary condition on the derivative or gradient of a field. Also called a Neumann boundary condition.
Oceananigans.Neumann
— TypeNeumann
An alias for the Gradient
boundary condition type.
Oceananigans.Periodic
— TypePeriodic
A type specifying a periodic boundary condition.
A condition may not be specified with a Periodic
boundary condition.
Oceananigans.Value
— TypeValue
A type specifying a boundary condition on the value of a field. Also called a Dirchlet boundary condition.
Oceananigans.ChannelBCs
— MethodChannelBCs(; north = BoundaryCondition(Flux, nothing),
south = BoundaryCondition(Flux, nothing),
top = BoundaryCondition(Flux, nothing),
bottom = BoundaryCondition(Flux, nothing))
Construct FieldBoundaryConditions
with Periodic
boundary conditions in the x direction and specified north
(+y), south
(-y), top
(+z) and bottom
(-z) boundary conditions for u, v, and tracer fields.
ChannelBCs
cannot be applied to the the vertical velocity w.
Oceananigans.ChannelSolutionBCs
— MethodChannelSolutionBCs(u=ChannelBCs(), ...)
Construct SolutionBoundaryConditions
for a reentrant channel model configuration with solution fields u
, v
, w
, T
, and S
specified by keyword arguments.
By default ChannelBCs
are applied to u
, v
, T
, and S
and ChannelBCs(top=NoPenetrationBC(), bottom=NoPenetrationBC())
is applied to w
.
Use ChannelBCs
when constructing non-default boundary conditions for u
, v
, w
, T
, S
.
Oceananigans.HorizontallyPeriodicBCs
— MethodHorizontallyPeriodicBCs(; top = BoundaryCondition(Flux, nothing),
bottom = BoundaryCondition(Flux, nothing))
Construct FieldBoundaryConditions
with Periodic
boundary conditions in the x and y directions and specified top
(+z) and bottom
(-z) boundary conditions for u, v, and tracer fields.
HorizontallyPeriodicBCs
cannot be applied to the the vertical velocity w.
Oceananigans.HorizontallyPeriodicSolutionBCs
— MethodHorizontallyPeriodicSolutionBCs(u=HorizontallyPeriodicBCs(), ...)
Construct SolutionBoundaryConditions
for a horizontally-periodic model configuration with solution fields u
, v
, w
, T
, and S
specified by keyword arguments.
By default HorizontallyPeriodicBCs
are applied to u
, v
, T
, and S
and HorizontallyPeriodicBCs(top=NoPenetrationBC(), bottom=NoPenetrationBC())
is applied to w
.
Use HorizontallyPeriodicBCs
when constructing non-default boundary conditions for u
, v
, w
, T
, S
.
Oceananigans.SolutionBoundaryConditions
— MethodSolutionBoundaryConditions(tracers, proposal_bcs)
Construct a NamedTuple
of FieldBoundaryConditions
for a model with fields u
, v
, w
, and tracers
from the proposal boundary conditions proposal_bcs
, which must contain the boundary conditions on u
, v
, and w
and may contain some or all of the boundary conditions on tracers
.
Buoyancy
Oceananigans.BuoyancyTracer
— TypeBuoyancyTracer <: AbstractBuoyancy{Nothing}
Type indicating that the tracer b
represents buoyancy.
Oceananigans.LinearEquationOfState
— TypeLinearEquationOfState{FT} <: AbstractEquationOfState
Linear equation of state for seawater.
Oceananigans.LinearEquationOfState
— TypeLinearEquationOfState([FT=Float64;] α=1.67e-4, β=7.80e-4)
Returns parameters for a linear equation of state for seawater with thermal expansion coefficient α
[K⁻¹] and haline contraction coefficient β
[ppt⁻¹]. The buoyancy perturbation associated with a linear equation of state is
Default constants are taken from Table 1.2 (page 33) of Vallis, "Atmospheric and Oceanic Fluid Dynamics: Fundamentals and Large-Scale Circulation" (2ed, 2017).
Oceananigans.SeawaterBuoyancy
— TypeSeawaterBuoyancy{G, EOS} <: AbstractBuoyancy{EOS}
Buoyancy model for temperature- and salt-stratified seawater.
Oceananigans.SeawaterBuoyancy
— TypeSeawaterBuoyancy([FT=Float64;] gravitational_acceleration = g_Earth,
equation_of_state = LinearEquationOfState(FT))
Returns parameters for a temperature- and salt-stratified seawater buoyancy model with a gravitational_acceleration
constant (typically called 'g'), and an equation_of_state
that related temperature and salinity (or conservative temperature and absolute salinity) to density anomalies and buoyancy.
Clock
Oceananigans.Clock
— TypeClock{T<:Number}
Clock{T}(time, iteration)
Keeps track of the current time
and iteration
number.
Coriolis
Oceananigans.BetaPlane
— TypeBetaPlane{T} <: AbstractRotation
A parameter object for meridionally increasing Coriolis parameter (f = f₀ + βy
).
Oceananigans.BetaPlane
— TypeBetaPlane([T=Float64;] f₀=nothing, β=nothing,
rotation_rate=nothing, latitude=nothing, radius=nothing)
A parameter object for meridionally increasing Coriolis parameter (f = f₀ + βy
).
The user may specify both f₀
and β
, or the three parameters rotation_rate
, latitude
, and radius
that specify the rotation rate and radius of a planet, and the central latitude at which the β
-plane approximation is to be made.
Oceananigans.FPlane
— TypeFPlane{FT} <: AbstractRotation
A parameter object for constant rotation around a vertical axis.
Oceananigans.FPlane
— TypeFPlane([FT=Float64;] f=nothing, rotation_rate=nothing, latitude=nothing)
Returns a parameter object for constant rotation at the angular frequency f/2
, and therefore with background vorticity f
, around a vertical axis. If f
is not specified, it is calculated from rotation_rate
and latitude
according to the relation `f = 2rotation_ratesind(latitude).
Also called FPlane
, after the "f-plane" approximation for the local effect of Earth's rotation in a planar coordinate system tangent to the Earth's surface.
Diagnostics
Fields
Oceananigans.Field
— TypeField{X, Y, Z, A, G} <: AbstractLocatedField{X, Y, Z, A, G}
A field defined at the location (X
, Y
, Z
) which can be either Cell
or Face
.
Oceananigans.Field
— MethodField(X, Y, Z, arch::AbstractArchitecture, grid)
Construct a Field
on architecture arch
and grid
at location X
, Y
, Z
, where each of X, Y, Z
is Cell
or Face
.
Oceananigans.Field
— MethodField(L::Tuple, data::AbstractArray, grid)
Construct a Field
on grid
using the array data
with location defined by the tuple L
of length 3 whose elements are Cell
or Face
.
Oceananigans.Field
— MethodField(L::Tuple, arch::AbstractArchitecture, grid)
Construct a Field
on architecture arch
and grid
at location L
, where L
is a tuple of Cell
or Face
types.
Oceananigans.CellField
— MethodCellField([T=eltype(grid)], arch, grid)
Return a Field{Cell, Cell, Cell}
on architecture arch
and grid
. Used for tracers and pressure fields.
Oceananigans.FaceFieldX
— MethodFaceFieldX([T=eltype(grid)], arch, grid)
Return a Field{Face, Cell, Cell}
on architecture arch
and grid
. Used for the x-velocity field.
Oceananigans.FaceFieldY
— MethodFaceFieldY([T=eltype(grid)], arch, grid)
Return a Field{Cell, Face, Cell}
on architecture arch
and grid
. Used for the y-velocity field.
Oceananigans.FaceFieldZ
— MethodFaceFieldZ([T=eltype(grid)], arch, grid)
Return a Field{Cell, Cell, Face}
on architecture arch
and grid
. Used for the z-velocity field.
Oceananigans.interior
— MethodReturns a view over the interior points of the field.data
.
Oceananigans.set!
— MethodSet the CPU field u
to the array v
.
Oceananigans.set!
— MethodSet the CPU field u
data to the function f(x, y, z)
.
Oceananigans.set!
— Methodset!(model; kwargs...)
Set velocity and tracer fields of model
. The keyword arguments kwargs...
take the form name=data
, where name
refers to one of the fields of model.velocities
or model.tracers
, and the data
may be an array, a function with arguments (x, y, z)
, or any data type for which a set!(ϕ::AbstractField, data)
function exists.
Example
model = Model(grid=RegularCartesianGrid(size=(32, 32, 32), length=(1, 1, 1))
# Set u to a parabolic function of z, v to random numbers damped
# at top and bottom, and T to some silly array of half zeros,
# half random numbers.
u₀(x, y, z) = z/model.grid.Lz * (1 + z/model.grid.Lz)
v₀(x, y, z) = 1e-3 * rand() * u₀(x, y, z)
T₀ = rand(size(model.grid)...)
T₀[T₀ .< 0.5] .= 0
set!(model, u=u₀, v=v₀, T=T₀)
Forcing
Oceananigans.SimpleForcing
— TypeSimpleForcing{X, Y, Z, F, P}
Callable object for specifying 'simple' forcings of x, y, z, t
and optionally parameters
of type P
at location X, Y, Z
.
Oceananigans.SimpleForcing
— MethodSimpleForcing([location=(Cell, Cell, Cell),] forcing; parameters=nothing)
Construct forcing for a field at location
using forcing::Function
, and optionally with parameters
. If parameters=nothing
, forcing
must have the signature
`forcing(x, y, z, t)`;
otherwise it must have the signature
`forcing(x, y, z, t, parameters)`.
Examples
julia> const a = 2.1
julia> fun_forcing(x, y, z, t) = a * exp(z) * cos(t)
julia> u_forcing = SimpleForcing(fun_forcing)
julia> parameterized_forcing(x, y, z, t, p) = p.μ * exp(z/p.λ) * cos(p.ω*t)
julia> v_forcing = SimpleForcing(parameterized_forcing, parameters=(μ=42, λ=0.1, ω=π))
Oceananigans.ModelForcing
— MethodModelForcing(; u=zeroforcing, v=zeroforcing, w=zeroforcing, tracer_forcings...)
Return a named tuple of forcing functions for each solution field.
Example
julia> u_forcing = SimpleForcing((x, y, z, t) -> exp(z) * cos(t))
julia> model = Model(forcing=ModelForcing(u=u_forcing))
Grids
Models
Oceananigans.Model
— MethodModel(; grid, kwargs...)
Construct an Oceananigans.jl
model on grid
.
Keyword arguments
grid
: (required) The resolution and discrete geometry on whichmodel
is solved. Currently the only option isRegularCartesianGrid
.architecture
:CPU()
orGPU()
. The computer architecture used to time-stepmodel
.float_type
:Float32
orFloat64
. The floating point type used formodel
data.closure
: The turbulence closure formodel
. SeeTurbulenceClosures
.buoyancy
: Buoyancy model parameters.coriolis
: Parameters for the background rotation rate of the model.forcing
: User-defined forcing functions that contribute to solution tendencies.boundary_conditions
: User-defined boundary conditions for model fields. Can be eitherSolutionBoundaryConditions
orModelBoundaryConditions
. SeeBoundaryConditions
,HorizontallyPeriodicSolutionBCs
, andChannelSolutionBCs
.parameters
: User-defined parameters for use in user-defined forcing functions and boundary condition functions.
Oceananigans.ChannelModel
— MethodChannelModel(; kwargs...)
Construct a Model
with walls in the y-direction. This is done by imposing FreeSlip
boundary conditions in the y-direction instead of Periodic
.
kwargs are passed to the regular Model
constructor.
Oceananigans.NonDimensionalModel
— MethodNonDimensionalModel(; N, L, Re, Pr=0.7, Ro=Inf, float_type=Float64, kwargs...)
Construct a "Non-dimensional" Model
with resolution N
, domain extent L
, precision float_type
, and the four non-dimensional numbers:
* `Re = U λ / ν` (Reynolds number)
* `Pr = U λ / κ` (Prandtl number)
* `Ro = U / f λ` (Rossby number)
for characteristic velocity scale U
, length-scale λ
, viscosity ν
, tracer diffusivity κ
, and Coriolis parameter f
. Buoyancy is scaled with λ U²
, so that the Richardson number is Ri=B
, where B
is a non-dimensional buoyancy scale set by the user via initial conditions or forcing.
Note that N
, L
, and Re
are required.
Additional kwargs
are passed to the regular Model
constructor.
Output writers
Time steppers
Tubrulence closures
Utilities
Oceananigans.GiB
— ConstantGiB
A Float64
constant equal to 1024MiB
. Useful for increasing the clarity of scripts, e.g. max_filesize = 50GiB
.
Oceananigans.KiB
— ConstantKiB
A Float64
constant equal to 1024.0. Useful for increasing the clarity of scripts, e.g. max_filesize = 250KiB
.
Oceananigans.MiB
— ConstantMiB
A Float64
constant equal to 1024KiB
. Useful for increasing the clarity of scripts, e.g. max_filesize = 100MiB
.
Oceananigans.TiB
— ConstantTiB
A Float64
constant equal to 1024GiB
. Useful for increasing the clarity of scripts, e.g. max_filesize = 2TiB
.
Oceananigans.day
— Constantday
A Float64
constant equal to 24hour
. Useful for increasing the clarity of scripts, e.g. Δt = 0.5day
.
Oceananigans.hour
— Constanthour
A Float64
constant equal to 60minute
. Useful for increasing the clarity of scripts, e.g. Δt = 3hour
.
Oceananigans.minute
— Constantminute
A Float64
constant equal to 60second
. Useful for increasing the clarity of scripts, e.g. Δt = 15minute
.
Oceananigans.second
— Constantsecond
A Float64
constant equal to 1.0. Useful for increasing the clarity of scripts, e.g. Δt = 1second
.
Oceananigans.pretty_filesize
— Functionpretty_filesize(s, suffix="B")
Convert a floating point value s
representing a file size to a more human-friendly formatted string with one decimal places with a suffix
defaulting to "B". Depending on the value of s
the string will be formatted to show s
using an SI prefix from bytes, kiB (1024 bytes), MiB (1024² bytes), and so on up to YiB (1024⁸ bytes).
Oceananigans.prettytime
— Methodprettytime(t)
Convert a floating point value t
representing an amount of time in seconds to a more human-friendly formatted string with three decimal places. Depending on the value of t
the string will be formatted to show t
in nanoseconds (ns), microseconds (μs), milliseconds (ms), seconds (s), minutes (min), hours (hr), or days (day).
Oceananigans.update_Δt!
— Methodupdate_Δt!(wizard, model)
Compute wizard.Δt
given the velocities and diffusivities of model
, and the parameters of wizard
.