Library

Documenting the public user interface.

Advection

Oceananigans.Advection.div_UcMethod
div_uc(i, j, k, grid, advection, U, c)

Calculates the divergence of the flux of a tracer quantity c being advected by a velocity field U = (u, v, w), ∇·(Uc),

1/V * [δxᶜᵃᵃ(Ax * u * ℑxᶠᵃᵃ(c)) + δyᵃᶜᵃ(Ay * v * ℑyᵃᶠᵃ(c)) + δzᵃᵃᶜ(Az * w * ℑzᵃᵃᶠ(c))]

which will end up at the location ccc.

source
Oceananigans.Advection.div_UuMethod
div_Uu(i, j, k, grid, advection, U, u)

Calculate the advection of momentum in the x-direction using the conservative form, ∇·(Uu)

1/Vᵘ * [δxᶠᵃᵃ(ℑxᶜᵃᵃ(Ax * u) * ℑxᶜᵃᵃ(u)) + δy_fca(ℑxᶠᵃᵃ(Ay * v) * ℑyᵃᶠᵃ(u)) + δz_fac(ℑxᶠᵃᵃ(Az * w) * ℑzᵃᵃᶠ(u))]

which will end up at the location fcc.

source
Oceananigans.Advection.div_UvMethod
div_Uv(i, j, k, grid, advection, U, v)

Calculate the advection of momentum in the y-direction using the conservative form, ∇·(Uv)

1/Vʸ * [δx_cfa(ℑyᵃᶠᵃ(Ax * u) * ℑxᶠᵃᵃ(v)) + δyᵃᶠᵃ(ℑyᵃᶜᵃ(Ay * v) * ℑyᵃᶜᵃ(v)) + δz_afc(ℑxᶠᵃᵃ(Az * w) * ℑzᵃᵃᶠ(w))]

which will end up at the location cfc.

source
Oceananigans.Advection.div_UwMethod
div_Uw(i, j, k, grid, advection, U, w)

Calculate the advection of momentum in the z-direction using the conservative form, ∇·(Uw)

1/Vʷ * [δx_caf(ℑzᵃᵃᶠ(Ax * u) * ℑxᶠᵃᵃ(w)) + δy_acf(ℑzᵃᵃᶠ(Ay * v) * ℑyᵃᶠᵃ(w)) + δzᵃᵃᶠ(ℑzᵃᵃᶜ(Az * w) * ℑzᵃᵃᶜ(w))]

which will end up at the location ccf.

source

Architectures

Oceananigans.Architectures.CPUType
CPU <: AbstractArchitecture

Run Oceananigans on one CPU node. Uses multiple threads if the environment variable JULIA_NUM_THREADS is set.

source
Oceananigans.Architectures.@hascudaMacro
@hascuda expr

A macro to compile and execute expr only if CUDA is installed and available. Generally used to wrap expressions that can only be compiled if CuArrays and CUDAnative can be loaded.

source

Boundary conditions

Oceananigans.BoundaryConditions.FluxType
Flux

A type specifying a boundary condition on the flux of a field.

The sign convention is such that a positive flux represents the flux of a quantity in the positive direction. For example, a positive vertical flux implies a quantity is fluxed upwards, in the +z direction.

Due to this convention, a positive flux applied to the top boundary specifies that a quantity is fluxed upwards across the top boundary and thus out of the domain. As a result, a positive flux applied to a top boundary leads to a reduction of that quantity in the interior of the domain; for example, a positive, upwards flux of heat at the top of the domain acts to cool the interior of the domain. Conversely, a positive flux applied to the bottom boundary leads to an increase of the quantity in the interior of the domain. The same logic holds for east, west, north, and south boundaries.

source
Oceananigans.BoundaryConditions.NormalFlowType
NormalFlow

A type specifying the component of a velocity field normal to a boundary.

Thus NormalFlow can only be applied to u along x, v along y, or w along z. For all other cases –- fields located at (Cell, Cell, Cell), or u, v, and w in (y, z), (x, z), and (x, y), respectively, either Value, Gradient, or Flux conditions must be used.

Note that NormalFlow differs from a zero Value boundary condition: Value imposes values at cell centers, while NormalFlow imposes values on a boundary, at cell faces. Only wall-normal components of the velocity field are defined on cell faces with respect to the wall-normal direction, and therefore only wall-normal components of the velocity field are defined on boundaries. Both tracers and wall-tangential components of velocity fields are defined at cell centers with respect to the wall-normal direction.

source
Oceananigans.BoundaryConditions.BoundaryConditionMethod
BoundaryCondition(BC, condition::Function; parameters=nothing, discrete_form=false)

Construct a boundary condition of type BC with a function boundary condition.

By default, the function boudnary condition is assumed to have the 'continuous form' condition(ξ, η, t), where t is time and ξ and η vary along the boundary. In particular:

  • On x-boundaries, condition(y, z, t).
  • On y-boundaries, condition(x, z, t).
  • On z-boundaries, condition(x, y, t).

If parameters is not nothing, then function boundary conditions have the form func(ξ, η, t, parameters), where ξ and η are spatial coordinates varying along the boundary as explained above.

If discrete_form=true, the function condition is assumed to have the "discrete form",

`condition(i, j, grid, clock, model_fields)`,

where i, and j are indices that vary along the boundary. If discrete_form=true and parameters is not nothing, the function condition is called with

`condition(i, j, grid, clock, model_fields, parameters)`.
source
Oceananigans.BoundaryConditions.CoordinateBoundaryConditionsType
CoordinateBoundaryConditions(left, right)

A set of two BoundaryConditions 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.

source
Oceananigans.BoundaryConditions.FieldBoundaryConditionsMethod
FieldBoundaryConditions(grid, loc;   east = DefaultBoundaryCondition(topology(grid, 1), loc[1]),
                                     west = DefaultBoundaryCondition(topology(grid, 1], loc[1]),
                                    south = DefaultBoundaryCondition(topology(grid, 2), loc[2]),
                                    north = DefaultBoundaryCondition(topology(grid, 2), loc[2]),
                                   bottom = DefaultBoundaryCondition(topology(grid, 3), loc[3]),
                                      top = DefaultBoundaryCondition(topology(grid, 3), loc[3]))

Construct FieldBoundaryConditions for a field with location loc (a 3-tuple of Face or Cell) defined on grid.

Boundary conditions on x-, y-, and z-boundaries are specified via keyword arguments:

* `west` and `east` for the `-x` and `+x` boundary;
* `south` and `north` for the `-y` and `+y` boundary;
* `bottom` and `top` for the `-z` and `+z` boundary.

Default boundary conditions depend on topology(grid) and loc.

source
Oceananigans.BoundaryConditions.TracerBoundaryConditionsMethod
TracerBoundaryConditions(grid;   east = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                 west = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                south = DefaultBoundaryCondition(topology(grid, 2), Cell),
                                north = DefaultBoundaryCondition(topology(grid, 2), Cell),
                               bottom = DefaultBoundaryCondition(topology(grid, 3), Cell),
                                  top = DefaultBoundaryCondition(topology(grid, 3), Cell))

Construct FieldBoundaryConditions for a tracer field on grid. Boundary conditions on x-, y-, and z-boundaries are specified via keyword arguments:

* `west` and `east` for the `-x` and `+x` boundary;
* `south` and `north` for the `-y` and `+y` boundary;
* `bottom` and `top` for the `-z` and `+z` boundary.

Default boundary conditions depend on topology(grid). See DefaultBoundaryCondition.

source
Oceananigans.BoundaryConditions.UVelocityBoundaryConditionsMethod
UVelocityBoundaryConditions(grid;   east = DefaultBoundaryCondition(topology(grid, 1), Face),
                                    west = DefaultBoundaryCondition(topology(grid, 1), Face),
                                   south = DefaultBoundaryCondition(topology(grid, 2), Cell),
                                   north = DefaultBoundaryCondition(topology(grid, 2), Cell),
                                  bottom = DefaultBoundaryCondition(topology(grid, 3), Cell),
                                     top = DefaultBoundaryCondition(topology(grid, 3), Cell))

Construct FieldBoundaryConditions for the u-velocity field on grid. Boundary conditions on x-, y-, and z-boundaries are specified via keyword arguments:

* `west` and `east` for the `-x` and `+x` boundary;
* `south` and `north` for the `-y` and `+y` boundary;
* `bottom` and `top` for the `-z` and `+z` boundary.

Default boundary conditions depend on topology(grid). See DefaultBoundaryCondition.

source
Oceananigans.BoundaryConditions.VVelocityBoundaryConditionsMethod
VVelocityBoundaryConditions(grid;   east = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                    west = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                   south = DefaultBoundaryCondition(topology(grid, 2), Face),
                                   north = DefaultBoundaryCondition(topology(grid, 2), Face),
                                  bottom = DefaultBoundaryCondition(topology(grid, 3), Cell),
                                     top = DefaultBoundaryCondition(topology(grid, 3), Cell))

Construct FieldBoundaryConditions for the v-velocity field on grid. Boundary conditions on x-, y-, and z-boundaries are specified via keyword arguments:

* `west` and `east` for the `-x` and `+x` boundary;
* `south` and `north` for the `-y` and `+y` boundary;
* `bottom` and `top` for the `-z` and `+z` boundary.

Default boundary conditions depend on topology(grid). See DefaultBoundaryCondition.

source
Oceananigans.BoundaryConditions.WVelocityBoundaryConditionsMethod
WVelocityBoundaryConditions(grid;   east = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                    west = DefaultBoundaryCondition(topology(grid, 1), Cell),
                                   south = DefaultBoundaryCondition(topology(grid, 2), Cell),
                                   north = DefaultBoundaryCondition(topology(grid, 2), Cell),
                                  bottom = DefaultBoundaryCondition(topology(grid, 3), Face),
                                     top = DefaultBoundaryCondition(topology(grid, 3), Face))

