Library

Documenting the public user interface.

Architectures

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.NoPenetrationType
NoPenetration

A type specifying a no-penetration boundary condition for a velocity component that is normal to a wall.

Thus NoPenetration 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.

A condition may not be specified with a NoPenetration boundary condition.

Note that this differs from a zero Value boundary condition as Value imposes values at the cell centers (and could apply to tracers) while a no-penetration boundary condition only applies to normal velocity components at a wall, where the velocity at the cell face collocated at the wall is known and set to zero.

source
Oceananigans.BoundaryConditions.BoundaryConditionType
BoundaryCondition{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.

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 (the grid's topology is what defined the default boundary conditions that are imposed).

Specific boundary conditions can be applied along the x dimension with the west and east kwargs, along the y-dimension with the south and north kwargs, and along the z-dimension with the bottom and top kwargs.

source
Oceananigans.BoundaryConditions.BoundaryFunctionType
    BoundaryFunction{B, X1, X2}(func, parameters=nothing)

A wrapper for the user-defined boundary condition function func, on the boundary specified by symbol B and at location (X1, X2), and with parameters.

Example

julia> using Oceananigans, Oceananigans.BoundaryConditions, Oceananigans.Fields

julia> toptracerflux = BoundaryFunction{:z, Cell, Cell}((x, y, t) -> cos(2π*x) * cos(t)) (::BoundaryFunction{:z,Cell,Cell,var"#7#8",Nothing}) (generic function with 1 method)

julia> toptracerbc = BoundaryCondition(Flux, toptracerflux);

julia> flux_func(x, y, t, p) = cos(p.k * x) * cos(p.ω * t); # function with parameters

julia> parameterizeduvelocityflux = BoundaryFunction{:z, Face, Cell}(fluxfunc, (k=4π, ω=3.0)) (::BoundaryFunction{:z,Face,Cell,typeof(flux_func),NamedTuple{(:k, :ω),Tuple{Float64,Float64}}}) (generic function with 1 method)

julia> topubc = BoundaryCondition(Flux, parameterizeduvelocity_flux);

source
Oceananigans.BoundaryConditions.TracerBoundaryConditionMethod
TracerBoundaryCondition(bctype, B, args...)

Returns a BoundaryCondition of type bctype, that applies the function func to a tracer on the boundary B, which is one of :x, :y, :z. The boundary function has the signature

`func(ξ, η, t)`

where t is time, and ξ and η are coordinates along the boundary, eg: (y, z) for B = :x, (x, z) for B = :y, or (x, y) for B = :z.

source
Oceananigans.BoundaryConditions.UVelocityBoundaryConditionMethod
UVelocityBoundaryCondition(bctype, B, args...)

Returns a BoundaryCondition of type bctype, that applies the function func to u, the x-velocity field, on the boundary B, which is one of :x, :y, :z. The boundary function has the signature

`func(ξ, η, t)`

where t is time, and ξ and η are coordinates along the boundary, eg: (y, z) for B = :x, (x, z) for B = :y, or (x, y) for B = :z.

source
Oceananigans.BoundaryConditions.VVelocityBoundaryConditionMethod
VVelocityBoundaryCondition(bctype, B, args...)

Returns a BoundaryCondition of type bctype, that applies the function func to v, the y-velocity field, on the boundary B, which is one of :x, :y, :z. The boundary function has the signature

`func(ξ, η, t)`

where t is time, and ξ and η are coordinates along the boundary, eg: (y, z) for B = :x, (x, z) for B = :y, or (x, y) for B = :z.

source
Oceananigans.BoundaryConditions.WVelocityBoundaryConditionMethod
VVelocityBoundaryCondition(bctype, B, args...)

Returns a BoundaryCondition of type bctype, that applies the function func to w, the z-velocity field, on the boundary B, which is one of :x, :y, :z. The boundary function has the signature

`func(ξ, η, t)`

where t is time, and ξ and η are coordinates along the boundary, eg: (y, z) for B = :x, (x, z) for B = :y, or (x, y) for B = :z.

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. If either temperature or salinity are specified, buoyancy is calculated

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
Oceananigans.Buoyancy.RoquetIdealizedNonlinearEquationOfStateType
RoquetIdealizedNonlinearEquationOfState{F, C, T} <: AbstractNonlinearEquationOfState

Parameters associated with the idealized nonlinear equation of state proposed by Roquet et al., "Defining a Simplified yet 'Realistic' Equation of State for Seawater", Journal of Physical Oceanography (2015).

source
Oceananigans.Buoyancy.RoquetIdealizedNonlinearEquationOfStateType
RoquetIdealizedNonlinearEquationOfState([FT=Float64,] flavor, ρ₀=1024.6,
                                        polynomial_coeffs=optimized_roquet_coeffs[flavor])

Returns parameters for the idealized polynomial nonlinear equation of state with reference density ρ₀ and polynomial_coeffs proposed by Roquet et al., "Defining a Simplified yet 'Realistic' Equation of State for Seawater", Journal of Physical Oceanography (2015). The default reference density is ρ₀ = 1024.6 kg m⁻³, the average surface density of seawater in the world ocean.

The flavor of the nonlinear equation of state is a symbol that selects a set of optimized polynomial coefficients defined in Table 2 of Roquet et al., "Defining a Simplified yet 'Realistic' Equation of State for Seawater", Journal of Physical Oceanography (2015), and further discussed in the text surrounding equations (12)–(15). The optimization minimizes errors in estimated horizontal density gradient estiamted from climatological temperature and salinity distributions between the 5 simplified forms chosen by Roquet et. al and the full-fledged TEOS-10 equation of state.

The equations of state define the density anomaly ρ′, and have the polynomial form

`ρ′(T, S, D) = Σᵢⱼₐ Rᵢⱼₐ Tⁱ Sʲ Dᵃ`,

where T is conservative temperature, S is absolute salinity, and D is the geopotential depth, currently just D = -z. The Rᵢⱼₐ are constant coefficients.

Flavors of idealized nonlinear equations of state

- `:linear`: a linear equation of state, `ρ′ = R₁₀₀ * T + R₀₁₀ * S`

- `:cabbeling`: includes quadratic temperature term,
                `ρ′ = R₁₀₀ * T + R₀₁₀ * S + R₀₂₀ * T^2`

- `:cabbeling_thermobaricity`: includes 'thermobaricity' term,
                               `ρ′ = R₁₀₀ * T + R₀₁₀ * S + R₀₂₀ * T^2 + R₀₁₁ * T * D`

- `:freezing`: same as `:cabbeling_thermobaricity` with modified constants to increase
               accuracy near freezing

- `:second_order`: includes quadratic salinity, halibaricity, and thermohaline term,
                   `ρ′ = R₁₀₀ * T + R₀₁₀ * S + R₀₂₀ * T^2 + R₀₁₁ * T * D`
                         + R₂₀₀ * S^2 + R₁₀₁ * S * D + R₁₁₀ * S * T`

Example

julia> using Oceananigans

julia> eos = Oceananigans.RoquetIdealizedNonlinearEquationOfState(:cabbeling);

julia> eos.polynomial_coeffs (R₀₁₀ = -0.0844, R₁₀₀ = 0.7718, R₀₂₀ = -0.004561, R₀₁₁ = 0.0, R₂₀₀ = 0.0, R₁₀₁ = 0.0, R₁₁₀ = 0.0)

References

- Roquet et al., "Defining a Simplified yet 'Realistic' Equation of State for
  Seawater", Journal of Physical Oceanography (2015).

- "Thermodynamic Equation of State for Seawater" (TEOS-10), http://www.teos-10.org
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.NonTraditionalFPlaneType
NonTraditionalFPlane{FT} <: AbstractRotation

A Coriolis implementation that facilitates non-traditional Coriolis terms in the zonal and vertical momentum equations along with the traditional Coriolis terms.

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.BetaPlaneType
BetaPlane([T=Float64;] f₀=nothing, β=nothing,
                       rotation_rate=Ω_Earth, latitude=nothing, radius=R_Earth)

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.

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

source

Diagnostics

Oceananigans.Diagnostics.HorizontalAverageMethod
HorizontalAverage(model, field; frequency=nothing, interval=nothing, return_type=Array)

Construct a HorizontalAverage of field.

After the horizontal average is computed it will be stored in the result property.

The HorizontalAverage can be used as a callable object that computes and returns the horizontal average.

A frequency or interval (or both) can be passed to indicate how often to run this diagnostic if it is part of model.diagnostics. frequency is a number of iterations while interval is a time interval in units of model.clock.time.

A return_type can be used to specify the type returned when the HorizontalAverage is used as a callable object. The default return_type=Array is useful when running a GPU model and you want to save the output to disk by passing it to an output writer.

source
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
Oceananigans.Diagnostics.NaNCheckerMethod
NaNChecker(model; frequency, fields)

Construct a NaNChecker for model. fields should be a Dict{Symbol,Field}. A frequency should be passed to indicate how often to check for NaNs (in number of iterations).

source

Fields

Oceananigans.Fields.AbstractFieldType
AbstractField{X, Y, Z, A, G}

Abstract supertype for fields located at (X, Y, Z) with data stored in a container of type A. The field is defined on a grid G.

source
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 = zeros(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 = zeros(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 = zeros(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 = zeros(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 = zeros(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

Forcing

Oceananigans.Forcing.SimpleForcingMethod
SimpleForcing([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, ω=π))
source
Oceananigans.Forcing.ModelForcingMethod
ModelForcing(; 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 = IncompressibleModel(forcing=ModelForcing(u=u_forcing))

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)

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

See xnodes, ynodes, and znodes.

source
Oceananigans.Grids.xnodesMethod
xnodes(loc, grid)

Returns a view over the interior loc=Cell or loc=Facenodes ongridin the x-direction. ForBoundeddirections,Face` nodes include the boundary points.

See znodes for examples.

source
Oceananigans.Grids.ynodesMethod
ynodes(loc, grid)

Returns a view over the interior loc=Cell or loc=Facenodes ongridin the y-direction. ForBoundeddirections,Face` nodes include the boundary points.

See znodes for examples.

source
Oceananigans.Grids.znodesMethod
znodes(loc, grid)

Returns a view over the interior loc=Cell or loc=Facenodes ongridin the z-direction. ForBoundeddirections,Face` nodes include the boundary points.

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)
1×1×3 view(OffsetArray(reshape(::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}, 1, 1, 5), 1:1, 1:1, 0:4), :, :, 1:3) with eltype Float64 with indices 1:1×1:1×Base.OneTo(3):
[:, :, 1] =
 -0.8333333333333331