Construct FieldBoundaryConditions for the w-velocity field on grid. Boundary conditions on x-, y-, and z-boundaries are specified via keyword arguments:

* `west` and `east` for the `-x` and `+x` boundary;
* `south` and `north` for the `-y` and `+y` boundary;
* `bottom` and `top` for the `-z` and `+z` boundary.

Default boundary conditions depend on topology(grid). See DefaultBoundaryCondition.

source

Buoyancy

Oceananigans.Buoyancy.SeawaterBuoyancyType
SeawaterBuoyancy{FT, EOS, T, S} <: AbstractBuoyancy{EOS}

Buoyancy model for seawater. T and S are either nothing if both temperature and salinity are active, or of type FT if temperature or salinity are constant, respectively.

source
Oceananigans.Buoyancy.SeawaterBuoyancyType
SeawaterBuoyancy([FT=Float64;] gravitational_acceleration = g_Earth,
                              equation_of_state = LinearEquationOfState(FT), 
                              constant_temperature = false, constant_salinity = false)

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.

constant_temperature indicates that buoyancy depends only on salinity. For a nonlinear equation of state, constant_temperature is used as the temperature of the system. true. The same logic with the role of salinity and temperature reversed holds when constant_salinity is provided.

For a linear equation of state, the values of constant_temperature or constant_salinity are irrelevant; in this case, constant_temperature=true (and similar for constant_salinity) is valid input.

source
Oceananigans.Buoyancy.∂x_bMethod
∂x_b(i, j, k, grid, b::SeawaterBuoyancy, C)

Returns the x-derivative of buoyancy for temperature and salt-stratified water,

\[∂_x b = g ( α ∂_x Θ - β ∂_x sᴬ ) ,\]

where g is gravitational acceleration, α is the thermal expansion coefficient, β is the haline contraction coefficient, Θ is conservative temperature, and sᴬ is absolute salinity.

Note: In Oceananigans, model.tracers.T is conservative temperature and model.tracers.S is absolute salinity.

Note that ∂x_Θ, ∂x_sᴬ, α, and β are all evaluated at cell interfaces in x and cell centers in y and z.

source
Oceananigans.Buoyancy.∂y_bMethod
∂y_b(i, j, k, grid, b::SeawaterBuoyancy, C)

Returns the y-derivative of buoyancy for temperature and salt-stratified water,

\[∂_y b = g ( α ∂_y Θ - β ∂_y sᴬ ) ,\]

where g is gravitational acceleration, α is the thermal expansion coefficient, β is the haline contraction coefficient, Θ is conservative temperature, and sᴬ is absolute salinity.

Note: In Oceananigans, model.tracers.T is conservative temperature and model.tracers.S is absolute salinity.

Note that ∂y_Θ, ∂y_sᴬ, α, and β are all evaluated at cell interfaces in y and cell centers in x and z.

source
Oceananigans.Buoyancy.∂z_bMethod
∂z_b(i, j, k, grid, b::SeawaterBuoyancy, C)

Returns the vertical derivative of buoyancy for temperature and salt-stratified water,

\[∂_z b = N^2 = g ( α ∂_z Θ - β ∂_z sᴬ ) ,\]

where g is gravitational acceleration, α is the thermal expansion coefficient, β is the haline contraction coefficient, Θ is conservative temperature, and sᴬ is absolute salinity.

Note: In Oceananigans, model.tracers.T is conservative temperature and model.tracers.S is absolute salinity.

Note that ∂z_Θ, ∂z_sᴬ, α, and β are all evaluated at cell interfaces in z and cell centers in x and y.