[:, :, 2] =
 -0.4999999999999999

[:, :, 3] =
 -0.16666666666666652

julia> zF = znodes(Face, horz_periodic_grid)
1×1×4 view(OffsetArray(reshape(::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}, 1, 1, 6), 1:1, 1:1, 0:5), :, :, 1:4) with eltype Float64 with indices 1:1×1:1×Base.OneTo(4):
[:, :, 1] =
 -1.0

[:, :, 2] =
 -0.6666666666666666

[:, :, 3] =
 -0.33333333333333337

[:, :, 4] =
 -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 3-tuple (Nx, Ny, Nz) prescribing the number of grid points in x, y, z

  • extent: A 3-tuple (Lx, Ly, Lz) prescribing the physical extent of the grid. The origin 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.

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

  • 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. In a Flat direction, derivatives are zero. The default is (Periodic, Periodic, Bounded).

  • halo: A 3-tuple of integers that specifies the size of the halo region of cells surrounding the physical interior in x, y, and z.

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, Flat)` 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, reshaped for broadcasting.

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

Examples

julia> grid = RegularCartesianGrid(size=(32, 32, 32), extent=(1, 2, 3))
RegularCartesianGrid{Float64}
domain: x ∈ [0.0, 1.0], y ∈ [0.0, 2.0], z ∈ [0.0, -3.0]
  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)
julia> grid = RegularCartesianGrid(Float32; size=(32, 32, 16), x=(0, 8), y=(-10, 10), z=(-π, π))
RegularCartesianGrid{Float32}
domain: x ∈ [0.0, 8.0], y ∈ [-10.0, 10.0], z ∈ [3.141592653589793, -3.141592653589793]
  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)
source

Models

Oceananigans.Models.stateMethod
state(model)

Returns a NamedTuple with fields velocities, tracers, diffusivities, tendencies corresponding to NamedTuples of OffsetArrays that reference each of the field's data.

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

Clock{T}(time, iteration)

Keeps track of the current time and iteration number. The time::T can be either a number of a DateTime object.

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.FieldOutputType
FieldOutput([return_type=Array], field)

Returns a FieldOutput type intended for use with the JLD2OutputWriter. Calling FieldOutput(model) returns return_type(field.data.parent).

source
Oceananigans.OutputWriters.JLD2OutputWriterMethod
JLD2OutputWriter(model, outputs; interval=nothing, frequency=nothing, dir=".",
                 prefix="", init=noinit, including=[:grid, :coriolis, :buoyancy, :closure],
                 part=1, max_filesize=Inf, force=false, async=false, verbose=false, jld2_kw=Dict{Symbol, Any}())

Construct a JLD2OutputWriter that writes label, func pairs in outputs (which can be a Dict or NamedTuple) to a JLD2 file, where label is a symbol that labels the output and func is a function of the form func(model) that returns the data to be saved.

Keyword arguments

  • frequency::Int : Save output every n model iterations.
  • interval::Int : Save output every t units of model clock time.
  • dir::String : Directory to save output to. Default: "." (current working directory).
  • prefix::String : Descriptive filename prefixed to all output files. Default: "".
  • init::Function : A function of the form init(file, model) that runs when a JLD2 output file is initialized. Default: noinit(args...) = nothing.
  • including::Array : List of model properties to save with every file. By default, the grid, equation of state, Coriolis parameters, buoyancy parameters, and turbulence closure parameters are saved.
  • part::Int : The starting part number used if max_filesize is finite. Default: 1.
  • max_filesize::Int: 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::Bool : Remove existing files if their filenames conflict. Default: false.
  • async::Bool : Write output asynchronously. Default: false.
  • verbose::Bool : Log what the output writer is doing with statistics on compute/write times and file sizes. Default: false.
  • jld2_kw::Dict : Dict of kwargs to be passed to jldopen when data is written.
source
Oceananigans.OutputWriters.FieldOutputsMethod
FieldOutputs(fields)

Returns a dictionary of FieldOutput objects with key, value pairs corresponding to each name and value in the NamedTuple fields. Intended for use with JLD2OutputWriter.

Examples

julia> output_writer = JLD2OutputWriter(model, FieldOutputs(model.velocities), frequency=1);

julia> keys(output_writer.outputs)
Base.KeySet for a Dict{Symbol,FieldOutput{UnionAll,F} where F} with 3 entries. Keys:
  :w
  :v
  :u
source
Oceananigans.OutputWriters.write_grid_and_attributesMethod
write_grid_and_attributes(model; filename="grid.nc", mode="c",
                          compression=0, attributes=Dict(), slice_kw...)

Writes grid and global attributes to filename. By default writes information to a standalone grid.nc file.

Keyword arguments

  • filename : File name to be saved under.
  • mode: NetCDF file is opened in either clobber ("c") or append ("a") mode. Default: "c".
  • compression: Defines the compression level of data from 0-9. Default: 0.
  • attributes: Global attributes. Default: Dict().
source
Oceananigans.OutputWriters.CheckpointerMethod
Checkpointer(model; frequency=nothing, interval=nothing, 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 frequency or interval. 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 (at least not 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

  • frequency::Int : Save output every n model iterations.
  • interval::Int : Save output every t units of model clock time.
  • dir::String : Directory to save output to. Default: "." (current working directory).
  • prefix::String : Descriptive filename prefixed to all output files. Default: "checkpoint".
  • force::Bool : Remove existing files if their filenames conflict. Default: false.
  • verbose::Bool : Log what the output writer is doing with statistics on compute/write times and file sizes. Default: false.
  • properties::Array: List of model properties to checkpoint.
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

Oceananigans.TimeSteppers.AdamsBashforthTimeStepperType
AdamsBashforthTimeStepper(float_type, arch, grid, tracers, χ=0.125;
                          Gⁿ = TendencyFields(arch, grid, tracers),
                          G⁻ = TendencyFields(arch, grid, tracers))

Return an AdamsBashforthTimeStepper object with tendency fields on arch and grid with AB2 parameter χ. The tendency fields can be specified via optional kwargs.

source
Oceananigans.TimeSteppers.time_step!Method
time_step!(model::IncompressibleModel{<:AdamsBashforthTimeStepper}, Δt; euler=false)

Step forward model one time step Δt with a 2nd-order Adams-Bashforth method and pressure-correction substep. Setting euler=true will take a forward Euler time step.

source

Simulations

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,
progress_frequency = 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 progress_frequency iterations. Useful for logging simulation health.
  • progress_frequency: 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.IsotropicViscosityType
IsotropicViscosity{FT} <: TurbulenceClosure{FT}

Abstract supertype for turbulence closures that are defined by an isotropic viscosity and isotropic diffusivities with model parameters stored as properties of type FT.

source
Oceananigans.TurbulenceClosures.AnisotropicBiharmonicDiffusivityType
AnisotropicBiharmonicDiffusivity(; νh, νv, κh, κv)

Returns parameters for a fourth-order, anisotropic biharmonic diffusivity closure with constant horizontal and vertical biharmonic viscosities νh, νv and constant horizontal and vertical tracer biharmonic diffusivities κh, κv. κh and κv may be NamedTuples with fields corresponding to each tracer, or a single number to be a applied to all tracers. The tracer flux divergence associated with an anisotropic biharmonic diffusivity is, for example

\[ ∂ᵢ κᵢⱼ ∂ⱼc = (κh ∇⁴h + κv ∂⁴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, 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.ConstantIsotropicDiffusivityType
ConstantIsotropicDiffusivity([FT=Float64;] ν, κ)

Returns parameters for a constant isotropic diffusivity model with constant viscosity ν and constant thermal diffusivities κ for each tracer field in tracers ν and the fields of κ 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, 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.ConstantAnisotropicDiffusivityType
ConstantAnisotropicDiffusivity(; νh, νv, κh, κv)

Returns parameters for a constant anisotropic diffusivity closure with constant horizontal and vertical viscosities νh, νv and constant horizontal and vertical tracer diffusivities κh, κv. κh and κv may be NamedTuples with fields corresponding to each tracer, or a single number to be a applied to all tracers.

By default, a viscosity of ν = 1.05×10⁻⁶ m² s⁻¹ is used for both the horizontal and vertical viscosity, and a diffusivity of κ = 1.46×10⁻⁷ m² s⁻¹ is used for the horizontal and vertical diffusivities applied to 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, 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 24hour. Useful for increasing the clarity of scripts, e.g. Δt = 0.5day.

source
Oceananigans.Utils.hourConstant
hour

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

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 (s), minutes (min), hours (hr), or days (day).

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> squareit(x) = x^2 squareit (generic function with 1 method)

julia> @unary squareit 7-element Array{Any,1}: :sqrt :sin :cos :exp :tanh :- :squareit

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

julia> square_it(c) UnaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}} │ ├── size: (1, 1, 16) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [0.0, -1.0] └── tree:

square_it at (Cell, Cell, Cell) via identity └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}

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> 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((1, 1, 16), (1, 1, 1))) for i = 1:2);

julia> plusortimes(c, d) BinaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}} │ ├── size: (1, 1, 16) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [0.0, -1.0] └── tree:

plusortimes at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity ├── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}} └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}

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> harmonicplus(a, b, c) = 1/3 * (1/a + 1/b + 1/c) harmonicplus(generic function with 1 method)

julia> @multiary harmonicplus 3-element Array{Any,1}: :+ :* :harmonicplus

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

julia> harmonic_plus(c, d, e) # this calls the original function, which in turn returns a (correct) operation tree BinaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}} │ ├── size: (1, 1, 16) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [0.0, -1.0] └── tree:

  • at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity

├── 0.3333333333333333 └── + at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity    ├── + at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity    │   ├── / at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity    │   │   ├── 1    │   │   └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}    │   └── / at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity       │   ├── 1       │   └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}    └── / at (Cell, Cell, Cell) via Oceananigans.AbstractOperations.identity       ├── 1       └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}

julia> @at (Cell, Cell, Cell) harmonic_plus(c, d, e) # this returns a MultiaryOperation as expected MultiaryOperation at (Cell, Cell, Cell) ├── grid: RegularCartesianGrid{Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}} │ ├── size: (1, 1, 16) │ └── domain: x ∈ [0.0, 1.0], y ∈ [0.0, 1.0], z ∈ [0.0, -1.0] └── tree:

harmonic_plus at (Cell, Cell, Cell) ├── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}} ├── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}} └── OffsetArrays.OffsetArray{Float64,3,Array{Float64,3}}

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
Oceananigans.AbstractOperations.ComputationMethod
Computation(operation, result; return_type=Array)

Returns a Computation representing an operation performed over the elements of operation.grid and stored in result. return_type specifies the output type when the Computation instances is called as a function.

source