source
Oceananigans.Buoyancy.LinearEquationOfStateType
LinearEquationOfState([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

\[ b = g (α T - β S)\]

Default constants are taken from Table 1.2 (page 33) of Vallis, "Atmospheric and Oceanic Fluid Dynamics: Fundamentals and Large-Scale Circulation" (2ed, 2017).

source

Coriolis

Oceananigans.Coriolis.FPlaneType
FPlane([FT=Float64;] f=nothing, rotation_rate=Ω_Earth, 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).

By default, rotation_rate is assumed to be Earth's.

Also called FPlane, after the "f-plane" approximation for the local effect of a planet's rotation in a planar coordinate system tangent to the planet's surface.

source
Oceananigans.Coriolis.BetaPlaneType
BetaPlane{T} <: AbstractRotation

A parameter object for meridionally increasing Coriolis parameter (f = f₀ + βy) that accounts for the variation of the locally vertical component of the rotation vector with latitude.

source
Oceananigans.Coriolis.BetaPlaneType
BetaPlane([T=Float64;] f₀=nothing, β=nothing,
                       rotation_rate=Ω_Earth, latitude=nothing, radius=R_Earth)

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 (where y = 0) at which the β-plane approximation is to be made.

By default, the rotation_rate and planet radius is assumed to be Earth's.

source
Oceananigans.Coriolis.NonTraditionalFPlaneType
NonTraditionalFPlane{FT} <: AbstractRotation

A Coriolis implementation that accounts for both the locally vertical and the locally horizontal components of the rotation vector. Traditionally (see FPlane) only the locally vertical component is accounted for.

source
Oceananigans.Coriolis.NonTraditionalFPlaneType
NonTraditionalFPlane([FT=Float64;] fz=nothing, fy=nothing,
                                   rotation_rate=Ω_Earth, latitude=nothing)

Returns a parameter object for constant rotation about an axis in the y-z plane with y- and z-components fy/2 and fz/2, and the background vorticity is (0, fy, fz).

In oceanography fz and fy represent the components of planetary voriticity which are perpendicular and parallel to the ocean surface in a domain in which x, y, z correspond to the directions east, north, and up.

If fz and fy are not specified, they are calculated from rotation_rate and latitude according to the relations fz = 2*rotation_rate*sind(latitude) and fy = 2*rotation_rate*cosd(latitude), respectively. By default, rotation_rate is assumed to be Earth's.

source
Oceananigans.Coriolis.NonTraditionalBetaPlaneType
NonTraditionalBetaPlane{FT} <: AbstractRotation

A Coriolis implementation that accounts for the latitudinal variation of both the locally vertical and the locally horizontal components of the rotation vector. The "traditional" approximation in ocean models accounts for only the locally vertical component of the rotation vector (see BetaPlane).

This implementation is based off of section 5 of Dellar (2011). It conserve energy, angular momentum, and potential vorticity.

References

Dellar, P. (2011). Variations on a beta-plane: Derivation of non-traditional beta-plane equations from Hamilton's principle on a sphere. Journal of Fluid Mechanics, 674, 174-195. doi:10.1017/S0022112010006464

source
Oceananigans.Coriolis.NonTraditionalBetaPlaneType
NonTraditionalBetaPlane(FT=Float64;
    fz=nothing, fy=nothing, β=nothing, γ=nothing,
    rotation_rate=Ω_Earth, latitude=nothing, radius=R_Earth)

The user may directly specify fz, fy, β, γ, and radius or the three parameters rotation_rate, latitude, and radius that specify the rotation rate and radius of a planet, and the central latitude (where y = 0) at which the non-traditional β-plane approximation is to be made.

By default, the rotation_rate and planet radius is assumed to be Earth's.

source

Diagnostics

Oceananigans.Diagnostics.CFLMethod
CFL(Δt [, timescale=Oceananigans.cell_advection_timescale])

Returns an object for computing the Courant-Freidrichs-Lewy (CFL) number associated with time step or TimeStepWizard Δt and timescale.

See also AdvectiveCFL and DiffusiveCFL.

source
Oceananigans.Diagnostics.AdvectiveCFLMethod
AdvectiveCFL(Δt)

Returns an object for computing the Courant-Freidrichs-Lewy (CFL) number associated with time step or TimeStepWizard Δt and the time scale for advection across a cell.

Example

julia> model = IncompressibleModel(grid=RegularCartesianGrid(size=(16, 16, 16), length=(8, 8, 8)));

julia> cfl = AdvectiveCFL(1.0);

julia> data(model.velocities.u) .= π;

julia> cfl(model)
6.283185307179586
source
Oceananigans.Diagnostics.DiffusiveCFLMethod
DiffusiveCFL(Δt)

Returns an object for computing the diffusive Courant-Freidrichs-Lewy (CFL) number associated with time step or TimeStepWizard Δt and the time scale for diffusion across a cell associated with model.closure.

The maximum diffusive CFL number among viscosity and all tracer diffusivities is returned.

Example

julia> model = IncompressibleModel(grid=RegularCartesianGrid(size=(16, 16, 16), length=(1, 1, 1)));

julia> dcfl = DiffusiveCFL(0.1);

julia> dcfl(model)
2.688e-5
source
Oceananigans.Diagnostics.FieldMaximumType
FieldMaximum(mapping, field)

An object for calculating the maximum of a mapping function applied element-wise to field.

Examples

julia> model = IncompressibleModel(grid=RegularCartesianGrid(size=(16, 16, 16), length=(1, 1, 1)));

julia> max_abs_u = FieldMaximum(abs, model.velocities.u);

julia> max_w² = FieldMaximum(x->x^2, model.velocities.w);
source

Fields

Oceananigans.Fields.FieldType
Field{X, Y, Z, A, G, B} <: AbstractField{X, Y, Z, A, G}

A field defined at the location (X, Y, Z), each of which can be either Cell or Face, and with data stored in a container of type A (typically an array). The field is defined on a grid G and has field boundary conditions B.

source
Oceananigans.Fields.FieldType
Field(X, Y, Z, arch, grid, [  bcs = FieldBoundaryConditions(grid, (X, Y, Z)),
                             data = new_data(arch, grid, (X, Y, Z)) ] )

Construct a Field on grid with data on architecture arch with boundary conditions bcs. Each of (X, Y, Z) is either Cell or Face and determines the field's location in (x, y, z).

Example

julia> ω = Field(Face, Face, Cell, CPU(), RegularCartesianmodel.grid)

source
Oceananigans.Fields.FieldMethod
Field(L::Tuple, arch, grid, data, bcs)

Construct a Field at the location defined by the 3-tuple L, whose elements are Cell or Face.

source
Oceananigans.Fields.CellFieldFunction
CellField([ FT=eltype(grid) ], arch::AbstractArchitecture, grid,
          [  bcs = TracerBoundaryConditions(grid),
            data = new_data(FT, arch, grid, (Cell, Cell, Cell) ] )

Return a Field{Cell, Cell, Cell} on architecture arch and grid containing data with field boundary conditions bcs.

source
Oceananigans.Fields.XFaceFieldFunction
XFaceField([ FT=eltype(grid) ], arch::AbstractArchitecture, grid,
           [  bcs = UVelocityBoundaryConditions(grid),
             data = new_data(FT, arch, grid, (Face, Cell, Cell) ] )

Return a Field{Face, Cell, Cell} on architecture arch and grid containing data with field boundary conditions bcs.

source
Oceananigans.Fields.YFaceFieldFunction
YFaceField([ FT=eltype(grid) ], arch::AbstractArchitecture, grid,
           [  bcs = VVelocityBoundaryConditions(grid),
             data = new_data(FT, arch, grid, (Cell, Face, Cell)) ] )

Return a Field{Cell, Face, Cell} on architecture arch and grid containing data with field boundary conditions bcs.

source
Oceananigans.Fields.ZFaceFieldFunction
ZFaceField([ FT=eltype(grid) ], arch::AbstractArchitecture, grid,
           [  bcs = WVelocityBoundaryConditions(grid),
             data = new_data(FT, arch, grid, (Cell, Cell, Face)) ] )

Return a Field{Cell, Cell, Face} on architecture arch and grid containing data with field boundary conditions bcs.

source
Oceananigans.Fields.set!Method
set!(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 = IncompressibleModel(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₀)
source

Forcings

Oceananigans.Forcings.ContinuousForcingMethod
ContinuousForcing(func; parameters=nothing, field_dependencies=())

Construct a "continuous form" forcing with optional parameters and optional field_dependencies on other fields in a model.

If neither parameters nor field_dependencies are provided, then func must be callable with the signature

`func(x, y, z, t)`

where x, y, z are the east-west, north-south, and vertical spatial coordinates, and t is time.

If field_dependencies are provided, the signature of func must include them. For example, if field_dependencies=(:u, :S) (and parameters are not provided), then func must be callable with the signature

`func(x, y, z, t, u, S)`

where u is assumed to be the u-velocity component, and S is a tracer. Note that any field which does not have the name u, v, or w is assumed to be a tracer and must be present in model.tracers.

If parameters are provided, then the last argument to func must be parameters. For example, if func has no field_dependencies but does depend on parameters, then it must be callable with the signature

`func(x, y, z, t, parameters)`

With field_dependencies=(:u, :v, :w, :c) and parameters, then func must be callable with the signature

`func(x, y, z, t, u, v, w, c, parameters)`
source
Oceananigans.Forcings.DiscreteForcingMethod
DiscreteForcing(func; parameters=nothing)

Construct a "discrete form" forcing function with optional parameters. The forcing function is applied at grid point i, j, k.

When parameters are not specified, func must be callable with the signature

`func(i, j, k, grid, clock, model_fields)`

where grid is model.grid, clock.time is the current simulation time and clock.iteration is the current model iteration, and model_fields is a NamedTuple with u, v, w and the fields in model.tracers.

Note that the index end does not access the final physical grid point of a model field in any direction. The final grid point must be explicitly specified, as in model_fields.u[i, j, grid.Nz].

When parameters is specified, func must be callable with the signature.

`func(i, j, k, grid, clock, model_fields, parameters)`

parameters is arbitrary in principle, however GPU compilation can place constraints on typeof(parameters).

source
Oceananigans.Forcings.ForcingMethod
Forcing(func; parameters=nothing, field_dependencies=(), discrete_form=false)

Returns a forcing function added to the tendency of an Oceananigans model field.

If discrete_form=false (the default), and neither parameters nor field_dependencies are provided, then func must be callable with the signature

`func(x, y, z, t)`

where x, y, z are the east-west, north-south, and vertical spatial coordinates, and t is time. Note that this form is also default in the constructor for IncompressibleModel, so that Forcing is not needed.

If discrete_form=false (the default), and field_dependencies are provided, the signature of func must include them. For example, if field_dependencies=(:u, :S) (and parameters are not provided), then func must be callable with the signature

`func(x, y, z, t, u, S)`

where u is assumed to be the u-velocity component, and S is a tracer. Note that any field which does not have the name u, v, or w is assumed to be a tracer and must be present in model.tracers.

If discrete_form=false (the default) and parameters are provided, then the last argument to func must be parameters. For example, if func has no field_dependencies but does depend on parameters, then it must be callable with the signature

`func(x, y, z, t, parameters)`

The object parameters is arbitrary in principle, however GPU compilation can place constraints on typeof(parameters).

With field_dependencies=(:u, :v, :w, :c) and parameters, then func must be callable with the signature

`func(x, y, z, t, u, v, w, c, parameters)`

If discrete_form=true then func must be callable with the "discrete form"

`func(i, j, k, grid, clock, model_fields)`

where i, j, k is the grid point at which the forcing is applied, grid is model.grid, clock.time is the current simulation time and clock.iteration is the current model iteration, and model_fields is a NamedTuple with u, v, w, the fields in model.tracers, and the fields in model.diffusivities, each of which is an OffsetArrays (or NamedTuples of OffsetArrays depending on the turbulence closure) of field data.

When discrete_form=true and parameters is specified, func must be callable with the signature

`func(i, j, k, grid, clock, model_fields, parameters)`

Examples

using Oceananigans

# Parameterized forcing
parameterized_func(x, y, z, t, p) = p.μ * exp(z / p.λ) * cos(p.ω * t)

v_forcing = Forcing(parameterized_func, parameters = (μ=42, λ=0.1, ω=π))

# output
ContinuousForcing{NamedTuple{(:μ, :λ, :ω),Tuple{Int64,Float64,Irrational{:π}}}}
├── func: parameterized_func
├── parameters: (μ = 42, λ = 0.1, ω = π)
└── field dependencies: ()

Note that because forcing locations are regularized within the IncompressibleModel constructor:

grid = RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1))
model = IncompressibleModel(grid=grid, forcing=(v=v_forcing,))

model.forcing.v

# output
ContinuousForcing{NamedTuple{(:μ, :λ, :ω),Tuple{Int64,Float64,Irrational{:π}}}} at (Cell, Face, Cell)
├── func: parameterized_func
├── parameters: (μ = 42, λ = 0.1, ω = π)
└── field dependencies: ()

After passing through the constructor for IncompressibleModel, the v-forcing location information is available and set to Cell, Face, Cell.

# Field-dependent forcing
growth_in_sunlight(x, y, z, t, P) = exp(z) * P

plankton_forcing = Forcing(growth_in_sunlight, field_dependencies=:P)

# output
ContinuousForcing{Nothing}
├── func: growth_in_sunlight
├── parameters: nothing
└── field dependencies: (:P,)
# Parameterized, field-dependent forcing
tracer_relaxation(x, y, z, t, c, p) = p.μ * exp((z + p.H) / p.λ) * (p.dCdz * z - c) 

c_forcing = Forcing(tracer_relaxation,
                    field_dependencies = :c,
                            parameters = (μ=1/60, λ=10, H=1000, dCdz=1))

# output
ContinuousForcing{NamedTuple{(:μ, :λ, :H, :dCdz),Tuple{Float64,Int64,Int64,Int64}}}
├── func: tracer_relaxation
├── parameters: (μ = 0.016666666666666666, λ = 10, H = 1000, dCdz = 1)
└── field dependencies: (:c,)
# Unparameterized discrete-form forcing function
filtered_relaxation(i, j, k, grid, clock, model_fields) =
    @inbounds - (model_fields.c[i-1, j, k] + model_fields.c[i, j, k] + model_fields.c[i+1, j, k]) / 3

filtered_forcing = Forcing(filtered_relaxation, discrete_form=true)

# output
DiscreteForcing{Nothing}
├── func: filtered_relaxation
└── parameters: nothing
# Discrete-form forcing function with parameters
masked_damping(i, j, k, grid, clock, model_fields, parameters) = 
    @inbounds - parameters.μ * exp(grid.zC[k] / parameters.λ) * model_fields.u[i, j, k]

masked_damping_forcing = Forcing(masked_damping, parameters=(μ=42, λ=π), discrete_form=true)

# output
DiscreteForcing{NamedTuple{(:μ, :λ),Tuple{Int64,Irrational{:π}}}}
├── func: masked_damping
└── parameters: (μ = 42, λ = π)
source
Oceananigans.Forcings.GaussianMaskType
GaussianMask{D}(center, width)

Callable object that returns a Gaussian masking function centered on center, with width, and varying along direction D.

Examples

  • Create a Gaussian mask centered on z=0 with width 1 meter.
julia> mask = GaussianMask{:z}(center=0, width=1)
source
Oceananigans.Forcings.LinearTargetType
LinearTarget{D}(intercept, gradient)

Callable object that returns a Linear target function with intercept and gradient, and varying along direction D.

Examples

  • Create a linear target function varying in z, equal to 0 at z=0 and with gradient 10⁻⁶:
julia> target = LinearTarget{:z}(intercept=0, gradient=1e-6)
source
Oceananigans.Forcings.RelaxationMethod
Relaxation(; rate, mask=onefunction, target=zerofunction)

Returns a Forcing that restores a field to target(x, y, z, t) at the specified rate, in the region mask(x, y, z).

The functions onefunction and zerofunction always return 1 and 0, respectively. Thus the default mask leaves the whole domain uncovered, and the default target is zero.

Example

  • Restore a field to zero on a timescale of "3600" (equal to one hour if the time units of the simulation are seconds).
using Oceananigans

damping = Relaxation(rate = 1/3600)

# output
Relaxation{Float64, typeof(Oceananigans.Forcings.onefunction), typeof(Oceananigans.Forcings.zerofunction)}
├── rate: 0.0002777777777777778
├── mask: 1
└── target: 0
  • Restore a field to a linear z-gradient within the bottom 1/4 of a domain on a timescale of "60" (equal to one minute if the time units of the simulation are seconds).
dTdz = 0.001 # ⁰C m⁻¹, temperature gradient

T₀ = 20 # ⁰C, surface temperature at z=0

Lz = 100 # m, depth of domain

bottom_sponge_layer = Relaxation(; rate = 1/60,
                                   target = LinearTarget{:z}(intercept=T₀, gradient=dTdz),
                                   mask = GaussianMask{:z}(center=-Lz, width=Lz/4))

# output
Relaxation{Float64, GaussianMask{:z,Float64}, LinearTarget{:z,Float64}}
├── rate: 0.016666666666666666
├── mask: exp(-(z + 100.0)^2 / (2 * 25.0^2))
└── target: 20.0 + 0.001 * z
source

Grids

Oceananigans.Grids.FlatType
Flat

Grid topology for flat dimensions, generally with one grid point, along which the solution is uniform and does not vary.

source
Oceananigans.Grids.nodesMethod
nodes(loc, grid; reshape=false)

Returns a 3-tuple of views over the interior nodes at the locations in loc in x, y, z.

If reshape=true, the views are reshaped to 3D arrays with non-singleton dimensions 1, 2, 3 for x, y, z, respectively. These reshaped arrays can then be used in broadcast operations with 3D fields or arrays.

See xnodes, ynodes, and znodes.

source
Oceananigans.Grids.xnodesMethod
xnodes(loc, grid, reshape=false)

Returns a view over the interior loc=Cell or loc=Face nodes on grid in the x-direction. For Bounded directions, Face nodes include the boundary points. reshape=false will return a 1D array while reshape=true will return a 3D array with size Nx×1×1.

See znodes for examples.

source
Oceananigans.Grids.ynodesMethod
ynodes(loc, grid, reshape=false)

Returns a view over the interior loc=Cell or loc=Face nodes on grid in the y-direction. For Bounded directions, Face nodes include the boundary points. reshape=false will return a 1D array while reshape=true will return a 3D array with size 1×Ny×1.

See znodes for examples.

source
Oceananigans.Grids.znodesMethod
znodes(loc, grid, reshape=false)

Returns a view over the interior loc=Cell or loc=Face nodes on grid in the z-direction. For Bounded directions, Face nodes include the boundary points. reshape=false will return a 1D array while reshape=true will return a 3D array with size 1×1×Nz.

Examples

julia> using Oceananigans, Oceananigans.Grids

julia> horz_periodic_grid = RegularCartesianGrid(size=(3, 3, 3), extent=(2π, 2π, 1),
                                                 topology=(Periodic, Periodic, Bounded));

julia> zC = znodes(Cell, horz_periodic_grid)
3-element view(OffsetArray(::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}, 0:4), 1:3) with eltype Float64:
 -0.8333333333333331
 -0.4999999999999999
 -0.16666666666666652
julia> zF = znodes(Face, horz_periodic_grid)
4-element view(OffsetArray(::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}, 0:5), 1:4) with eltype Float64:
 -1.0
 -0.6666666666666666
 -0.33333333333333337
 -4.44089209850063e-17
source
Oceananigans.Grids.RegularCartesianGridType
RegularCartesianGrid{FT, TX, TY, TZ, R} <: AbstractGrid{FT, TX, TY, TZ}

A Cartesian grid with with constant grid spacings Δx, Δy, and Δz between cell centers and cell faces, elements of type FT, topology {TX, TY, TZ}, and coordinate ranges of type R.

source
Oceananigans.Grids.RegularCartesianGridType
RegularCartesianGrid([FT=Float64]; size,
                     extent = nothing, x = nothing, y = nothing, z = nothing,
                     topology = (Periodic, Periodic, Bounded), halo = (1, 1, 1))

Creates a RegularCartesianGrid with size = (Nx, Ny, Nz) grid points.

Keyword arguments

  • size (required): A tuple prescribing the number of grid points in non-Flat directions. size is a 3-tuple for 3D models, a 2-tuple for 2D models, and either a scalar or 1-tuple for 1D models.

  • topology: A 3-tuple (Tx, Ty, Tz) specifying the topology of the domain. Tx, Ty, and Tz specify whether the x-, y-, and z directions are Periodic, Bounded, or Flat. The topology Flat indicates that a model does not vary in that directions so that derivatives and interpolation are zero. The default is topology=(Periodic, Periodic, Bounded).

  • extent: A tuple prescribing the physical extent of the grid in non-Flat directions. The origin for three-dimensional domains is the oceanic default (0, 0, -Lz).

  • x, y, and z: Each of x, y, z are 2-tuples that specify the end points of the domain in their respect directions. Scalar values may be used in Flat directions.

Note: Either extent, or all of x, y, and z must be specified.

  • halo: A tuple of integers that specifies the size of the halo region of cells surrounding the physical interior for each non-Flat direction.

The physical extent of the domain can be specified via x, y, and z keyword arguments indicating the left and right endpoints of each dimensions, e.g. x=(-π, π) or via the extent argument, e.g. extent=(Lx, Ly, Lz) which specifies the extent of each dimension in which case 0 ≤ x ≤ Lx, 0 ≤ y ≤ Ly, and -Lz ≤ z ≤ 0.

A grid topology may be specified via a tuple assigning one of Periodic, Bounded, andFlatto each dimension. By default, a horizontally periodic grid topology(Periodic, Periodic, Bounded)` is assumed.

Constants are stored using floating point values of type FT. By default this is Float64. Make sure to specify the desired FT if not using Float64.

Grid properties

  • (Nx, Ny, Nz)::Int: Number of physical points in the (x, y, z)-direction

  • (Hx, Hy, Hz)::Int: Number of halo points in the (x, y, z)-direction

  • (Lx, Ly, Lz)::FT: Physical extent of the grid in the (x, y, z)-direction

  • (Δx, Δy, Δz)::FT: Cell width in the (x, y, z)-direction

  • (xC, yC, zC): (x, y, z) coordinates of cell centers.

  • (xF, yF, zF): (x, y, z) coordinates of cell faces.

Examples

  • A default grid with Float64 type:
julia> using Oceananigans

julia> grid = RegularCartesianGrid(size=(32, 32, 32), extent=(1, 2, 3))
RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}
                   domain: x ∈ [0.0, 1.0], y ∈ [0.0, 2.0], z ∈ [-3.0, 0.0]
                 topology: (Periodic, Periodic, Bounded)
  resolution (Nx, Ny, Nz): (32, 32, 32)
   halo size (Hx, Hy, Hz): (1, 1, 1)
grid spacing (Δx, Δy, Δz): (0.03125, 0.0625, 0.09375)
  • A default grid with Float32 type:
julia> using Oceananigans

julia> grid = RegularCartesianGrid(Float32; size=(32, 32, 16), x=(0, 8), y=(-10, 10), z=(-π, π))
RegularCartesianGrid{Float32, Periodic, Periodic, Bounded}
                   domain: x ∈ [0.0, 8.0], y ∈ [-10.0, 10.0], z ∈ [-3.1415927, 3.1415927]
                 topology: (Periodic, Periodic, Bounded)
  resolution (Nx, Ny, Nz): (32, 32, 16)
   halo size (Hx, Hy, Hz): (1, 1, 1)
grid spacing (Δx, Δy, Δz): (0.25f0, 0.625f0, 0.3926991f0)
  • A two-dimenisional, horizontally-periodic grid:
julia> using Oceananigans

julia> grid = RegularCartesianGrid(size=(32, 32), extent=(2π, 4π), topology=(Periodic, Periodic, Flat))
RegularCartesianGrid{Float64, Periodic, Periodic, Flat}
                   domain: x ∈ [0.0, 6.283185307179586], y ∈ [0.0, 12.566370614359172], z ∈ [0.0, 0.0]
                 topology: (Periodic, Periodic, Flat)
  resolution (Nx, Ny, Nz): (32, 32, 1)
   halo size (Hx, Hy, Hz): (1, 1, 0)
grid spacing (Δx, Δy, Δz): (0.19634954084936207, 0.39269908169872414, 0.0)
  • A one-dimensional "column" grid:
julia> using Oceananigans

julia> grid = RegularCartesianGrid(size=256, z=(-128, 0), topology=(Flat, Flat, Bounded))
RegularCartesianGrid{Float64, Flat, Flat, Bounded}
                   domain: x ∈ [0.0, 0.0], y ∈ [0.0, 0.0], z ∈ [-128.0, 0.0]
                 topology: (Flat, Flat, Bounded)
  resolution (Nx, Ny, Nz): (1, 1, 256)
   halo size (Hx, Hy, Hz): (0, 0, 1)
grid spacing (Δx, Δy, Δz): (0.0, 0.0, 0.5)
source

Logger

Oceananigans.Logger.OceananigansLoggerType
OceananigansLogger(stream::IO=stdout, level=Logging.Info; show_info_source=false)

Based on Logging.SimpleLogger, it tries to log all messages in the following format:

[yyyy/mm/dd HH:MM:SS.sss] log_level message [-@-> source_file:line_number]

where the source of the message between the square brackets is included only if show_info_source=true or if the message is not an info level message.

source

Models

Oceananigans.Models.ClockType
Clock{T<:Number}

Keeps track of the current time, iteration number, and time-stepping stage. stage is updated only for multi-stage time-stepping methods. The time::T can be either a number of a DateTime object.

source
Oceananigans.Models.ClockMethod
Clock(; time, iteration=0, stage=1)

Returns a Clock initialized to the zeroth iteration and first time step stage.

source
Oceananigans.Models.IncompressibleModelMethod
IncompressibleModel(;
               grid,
       architecture = CPU(),
         float_type = Float64,
              clock = Clock{float_type}(0, 0, 1),
          advection = CenteredSecondOrder(),
           buoyancy = SeawaterBuoyancy(float_type),
           coriolis = nothing,
      surface_waves = nothing,
            forcing = NamedTuple(),
            closure = IsotropicDiffusivity(float_type, ν=ν₀, κ=κ₀),
boundary_conditions = NamedTuple(),
            tracers = (:T, :S),
        timestepper = :QuasiAdamsBashforth2,
  background_fields = NamedTuple(),
         velocities = nothing,
          pressures = nothing,
      diffusivities = nothing,
    pressure_solver = nothing
)

Construct an incompressible Oceananigans.jl model on grid.

Keyword arguments

- `grid`: (required) The resolution and discrete geometry on which `model` is solved.
- `architecture`: `CPU()` or `GPU()`. The computer architecture used to time-step `model`.
- `float_type`: `Float32` or `Float64`. The floating point type used for `model` data.
- `advection`: The scheme that advects velocities and tracers. See `Oceananigans.Advection`.
- `buoyancy`: The buoyancy model. See `Oceananigans.Buoyancy`.
- `closure`: The turbulence closure for `model`. See `Oceananigans.TurbulenceClosures`.
- `coriolis`: Parameters for the background rotation rate of the model.
- `forcing`: `NamedTuple` of user-defined forcing functions that contribute to solution tendencies.
- `boundary_conditions`: `NamedTuple` containing field boundary conditions.
- `tracers`: A tuple of symbols defining the names of the modeled tracers, or a `NamedTuple` of
             preallocated `CellField`s.
- `timestepper`: A symbol that specifies the time-stepping method. Either `:QuasiAdamsBashforth2` or
                 `:RungeKutta3`.
source
Oceananigans.Models.NonDimensionalModelMethod
NonDimensionalModel(; 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 IncompressibleModel constructor.

source

Output writers

Oceananigans.OutputWriters.JLD2OutputWriterMethod
JLD2OutputWriter(model, outputs; prefix, schedule,
                          dir = ".",
                 field_slicer = FieldSlicer(),
                   array_type = Array{Float32},
                 max_filesize = Inf,
                        force = false,
                         init = noinit,
                    including = [:grid, :coriolis, :buoyancy, :closure],
                      verbose = false,
                         part = 1,
                      jld2_kw = Dict{Symbol, Any}())

Construct a JLD2OutputWriter for an Oceananigans model that writes label, output pairs in outputs to a JLD2 file.

The argument outputs may be a Dict or NamedTuple. The keys of outputs are symbols or strings that "name" output data. The values of outputs are either AbstractFields, objects that are called with the signature output(model), or WindowedTimeAverages of AbstractFieldss, functions, or callable objects.

Keyword arguments

## Filenaming

- `prefix` (required): Descriptive filename prefixed to all output files.

- `dir`: Directory to save output to.
         Default: "." (current working directory).

## Output frequency and time-averaging

- `schedule` (required): `AbstractSchedule` that determines when output is saved.

## Slicing and type conversion prior to output

- `field_slicer`: An object for slicing field output in ``(x, y, z)``, including omitting halos.
                  Has no effect on output that is not a field. `field_slicer = nothing` means
                  no slicing occurs, so that all field data, including halo regions, is saved.
                  Default: FieldSlicer(), which slices halo regions.

- `array_type`: The array type to which output arrays are converted to prior to saving.
                Default: Array{Float32}.

## File management

- `max_filesize`: The writer will stop writing to the output file once the file size exceeds `max_filesize`,
                  and write to a new one with a consistent naming scheme ending in `part1`, `part2`, etc.
                  Defaults to `Inf`.

- `force`: Remove existing files if their filenames conflict.
           Default: `false`.

## Output file metadata management

- `init`: A function of the form `init(file, model)` that runs when a JLD2 output file is initialized.
          Default: `noinit(args...) = nothing`.

- `including`: List of model properties to save with every file.
               Default: `[:grid, :coriolis, :buoyancy, :closure]`

## Miscellaneous keywords

- `verbose`: Log what the output writer is doing with statistics on compute/write times and file sizes.
             Default: `false`.

- `part`: The starting part number used if `max_filesize` is finite.
          Default: 1.

- `jld2_kw`: Dict of kwargs to be passed to `jldopen` when data is written.

Example

Write out 3D fields for w and T and a horizontal average:

using Oceananigans, Oceananigans.OutputWriters, Oceananigans.Fields
using Oceananigans.Utils: hour, minute

model = IncompressibleModel(grid=RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1)))
simulation = Simulation(model, Δt=12, stop_time=1hour)

function init_save_some_metadata!(file, model)
    file["author"] = "Chim Riggles"
    file["parameters/coriolis_parameter"] = 1e-4
    file["parameters/density"] = 1027
    return nothing
end

T_avg =  AveragedField(model.tracers.T, dims=(1, 2))

# Note that model.velocities is NamedTuple
simulation.output_writers[:velocities] = JLD2OutputWriter(model, model.velocities,
                                                          prefix = "some_data",
                                                          schedule = TimeInterval(20minute),
                                                          init = init_save_some_metadata!)

# output
JLD2OutputWriter scheduled on TimeInterval(20 minutes):
├── filepath: ./some_data.jld2
├── 3 outputs: (:u, :v, :w)
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB

and a time- and horizontal-average of temperature T every 1 hour of simulation time to a file called some_averaged_data.jld2

simulation.output_writers[:avg_T] = JLD2OutputWriter(model, (T=T_avg,),
                                                     prefix = "some_averaged_data",
                                                     schedule = AveragedTimeInterval(20minute, window=5minute))

# output
JLD2OutputWriter scheduled on TimeInterval(20 minutes):
├── filepath: ./some_averaged_data.jld2
├── 1 outputs: (:T,) averaged on AveragedTimeInterval(window=5 minutes, stride=1, interval=20 minutes)
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB
source
Oceananigans.OutputWriters.NetCDFOutputWriterMethod

function NetCDFOutputWriter(model, outputs; filepath, schedule arraytype = Array{Float32}, fieldslicer = FieldSlicer(), globalattributes = Dict(), outputattributes = Dict(), dimensions = Dict(), mode = "c", compression = 0, verbose = false)

Construct a NetCDFOutputWriter that writes (label, output) pairs in outputs (which should be a Dict) to a NetCDF file, where label is a string that labels the output and output is either a Field (e.g. model.velocities.u or an AveragedField) or a function f(model) that returns something to be written to disk. Custom output requires the spatial dimensions (a Dict) to be manually specified (see examples).

Keyword arguments

  • filepath (required): Filepath to save output to.

  • schedule (required): AbstractSchedule that determines when output is saved.

  • array_type: The array type to which output arrays are converted to prior to saving. Default: Array{Float32}.

  • field_slicer: An object for slicing field output in $(x, y, z)$, including omitting halos. Has no effect on output that is not a field. field_slicer = nothing means no slicing occurs, so that all field data, including halo regions, is saved. Default: FieldSlicer(), which slices halo regions.

  • global_attributes: Dict of model properties to save with every file (deafult: Dict())

  • output_attributes: Dict of attributes to be saved with each field variable (reasonable defaults are provided for velocities, buoyancy, temperature, and salinity; otherwise output_attributes must be user-provided).

  • dimensions: A Dict of dimension tuples to apply to outputs (required for function outputs)

  • with_halos: Include the halo regions in the grid coordinates and output fields (default: false).

  • mode: "a" (for append) and "c" (for clobber or create). Default: "c". See NCDatasets.jl documentation for more information on the mode option.

  • compression: Determines the compression level of data (0-9, default 0)

Examples

Saving the u velocity field and temperature fields, the full 3D fields and surface 2D slices to separate NetCDF files:

using Oceananigans, Oceananigans.OutputWriters

grid = RegularCartesianGrid(size=(16, 16, 16), extent=(1, 1, 1));

model = IncompressibleModel(grid=grid);

simulation = Simulation(model, Δt=12, stop_time=3600);

fields = Dict("u" => model.velocities.u, "T" => model.tracers.T);

simulation.output_writers[:field_writer] =
    NetCDFOutputWriter(model, fields, filepath="fields.nc", schedule=TimeInterval(60))

# output
NetCDFOutputWriter scheduled on TimeInterval(1 minute):
├── filepath: fields.nc
├── dimensions: zC(16), zF(17), xC(16), yF(16), xF(16), yC(16), time(0)
├── 2 outputs: ["T", "u"]
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
└── array type: Array{Float32}
simulation.output_writers[:surface_slice_writer] =
    NetCDFOutputWriter(model, fields, filepath="surface_xy_slice.nc",
                       schedule=TimeInterval(60), field_slicer=FieldSlicer(k=grid.Nz))

# output
NetCDFOutputWriter scheduled on TimeInterval(1 minute):
├── filepath: surface_xy_slice.nc
├── dimensions: zC(1), zF(1), xC(16), yF(16), xF(16), yC(16), time(0)
├── 2 outputs: ["T", "u"]
├── field slicer: FieldSlicer(:, :, 16, with_halos=false)
└── array type: Array{Float32}
simulation.output_writers[:averaged_profile_writer] =
    NetCDFOutputWriter(model, fields,
                       filepath = "averaged_z_profile.nc",
                       schedule = AveragedTimeInterval(60, window=20),
                       field_slicer = FieldSlicer(i=1, j=1))

# output
NetCDFOutputWriter scheduled on TimeInterval(1 minute):
├── filepath: averaged_z_profile.nc
├── dimensions: zC(16), zF(17), xC(1), yF(1), xF(1), yC(1), time(0)
├── 2 outputs: ["T", "u"] averaged on AveragedTimeInterval(window=20 seconds, stride=1, interval=1 minute)
├── field slicer: FieldSlicer(1, 1, :, with_halos=false)
└── array type: Array{Float32}

NetCDFOutputWriter also accepts output functions that write scalars and arrays to disk, provided that their dimensions are provided:

using Oceananigans, Oceananigans.OutputWriters

grid = RegularCartesianGrid(size=(16, 16, 16), extent=(1, 2, 3));

model = IncompressibleModel(grid=grid);

simulation = Simulation(model, Δt=1.25, stop_iteration=3);

f(model) = model.clock.time^2; # scalar output

g(model) = model.clock.time .* exp.(znodes(Cell, grid)); # vector/profile output

h(model) = model.clock.time .* (   sin.(xnodes(Cell, grid, reshape=true)[:, :, 1])
                            .*     cos.(ynodes(Face, grid, reshape=true)[:, :, 1])); # xy slice output

outputs = Dict("scalar" => f, "profile" => g, "slice" => h);

dims = Dict("scalar" => (), "profile" => ("zC",), "slice" => ("xC", "yC"));

output_attributes = Dict(
    "scalar"  => Dict("longname" => "Some scalar", "units" => "bananas"),
    "profile" => Dict("longname" => "Some vertical profile", "units" => "watermelons"),
    "slice"   => Dict("longname" => "Some slice", "units" => "mushrooms")
);

global_attributes = Dict("location" => "Bay of Fundy", "onions" => 7);

simulation.output_writers[:things] =
    NetCDFOutputWriter(model, outputs,
                       schedule=IterationInterval(1), filepath="things.nc", dimensions=dims, verbose=true,
                       global_attributes=global_attributes, output_attributes=output_attributes)

# output
NetCDFOutputWriter scheduled on IterationInterval(1):
├── filepath: things.nc
├── dimensions: zC(16), zF(17), xC(16), yF(16), xF(16), yC(16), time(0)
├── 3 outputs: ["profile", "slice", "scalar"]
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
└── array type: Array{Float32}
source
Oceananigans.OutputWriters.write_output!Method
write_output!(output_writer, model)

Writes output to netcdf file output_writer.filepath at specified intervals. Increments the time dimension every time an output is written to the file.

source
Oceananigans.OutputWriters.AveragedTimeIntervalMethod
AveragedTimeInterval(interval; window=interval, stride=1)

Returns a schedule that specifies periodic time-averaging of output. The time window specifies the extent of the time-average, which reoccurs every interval.

output is computed and accumulated into the average every stride iterations during the averaging window. For example, stride=1 computs output every iteration, whereas stride=2 computes output every other iteration. Time-averages with longer strides are faster to compute, but less accurate.

The time-average of $a$ is a left Riemann sum corresponding to

$⟨a⟩ = 1/T \int_{tᵢ-T}^T a \mathrm{d} t ,$

where $⟨a⟩$ is the time-average of $a$, $T$ is the time-window for averaging, and the $tᵢ$ are discrete times separated by the time interval. The $tᵢ$ specify both the end of the averaging window and the time at which output is written.

Example

using Oceananigans.OutputWriters: AveragedTimeInterval
using Oceananigans.Utils: year, years

schedule = AveragedTimeInterval(4years, window=1year)

# output
AveragedTimeInterval(window=1 year, stride=1, interval=4 years)

An AveragedTimeInterval schedule directs an output writer to time-average its outputs before writing them to disk:

using Oceananigans
using Oceananigans.OutputWriters: JLD2OutputWriter
using Oceananigans.Utils: minutes

model = IncompressibleModel(grid=RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1)))

simulation = Simulation(model, Δt=10minutes, stop_time=30years)

simulation.output_writers[:velocities] = JLD2OutputWriter(model, model.velocities,
                                                          prefix = "averaged_velocity_data",
                                                          schedule = AveragedTimeInterval(4years, window=1year, stride=2))

# output
JLD2OutputWriter scheduled on TimeInterval(4 years):
├── filepath: ./averaged_velocity_data.jld2
├── 3 outputs: (:u, :v, :w) averaged on AveragedTimeInterval(window=1 year, stride=2, interval=4 years)
├── field slicer: FieldSlicer(:, :, :, with_halos=false)
├── array type: Array{Float32}
├── including: [:grid, :coriolis, :buoyancy, :closure]
└── max filesize: Inf YiB
source
Oceananigans.OutputWriters.WindowedTimeAverageType
WindowedTimeAverage(operand, model=nothing; schedule, field_slicer=FieldSlicer())

Returns an object for computing running averages of operand over schedule.window and recurring on schedule.interval, where schedule is an AveragedTimeInterval. During the collection period, averages are computed every schedule.stride iteration.

operand may be a Oceananigans.Field or a function that returns an array or scalar.

Calling wta(model) for wta::WindowedTimeAverage object returns wta.result.

source
Oceananigans.OutputWriters.CheckpointerMethod
Checkpointer(model; schedule,
                    dir = ".",
                 prefix = "checkpoint",
                  force = false,
                verbose = false,
             properties = [:architecture, :boundary_conditions, :grid, :clock, :coriolis,
                           :buoyancy, :closure, :velocities, :tracers, :timestepper]
             )

Construct a Checkpointer that checkpoints the model to a JLD2 file every so often as specified by schedule. The model.clock.iteration is included in the filename to distinguish between multiple checkpoint files.

Note that extra model properties can be safely specified, but removing crucial properties such as :velocities will make restoring from the checkpoint impossible.

The checkpoint file is generated by serializing model properties to JLD2. However, functions cannot be serialized to disk with JLD2. So if a model property contains a reference somewhere in its hierarchy it will not be included in the checkpoint file (and you will have to manually restore them).

Keyword arguments

  • schedule (required): Schedule that determines when to checkpoint.

  • dir: Directory to save output to. Default: "." (current working directory).

  • prefix: Descriptive filename prefixed to all output files. Default: "checkpoint".

  • force: Remove existing files if their filenames conflict. Default: false.

  • verbose: Log what the output writer is doing with statistics on compute/write times and file sizes. Default: false.

  • properties: List of model properties to checkpoint. Some are required.

source
Oceananigans.OutputWriters.restore_from_checkpointMethod
restore_from_checkpoint(filepath; kwargs=Dict())

Restore a model from the checkpoint file stored at filepath. kwargs can be passed to the model constructor, which can be especially useful if you need to manually restore forcing functions or boundary conditions that rely on functions.

source

Time steppers

Simulations

Oceananigans.Simulations.TimeStepWizardMethod
TimeStepWizard(cfl=0.1, max_change=2.0, min_change=0.5, max_Δt=Inf)

A type for calculating adaptive time steps based on capping the CFL number at cfl.

On calling update_Δt!(wizard, model), the TimeStepWizard computes a time-step such that $cfl = max(u/Δx, v/Δy, w/Δz) Δt$, where $max(u/Δx, v/Δy, w/Δz)$ is the maximum ratio between model velocity and along-velocity grid spacing anywhere on the model grid. The new Δt is constrained to change by a multiplicative factor no more than max_change or no less than min_change from the previous Δt, and to be no greater in absolute magnitude than max_Δt.

source
Oceananigans.Simulations.SimulationMethod
Simulation(model; Δt,
     stop_criteria = Function[iteration_limit_exceeded, stop_time_exceeded, wall_time_limit_exceeded],
    stop_iteration = Inf,
         stop_time = Inf,
   wall_time_limit = Inf,
       diagnostics = OrderedDict{Symbol, AbstractDiagnostic}(),
    output_writers = OrderedDict{Symbol, AbstractOutputWriter}(),
          progress = nothing,
iteration_interval = 1,
        parameters = nothing)

Construct an Oceananigans.jl Simulation for a model with time step Δt.

Keyword arguments

  • Δt: Required keyword argument specifying the simulation time step. Can be a Number for constant time steps or a TimeStepWizard for adaptive time-stepping.
  • stop_criteria: A list of functions or callable objects (each taking a single argument, the simulation). If any of the functions return true when the stop criteria is evaluated the simulation will stop.
  • stop_iteration: Stop the simulation after this many iterations.
  • stop_time: Stop the simulation once this much model clock time has passed.
  • wall_time_limit: Stop the simulation if it's been running for longer than this many seconds of wall clock time.
  • progress: A function with a single argument, the simulation. Will be called every iteration_interval iterations. Useful for logging simulation health.
  • iteration_interval: How often to update the time step, check stop criteria, and call progress function (in number of iterations).
  • parameters: Parameters that can be accessed in the progress function.
source

Tubrulence closures

Oceananigans.TurbulenceClosures.AnisotropicBiharmonicDiffusivityType
AnisotropicBiharmonicDiffusivity(FT=Float64; νx=0, νy=0, νz=0, νh=nothing, κx=0, κy=0, κz=0, κh=nothing)

Returns parameters for a fourth-order, anisotropic biharmonic diffusivity closure with constant x-, y, and z-direction biharmonic viscosities νx, νy, and νz, and constant x-, y, and z-direction biharmonic diffusivities κx, κy, and κz, κx, κy, and κz may be NamedTuples with fields corresponding to each tracer, or a single number to be a applied to all tracers.

If νh or κh are provided, then νx = νy = νh or κx = κy = κh.

The tracer flux divergence associated with an anisotropic biharmonic diffusivity is, for example

\[ ∂ᵢ κᵢⱼ ∂ⱼc = - [κx ∂⁴x + κy ∂⁴y + κz ∂⁴z] c\]
source
Oceananigans.TurbulenceClosures.SmagorinskyLillyType
SmagorinskyLilly([FT=Float64;] C=0.23, Pr=1, ν=1.05e-6, κ=1.46e-7)

Return a SmagorinskyLilly type associated with the turbulence closure proposed by Lilly (1962) and Smagorinsky (1958, 1963), which has an eddy viscosity of the form

`νₑ = (C * Δᶠ)² * √(2Σ²) * √(1 - Cb * N² / Σ²) + ν`,

and an eddy diffusivity of the form

`κₑ = (νₑ - ν) / Pr + κ`

where Δᶠ is the filter width, Σ² = ΣᵢⱼΣᵢⱼ is the double dot product of the strain tensor Σᵢⱼ, Pr is the turbulent Prandtl number, and is the total buoyancy gradient, and Cb is a constant the multiplies the Richardson number modification to the eddy viscosity.

Keyword arguments

- `C`  : Model constant
- `Cb` : Buoyancy term multipler (`Cb = 0` turns it off, `Cb ≠ 0` turns it on.
         Typically `Cb=1/Pr`.)
- `Pr` : Turbulent Prandtl numbers for each tracer. Either a constant applied to every
         tracer, or a `NamedTuple` with fields for each tracer individually.
- `ν`  : Constant background viscosity for momentum
- `κ`  : Constant background diffusivity for tracer. Can either be a single number
         applied to all tracers, or `NamedTuple` of diffusivities corresponding to each
         tracer.

References

Smagorinsky, J. "On the numerical integration of the primitive equations of motion for baroclinic flow in a closed region." Monthly Weather Review (1958)

Lilly, D. K. "On the numerical simulation of buoyant convection." Tellus (1962)

Smagorinsky, J. "General circulation experiments with the primitive equations: I. The basic experiment." Monthly weather review (1963)

source
Oceananigans.TurbulenceClosures.∇_κ_∇cMethod
∇_κ_∇c(i, j, k, grid, clock, c, closure, diffusivities)

Return the diffusive flux divergence ∇ ⋅ (κ ∇ c) for the turbulence closure, where c is an array of scalar data located at cell centers.

source
Oceananigans.TurbulenceClosures.IsotropicDiffusivityType
IsotropicDiffusivity(; ν=ν₀, κ=κ₀)

Returns parameters for an isotropic diffusivity model with viscosity ν and thermal diffusivities κ for each tracer field in tracers ν and the fields of κ may be constants or functions of (x, y, z, t), and may represent molecular diffusivities in cases that all flow features are explicitly resovled, or turbulent eddy diffusivities that model the effect of unresolved, subgrid-scale turbulence. κ may be a NamedTuple with fields corresponding to each tracer, or a single number to be a applied to all tracers.

By default, a molecular viscosity of ν₀ = 1.05×10⁻⁶ m² s⁻¹ and a molecular thermal diffusivity of κ₀ = 1.46×10⁻⁷ m² s⁻¹ is used for each tracer. These molecular values are the approximate viscosity and thermal diffusivity for seawater at 20°C and 35 psu, according to Sharqawy et al., "Thermophysical properties of seawater: A review of existing correlations and data" (2010).

source
Oceananigans.TurbulenceClosures.VerstappenAnisotropicMinimumDissipationType
VerstappenAnisotropicMinimumDissipation(FT=Float64; C=1/12, Cν=nothing, Cκ=nothing,
                                        Cb=0.0, ν=ν₀, κ=κ₀)

Returns parameters of type FT for the VerstappenAnisotropicMinimumDissipation turbulence closure.

Keyword arguments

- `C`  : Poincaré constant for both eddy viscosity and eddy diffusivities. `C` is overridden
         for eddy viscosity or eddy diffusivity if `Cν` or `Cκ` are set, respecitvely.
- `Cν` : Poincaré constant for momentum eddy viscosity.
- `Cκ` : Poincaré constant for tracer eddy diffusivities. If one number or function, the same
         number or function is applied to all tracers. If a `NamedTuple`, it must possess
         a field specifying the Poncaré constant for every tracer.
- `Cb` : Buoyancy modification multiplier (`Cb = 0` turns it off, `Cb = 1` turns it on)
- `ν`  : Constant background viscosity for momentum.
- `κ`  : Constant background diffusivity for tracer. If a single number, the same background
         diffusivity is applied to all tracers. If a `NamedTuple`, it must possess a field
         specifying a background diffusivity for every tracer.

By default: C = Cν = Cκ = 1/12, which is appropriate for a finite-volume method employing a second-order advection scheme, Cb = 0, which terms off the buoyancy modification term, the molecular viscosity of seawater at 20 deg C and 35 psu is used for ν, and the molecular diffusivity of heat in seawater at 20 deg C and 35 psu is used for κ.

or may be constant numbers, or functions of x, y, z.

Example

julia> prettydiffusiveclosure = AnisotropicMinimumDissipation(C=1/2) VerstappenAnisotropicMinimumDissipation{Float64} turbulence closure with: Poincaré constant for momentum eddy viscosity Cν: 0.5 Poincaré constant for tracer(s) eddy diffusivit(ies) Cκ: 0.5 Buoyancy modification multiplier Cb: 0.0 Background diffusivit(ies) for tracer(s), κ: 1.46e-7 Background kinematic viscosity for momentum, ν: 1.05e-6

julia> const Δz = 0.5; # grid resolution at surface

julia> surfaceenhancedtracerC(x, y, z) = 1/12 * (1 + exp((z + Δz/2) / 8Δz)) surfaceenhancedtracerC (generic function with 1 method)

julia> fancyclosure = AnisotropicMinimumDissipation(Cκ=surfaceenhancedtracerC) VerstappenAnisotropicMinimumDissipation{Float64} turbulence closure with: Poincaré constant for momentum eddy viscosity Cν: 0.08333333333333333 Poincaré constant for tracer(s) eddy diffusivit(ies) Cκ: surfaceenhancedtracer_C Buoyancy modification multiplier Cb: 0.0 Background diffusivit(ies) for tracer(s), κ: 1.46e-7 Background kinematic viscosity for momentum, ν: 1.05e-6

julia> tracerspecificclosure = AnisotropicMinimumDissipation(Cκ=(c₁=1/12, c₂=1/6)) VerstappenAnisotropicMinimumDissipation{Float64} turbulence closure with: Poincaré constant for momentum eddy viscosity Cν: 0.08333333333333333 Poincaré constant for tracer(s) eddy diffusivit(ies) Cκ: (c₁ = 0.08333333333333333, c₂ = 0.16666666666666666) Buoyancy modification multiplier Cb: 0.0 Background diffusivit(ies) for tracer(s), κ: 1.46e-7 Background kinematic viscosity for momentum, ν: 1.05e-6

References

Vreugdenhil C., and Taylor J. (2018), "Large-eddy simulations of stratified plane Couette flow using the anisotropic minimum-dissipation model", Physics of Fluids 30, 085104.

Verstappen, R. (2018), "How much eddy dissipation is needed to counterbalance the nonlinear production of small, unresolved scales in a large-eddy simulation of turbulence?", Computers & Fluids 176, pp. 276-284.

source
Oceananigans.TurbulenceClosures.∇_κ_∇cMethod
∇_κ_∇c(i, j, k, grid, clock, c, tracer_index, closure, diffusivities)

Return the diffusive flux divergence ∇ ⋅ (κ ∇ c) for the turbulence closure, where c is an array of scalar data located at cell centers.

source
Oceananigans.TurbulenceClosures.BlasiusSmagorinskyType
BlasiusSmagorinsky(FT=Float64; Pr=1.0, ν=1.05e-6, κ=1.46e-7)

Returns a BlasiusSmagorinsky closure object of type FT.

Keyword arguments

- `Pr` : Turbulent Prandtl numbers for each tracer. Either a constant applied to every
         tracer, or a `NamedTuple` with fields for each tracer individually.
- `ν`  : Constant background viscosity for momentum
- `κ`  : Constant background diffusivity for tracer. Can either be a single number
         applied to all tracers, or `NamedTuple` of diffusivities corresponding to each
         tracer.

References

Polton, J. A., and Belcher, S. E. (2007), "Langmuir turbulence and deeply penetrating jets in an unstratified mixed layer." Journal of Geophysical Research: Oceans.

source
Oceananigans.TurbulenceClosures.AnisotropicDiffusivityType
AnisotropicDiffusivity(; νx=ν₀, νy=ν₀, νz=ν₀, κx=κ₀, κy=κ₀, κz=κ₀,
                         νh=nothing, κh=nothing)

Returns parameters for a closure with a diagonal diffusivity tensor with heterogeneous 'anisotropic' components labeled by x, y, z. Each component may be a number or function. The tracer diffusivities κx, κy, and κz may be NamedTuples with fields corresponding to each tracer, or a single number or function to be a applied to all tracers.

If νh or κh are provided, then νx = νy = νh, and κx = κy = κh, respectively.

By default, a viscosity of ν₀ = 1.05×10⁻⁶ m² s⁻¹ is used for all viscosity components and a diffusivity of κ₀ = 1.46×10⁻⁷ m² s⁻¹ is used for all diffusivity components for every tracer. These values are the approximate viscosity and thermal diffusivity for seawater at 20°C and 35 psu, according to Sharqawy et al., "Thermophysical properties of seawater: A review of existing correlations and data" (2010).

source
Oceananigans.TurbulenceClosures.RozemaAnisotropicMinimumDissipationType
RozemaAnisotropicMinimumDissipation(FT=Float64; C=0.33, ν=1.05e-6, κ=1.46e-7)

Returns a RozemaAnisotropicMinimumDissipation closure object of type FT with

* `C` : Poincaré constant
* `ν` : 'molecular' background viscosity
* `κ` : 'molecular' background diffusivity for each tracer

See Rozema et al., " (2015)

source
Oceananigans.TurbulenceClosures.TwoDimensionalLeithType
TwoDimensionalLeith([FT=Float64;] C=0.3, C_Redi=1, C_GM=1)

Return a TwoDimensionalLeith type associated with the turbulence closure proposed by Leith (1965) and Fox-Kemper & Menemenlis (2008) which has an eddy viscosity of the form

`νₑ = (C * Δᶠ)³ * √(ζ² + (∇h ∂z w)²)`

and an eddy diffusivity of the form...

where Δᶠ is the filter width, ζ² = (∂x v - ∂y u)² is the squared vertical vorticity, and C is a model constant.

Keyword arguments

- `C`      : Model constant
- `C_Redi` : Coefficient for down-gradient tracer diffusivity for each tracer.
             Either a constant applied to every tracer, or a `NamedTuple` with fields
             for each tracer individually.
- `C_GM`   : Coefficient for down-gradient tracer diffusivity for each tracer.
             Either a constant applied to every tracer, or a `NamedTuple` with fields
             for each tracer individually.

References

Leith, C. E. (1968). "Diffusion Approximation for Two‐Dimensional Turbulence", The Physics of Fluids 11, 671. doi: 10.1063/1.1691968

Fox‐Kemper, B., & D. Menemenlis (2008), "Can large eddy simulation techniques improve mesoscale rich ocean models?", in Ocean Modeling in an Eddying Regime, Geophys. Monogr. Ser., vol. 177, pp. 319–337. doi: 10.1029/177GM19

Pearson, B. et al. (2017) , "Evaluation of scale-aware subgrid mesoscale eddy models in a global eddy rich model", Ocean Modelling 115, 42-58. doi: 10.1016/j.ocemod.2017.05.007

source
Oceananigans.TurbulenceClosures.∇_κ_∇cMethod
∇_κ_∇c(i, j, k, grid, clock, c, closure, diffusivities)

Return the diffusive flux divergence ∇ ⋅ (κ ∇ c) for the turbulence closure, where c is an array of scalar data located at cell centers.

source

Utilities

Oceananigans.Utils.GiBConstant
GiB

A Float64 constant equal to 1024MiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 50GiB.

source
Oceananigans.Utils.KiBConstant
KiB

A Float64 constant equal to 1024.0. Useful for increasing the clarity of scripts, e.g. max_filesize = 250KiB.

source
Oceananigans.Utils.MiBConstant
MiB

A Float64 constant equal to 1024KiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 100MiB.

source
Oceananigans.Utils.TiBConstant
TiB

A Float64 constant equal to 1024GiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 2TiB.

source
Oceananigans.Utils.dayConstant
day

A Float64 constant equal to 24hours. Useful for increasing the clarity of scripts, e.g. stop_time = 1day.

source
Oceananigans.Utils.daysConstant
days

A Float64 constant equal to 24hours. Useful for increasing the clarity of scripts, e.g. stop_time = 7days.

source
Oceananigans.Utils.hourConstant
hour

A Float64 constant equal to 60minutes. Useful for increasing the clarity of scripts, e.g. Δt = 1hour.

source
Oceananigans.Utils.hoursConstant
hours

A Float64 constant equal to 60minutes. Useful for increasing the clarity of scripts, e.g. Δt = 3hours.

source
Oceananigans.Utils.yearConstant
year

A Float64 constant equal to 365days. Useful for increasing the clarity of scripts, e.g. stop_time = 1year.

source
Oceananigans.Utils.yearsConstant
years

A Float64 constant equal to 365days. Useful for increasing the clarity of scripts, e.g. stop_time = 100years.

source
Oceananigans.Utils.prettytimeMethod
prettytime(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, minutes, hours, days, or years.

source
Oceananigans.Utils.pretty_filesizeFunction
pretty_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).

source
Oceananigans.Utils.with_tracersMethod
with_tracers(tracer_names, initial_tuple, tracer_default)

Create a tuple corresponding to the solution variables u, v, w, and tracer_names. initial_tuple is a NamedTuple that at least has fields u, v, and w, and may have some fields corresponding to the names in tracer_names. tracer_default is a function that produces a default tuple value for each tracer if not included in initial_tuple.

source

Abstract operations

Oceananigans.AbstractOperations.@unaryMacro
@unary op1 op2 op3...

Turn each unary function in the list (op1, op2, op3...) into a unary operator on Oceananigans.Fields for use in AbstractOperations.

Note: a unary function is a function with one argument: for example, sin(x) is a unary function.

Also note: a unary function in Base must be imported to be extended: use import Base: op; @unary op.

Example

julia> using Oceananigans, Oceananigans.Grids, Oceananigans.AbstractOperations

julia> square_it(x) = x^2
square_it (generic function with 1 method)

julia> @unary square_it
Set{Any} with 7 elements:
  :sqrt
  :square_it
  :cos
  :exp
  :-
  :tanh
  :sin

julia> c = Field(Cell, Cell, Cell, CPU(), RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1)));

julia> square_it(c)
UnaryOperation at (Cell, Cell, Cell)
├── grid: RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}(Nx=1, Ny=1, Nz=1)
│   └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [-1.0, 0.0]
└── tree:
    square_it at (Cell, Cell, Cell) via identity
    └── Field located at (Cell, Cell, Cell)
source
Oceananigans.AbstractOperations.@binaryMacro
@binary op1 op2 op3...

Turn each binary function in the list (op1, op2, op3...) into a binary operator on Oceananigans.Fields for use in AbstractOperations.

Note: a binary function is a function with two arguments: for example, +(x, y) is a binary function.

Also note: a binary function in Base must be imported to be extended: use import Base: op; @binary op.

Example

```jldoctest julia> using Oceananigans, Oceananigans.AbstractOperations, Oceananigans.Grids

julia> plusortimes(x, y) = x < 0 ? x + y : x * y plusortimes (generic function with 1 method)

julia> @binary plusortimes 6-element Array{Any,1}: :+ :- :/ :^ :* :plusortimes

julia> c, d = (Field(Cell, Cell, Cell, CPU(), RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1))) for i = 1:2);

julia> plusortimes(c, d) BinaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}(Nx=1, Ny=1, Nz=1) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [-1.0, 0.0] └── tree: plusortimes at (Cell, Cell, Cell) via identity    ├── Field located at (Cell, Cell, Cell)    └── Field located at (Cell, Cell, Cell)

source
Oceananigans.AbstractOperations.@multiaryMacro
@multiary op1 op2 op3...

Turn each multiary operator in the list (op1, op2, op3...) into a multiary operator on Oceananigans.Fields for use in AbstractOperations.

Note that a multiary operator: * is a function with two or more arguments: for example, +(x, y, z) is a multiary function; * must be imported to be extended if part of Base: use import Base: op; @multiary op; * can only be called on Oceananigans.Fields if the "location" is noted explicitly; see example.

Example

```jldoctest julia> using Oceananigans, Oceananigans.Grids, Oceananigans.AbstractOperations

julia> harmonicplus(a, b, c) = 1/3 * (1/a + 1/b + 1/c) harmonicplus(generic function with 1 method)

julia> c, d, e = Tuple(Field(Cell, Cell, Cell, CPU(), RegularCartesianGrid(size=(1, 1, 1), extent=(1, 1, 1))) for i = 1:3);

julia> harmonic_plus(c, d, e) # before magic @multiary transformation BinaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}(Nx=1, Ny=1, Nz=1) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [-1.0, 0.0] └── tree: * at (Cell, Cell, Cell) via identity    ├── 0.3333333333333333    └── + at (Cell, Cell, Cell)       ├── / at (Cell, Cell, Cell) via identity       │   ├── 1       │   └── Field located at (Cell, Cell, Cell)       ├── / at (Cell, Cell, Cell) via identity       │   ├── 1       │   └── Field located at (Cell, Cell, Cell)       └── / at (Cell, Cell, Cell) via identity          ├── 1          └── Field located at (Cell, Cell, Cell)

julia> @multiary harmonicplus Set{Any} with 3 elements: :+ :harmonicplus :*

julia> harmonicplus(c, d, e) MultiaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64, Periodic, Periodic, Bounded}(Nx=1, Ny=1, Nz=1) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [-1.0, 0.0] └── tree: harmonicplus at (Cell, Cell, Cell)    ├── Field located at (Cell, Cell, Cell)    ├── Field located at (Cell, Cell, Cell)    └── Field located at (Cell, Cell, Cell)

source
Oceananigans.AbstractOperations.∂xMethod
∂x(L::Tuple, a::AbstractField)

Return an abstract representation of an x-derivative acting on a followed by interpolation to L, where L is a 3-tuple of Faces and Cells.

source
Oceananigans.AbstractOperations.∂yMethod
∂y(L::Tuple, a::AbstractField)

Return an abstract representation of a y-derivative acting on a followed by interpolation to L, where L is a 3-tuple of Faces and Cells.

source
Oceananigans.AbstractOperations.∂zMethod
∂z(L::Tuple, a::AbstractField)

Return an abstract representation of a z-derivative acting on a followed by interpolation to L, where L is a 3-tuple of Faces and Cells.

source