Library
Documenting the public user interface.
Boundary conditions
Oceananigans.BoundaryCondition — TypeBoundaryCondition{C<:BCType}(condition)Construct a boundary condition of type C with a condition that may be given by a number, an array, or a function with signature:
condition(i, j, grid, time, iteration, U, Φ, parameters) = # function definitionthat returns a number and where i and j are indices along the boundary.
Boundary condition types include Periodic, Flux, Value, Gradient, and NoPenetration.
Oceananigans.BoundaryFunction — TypeBoundaryFunction{B, X1, X2, F}A wrapper for user-defined boundary condition functions.
Oceananigans.CoordinateBoundaryConditions — TypeCoordinateBoundaryConditions(left, right)A set of two 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.
Oceananigans.Dirchlet — TypeDirchletAn alias for the Value boundary condition type.
Oceananigans.FieldBoundaryConditions — TypeFieldBoundaryConditionsAn alias for NamedTuple{(:x, :y, :z)} that represents a set of three CoordinateBoundaryConditions applied to a field along x, y, and z.
Oceananigans.FieldBoundaryConditions — MethodFieldBoundaryConditions(x, y, z)Construct a FieldBoundaryConditions using a CoordinateBoundaryCondition for each of the x, y, and z coordinates.
Oceananigans.Flux — TypeFluxA type specifying a boundary condition on the flux of a field.
Oceananigans.Gradient — TypeGradientA type specifying a boundary condition on the derivative or gradient of a field. Also called a Neumann boundary condition.
Oceananigans.Neumann — TypeNeumannAn alias for the Gradient boundary condition type.
Oceananigans.Periodic — TypePeriodicA type specifying a periodic boundary condition.
A condition may not be specified with a Periodic boundary condition.
Oceananigans.Value — TypeValueA type specifying a boundary condition on the value of a field. Also called a Dirchlet boundary condition.
Oceananigans.ChannelBCs — MethodChannelBCs(; north = BoundaryCondition(Flux, nothing),
south = BoundaryCondition(Flux, nothing),
top = BoundaryCondition(Flux, nothing),
bottom = BoundaryCondition(Flux, nothing))Construct FieldBoundaryConditions with Periodic boundary conditions in the x direction and specified north (+y), south (-y), top (+z) and bottom (-z) boundary conditions for u, v, and tracer fields.
ChannelBCs cannot be applied to the the vertical velocity w.
Oceananigans.ChannelSolutionBCs — MethodChannelSolutionBCs(u=ChannelBCs(), ...)Construct SolutionBoundaryConditions for a reentrant channel model configuration with solution fields u, v, w, T, and S specified by keyword arguments.
By default ChannelBCs are applied to u, v, T, and S and ChannelBCs(top=NoPenetrationBC(), bottom=NoPenetrationBC()) is applied to w.
Use ChannelBCs when constructing non-default boundary conditions for u, v, w, T, S.
Oceananigans.HorizontallyPeriodicBCs — MethodHorizontallyPeriodicBCs(; top = BoundaryCondition(Flux, nothing),
bottom = BoundaryCondition(Flux, nothing))Construct FieldBoundaryConditions with Periodic boundary conditions in the x and y directions and specified top (+z) and bottom (-z) boundary conditions for u, v, and tracer fields.
HorizontallyPeriodicBCs cannot be applied to the the vertical velocity w.
Oceananigans.HorizontallyPeriodicSolutionBCs — MethodHorizontallyPeriodicSolutionBCs(u=HorizontallyPeriodicBCs(), ...)Construct SolutionBoundaryConditions for a horizontally-periodic model configuration with solution fields u, v, w, T, and S specified by keyword arguments.
By default HorizontallyPeriodicBCs are applied to u, v, T, and S and HorizontallyPeriodicBCs(top=NoPenetrationBC(), bottom=NoPenetrationBC()) is applied to w.
Use HorizontallyPeriodicBCs when constructing non-default boundary conditions for u, v, w, T, S.
Oceananigans.SolutionBoundaryConditions — MethodSolutionBoundaryConditions(tracers, proposal_bcs)Construct a NamedTuple of FieldBoundaryConditions for a model with fields u, v, w, and tracers from the proposal boundary conditions proposal_bcs, which must contain the boundary conditions on u, v, and w and may contain some or all of the boundary conditions on tracers.
Buoyancy
Oceananigans.BuoyancyTracer — TypeBuoyancyTracer <: AbstractBuoyancy{Nothing}Type indicating that the tracer b represents buoyancy.
Oceananigans.LinearEquationOfState — TypeLinearEquationOfState{FT} <: AbstractEquationOfStateLinear equation of state for seawater.
Oceananigans.LinearEquationOfState — TypeLinearEquationOfState([FT=Float64;] α=1.67e-4, β=7.80e-4)Returns parameters for a linear equation of state for seawater with thermal expansion coefficient α [K⁻¹] and haline contraction coefficient β [ppt⁻¹]. The buoyancy perturbation associated with a linear equation of state is
Default constants are taken from Table 1.2 (page 33) of Vallis, "Atmospheric and Oceanic Fluid Dynamics: Fundamentals and Large-Scale Circulation" (2ed, 2017).
Oceananigans.RoquetIdealizedNonlinearEquationOfState — TypeRoquetIdealizedNonlinearEquationOfState{F, C, T} <: AbstractNonlinearEquationOfStateParameters 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).
Oceananigans.RoquetIdealizedNonlinearEquationOfState — TypeRoquetIdealizedNonlinearEquationOfState([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.orgOceananigans.SeawaterBuoyancy — TypeSeawaterBuoyancy{G, EOS} <: AbstractBuoyancy{EOS}Buoyancy model for temperature- and salt-stratified seawater.
Oceananigans.SeawaterBuoyancy — TypeSeawaterBuoyancy([FT=Float64;] gravitational_acceleration = g_Earth,
equation_of_state = LinearEquationOfState(FT))Returns parameters for a temperature- and salt-stratified seawater buoyancy model with a gravitational_acceleration constant (typically called 'g'), and an equation_of_state that related temperature and salinity (or conservative temperature and absolute salinity) to density anomalies and buoyancy.
Clock
Oceananigans.Clock — TypeClock{T<:Number}
Clock{T}(time, iteration)Keeps track of the current time and iteration number.
Coriolis
Oceananigans.BetaPlane — TypeBetaPlane{T} <: AbstractRotationA parameter object for meridionally increasing Coriolis parameter (f = f₀ + βy).
Oceananigans.BetaPlane — TypeBetaPlane([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.
Oceananigans.FPlane — TypeFPlane{FT} <: AbstractRotationA parameter object for constant rotation around a vertical axis.
Oceananigans.FPlane — TypeFPlane([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.
Diagnostics
Oceananigans.Diagnostics.HorizontalAverage — TypeHorizontalAverage{F, R, P, I, Ω} <: AbstractDiagnosticA diagnostic for computing horizontal average of a field.
Oceananigans.Diagnostics.HorizontalAverage — MethodHorizontalAverage(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.
Oceananigans.Diagnostics.Timeseries — TypeTimeseries{D, Ω, I, T, TT} <: AbstractDiagnosticA diagnostic for collecting and storing timeseries.
Oceananigans.Diagnostics.Timeseries — MethodTimeseries(diagnostic, model; frequency=nothing, interval=nothing)A Timeseries Diagnostic that records a time series of diagnostic(model).
Example
julia> model = Model(grid=RegularCartesianGrid(size=(16, 16, 16), length=(1, 1, 1)));
julia> max_u = Timeseries(FieldMaximum(abs, model.velocities.u), model; frequency=1)
julia> model.diagnostics[:max_u] = max_u; data(model.velocities.u) .= π; time_step!(model, Nt=3, Δt=1e-16)
julia> max_u.data
3-element Array{Float64,1}:
3.141592653589793
3.1415926025389127
3.1415925323439517Oceananigans.Diagnostics.Timeseries — MethodTimeseries(diagnostics::NamedTuple, model; frequency=nothing, interval=nothing)A Timeseries Diagnostic that records a NamedTuple of time series of diag(model) for each diag in diagnostics.
Example
julia> model = Model(grid=RegularCartesianGrid(size=(16, 16, 16), length=(1, 1, 1))); Δt = 1.0;
julia> cfl = Timeseries((adv=AdvectiveCFL(Δt), diff=DiffusiveCFL(Δt)), model; frequency=1);
julia> model.diagnostics[:cfl] = cfl; time_step!(model, Nt=3, Δt=Δt)
julia> cfl.data
(adv = [0.0, 0.0, 0.0, 0.0], diff = [0.0002688, 0.0002688, 0.0002688, 0.0002688])
julia> cfl.adv
4-element Array{Float64,1}:
0.0
0.0
0.0
0.0Oceananigans.Diagnostics.CFL — TypeCFL{D, S}An object for computing the Courant-Freidrichs-Lewy (CFL) number.
Oceananigans.Diagnostics.CFL — MethodCFL(Δ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.
Oceananigans.Diagnostics.AdvectiveCFL — MethodAdvectiveCFL(Δ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 = Model(grid=RegularCartesianGrid(size=(16, 16, 16), length=(8, 8, 8)));
julia> cfl = AdvectiveCFL(1.0);
julia> data(model.velocities.u) .= π;
julia> cfl(model)
6.283185307179586Oceananigans.Diagnostics.DiffusiveCFL — MethodDiffusiveCFL(Δ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 = Model(grid=RegularCartesianGrid(size=(16, 16, 16), length=(1, 1, 1)));
julia> dcfl = DiffusiveCFL(0.1);
julia> dcfl(model)
2.688e-5Oceananigans.Diagnostics.FieldMaximum — TypeFieldMaximum(mapping, field)An object for calculating the maximum of a mapping function applied element-wise to field.
Examples
julia> model = Model(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);Oceananigans.Diagnostics.NaNChecker — TypeNaNChecker{F} <: AbstractDiagnosticA diagnostic that checks for NaN values and aborts the simulation if any are found.
Oceananigans.Diagnostics.NaNChecker — MethodNaNChecker(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).
Fields
Oceananigans.Field — TypeField{X, Y, Z, A, G} <: AbstractLocatedField{X, Y, Z, A, G}A field defined at the location (X, Y, Z) which can be either Cell or Face.
Oceananigans.Field — MethodField(X, Y, Z, arch::AbstractArchitecture, grid)Construct a Field on architecture arch and grid at location X, Y, Z, where each of X, Y, Z is Cell or Face.
Oceananigans.Field — MethodField(L::Tuple, data::AbstractArray, grid)Construct a Field on grid using the array data with location defined by the tuple L of length 3 whose elements are Cell or Face.
Oceananigans.Field — MethodField(L::Tuple, arch::AbstractArchitecture, grid)Construct a Field on architecture arch and grid at location L, where L is a tuple of Cell or Face types.
Oceananigans.CellField — MethodCellField([T=eltype(grid)], arch, grid)Return a Field{Cell, Cell, Cell} on architecture arch and grid. Used for tracers and pressure fields.
Oceananigans.FaceFieldX — MethodFaceFieldX([T=eltype(grid)], arch, grid)Return a Field{Face, Cell, Cell} on architecture arch and grid. Used for the x-velocity field.
Oceananigans.FaceFieldY — MethodFaceFieldY([T=eltype(grid)], arch, grid)Return a Field{Cell, Face, Cell} on architecture arch and grid. Used for the y-velocity field.
Oceananigans.FaceFieldZ — MethodFaceFieldZ([T=eltype(grid)], arch, grid)Return a Field{Cell, Cell, Face} on architecture arch and grid. Used for the z-velocity field.
Oceananigans.interior — MethodReturns a view over the interior points of the field.data.
Oceananigans.set! — MethodSet the CPU field u to the array v.
Oceananigans.set! — MethodSet the CPU field u data to the function f(x, y, z).
Oceananigans.set! — Methodset!(model; kwargs...)Set velocity and tracer fields of model. The keyword arguments kwargs... take the form name=data, where name refers to one of the fields of model.velocities or model.tracers, and the data may be an array, a function with arguments (x, y, z), or any data type for which a set!(ϕ::AbstractField, data) function exists.
Example
model = Model(grid=RegularCartesianGrid(size=(32, 32, 32), length=(1, 1, 1))
# Set u to a parabolic function of z, v to random numbers damped
# at top and bottom, and T to some silly array of half zeros,
# half random numbers.
u₀(x, y, z) = z/model.grid.Lz * (1 + z/model.grid.Lz)
v₀(x, y, z) = 1e-3 * rand() * u₀(x, y, z)
T₀ = rand(size(model.grid)...)
T₀[T₀ .< 0.5] .= 0
set!(model, u=u₀, v=v₀, T=T₀)Forcing
Oceananigans.SimpleForcing — TypeSimpleForcing{X, Y, Z, F, P}Callable object for specifying 'simple' forcings of x, y, z, t and optionally parameters of type P at location X, Y, Z.
Oceananigans.SimpleForcing — MethodSimpleForcing([location=(Cell, Cell, Cell),] forcing; parameters=nothing)Construct forcing for a field at location using forcing::Function, and optionally with parameters. If parameters=nothing, forcing must have the signature
`forcing(x, y, z, t)`;otherwise it must have the signature
`forcing(x, y, z, t, parameters)`.Examples
julia> const a = 2.1
julia> fun_forcing(x, y, z, t) = a * exp(z) * cos(t)
julia> u_forcing = SimpleForcing(fun_forcing)
julia> parameterized_forcing(x, y, z, t, p) = p.μ * exp(z/p.λ) * cos(p.ω*t)
julia> v_forcing = SimpleForcing(parameterized_forcing, parameters=(μ=42, λ=0.1, ω=π))Oceananigans.ModelForcing — MethodModelForcing(; u=zeroforcing, v=zeroforcing, w=zeroforcing, tracer_forcings...)Return a named tuple of forcing functions for each solution field.
Example
julia> u_forcing = SimpleForcing((x, y, z, t) -> exp(z) * cos(t))
julia> model = Model(forcing=ModelForcing(u=u_forcing))
Grids
Oceananigans.Grids.RegularCartesianGrid — TypeRegularCartesianGrid{FT<:AbstractFloat, R<:AbstractRange} <: AbstractGrid{FT}A Cartesian grid with with constant grid spacings Δx, Δy, and Δz between cell centers and cell faces.
Oceananigans.Grids.RegularCartesianGrid — TypeRegularCartesianGrid([FT=Float64]; size, length, x, y, z)Creates a RegularCartesianGrid with size = (Nx, Ny, Nz) grid points.
The physical length 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 length argument, e.g. length=(Lx, Ly, Lz) which specifies the length of each dimension in which case 0 ≤ x ≤ Lx, 0 ≤ y ≤ Ly, and -Lz ≤ z ≤ 0.
Constants are stored using floating point values of type FT.
Grid properties
(xC, yC, zC)::AbstractRange: (x, y, z) coordinates of cell centers(xF, yF, zF)::AbstractRange: (x, y, z) coordinates of cell faces(Hx, Hy, Hz)::Int: Halo size in the (x, y, z)-direction(Tx, Ty, Tz)::Int: "Total" grid size (interior + halo points) in the (x, y, z)-direction
Examples
julia> grid = RegularCartesianGrid(size=(32, 32, 32), length=(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)Models
Oceananigans.Model — MethodModel(;
grid,
architecture = CPU(),
float_type = Float64,
tracers = (:T, :S),
closure = ConstantIsotropicDiffusivity(float_type, ν=ν₀, κ=κ₀),
clock = Clock{float_type}(0, 0),
buoyancy = SeawaterBuoyancy(float_type),
coriolis = nothing,
surface_waves = nothing,
forcing = ModelForcing(),
boundary_conditions = HorizontallyPeriodicSolutionBCs(),
output_writers = OrderedDict{Symbol, AbstractOutputWriter}(),
diagnostics = OrderedDict{Symbol, AbstractDiagnostic}(),
parameters = nothing,
velocities = VelocityFields(architecture, grid),
pressures = PressureFields(architecture, grid),
diffusivities = TurbulentDiffusivities(architecture, grid, tracernames(tracers), closure),
timestepper = :AdamsBashforth,
poisson_solver = PoissonSolver(architecture, PoissonBCs(boundary_conditions), grid)
)Construct an Oceananigans.jl model on grid.
Keyword arguments
grid: (required) The resolution and discrete geometry on whichmodelis solved.architecture:CPU()orGPU(). The computer architecture used to time-stepmodel.float_type:Float32orFloat64. The floating point type used formodeldata.closure: The turbulence closure formodel. SeeTurbulenceClosures.buoyancy: Buoyancy model parameters.coriolis: Parameters for the background rotation rate of the model.forcing: User-defined forcing functions that contribute to solution tendencies.boundary_conditions: User-defined boundary conditions for model fields. Can be eitherSolutionBoundaryConditionsorModelBoundaryConditions. SeeBoundaryConditions,HorizontallyPeriodicSolutionBCs, andChannelSolutionBCs.parameters: User-defined parameters for use in user-defined forcing functions and boundary condition functions.
Oceananigans.ChannelModel — MethodChannelModel(; kwargs...)Construct a Model with walls in the y-direction. This is done by imposing FreeSlip boundary conditions in the y-direction instead of Periodic.
kwargs are passed to the regular Model constructor.
Oceananigans.NonDimensionalModel — MethodNonDimensionalModel(; N, L, Re, Pr=0.7, Ro=Inf, float_type=Float64, kwargs...)Construct a "Non-dimensional" Model with resolution N, domain extent L, precision float_type, and the four non-dimensional numbers:
* `Re = U λ / ν` (Reynolds number)
* `Pr = U λ / κ` (Prandtl number)
* `Ro = U / f λ` (Rossby number)for characteristic velocity scale U, length-scale λ, viscosity ν, tracer diffusivity κ, and Coriolis parameter f. Buoyancy is scaled with λ U², so that the Richardson number is Ri=B, where B is a non-dimensional buoyancy scale set by the user via initial conditions or forcing.
Note that N, L, and Re are required.
Additional kwargs are passed to the regular Model constructor.
Output writers
Oceananigans.OutputWriters.FieldOutput — TypeFieldOutput([return_type=Array], field)Returns a FieldOutput type intended for use with the JLD2OutputWriter. Calling FieldOutput(model) returns return_type(field.data.parent).
Oceananigans.OutputWriters.JLD2OutputWriter — TypeJLD2OutputWriter{F, I, O, IF, IN, KW} <: AbstractOutputWriterAn output writer for writing to JLD2 files.
Oceananigans.OutputWriters.JLD2OutputWriter — MethodJLD2OutputWriter(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 everynmodel iterations.interval::Int: Save output everytunits 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 forminit(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 ifmax_filesizeis finite. Default: 1.max_filesize::Int: The writer will stop writing to the output file once the file size exceedsmax_filesize, and write to a new one with a consistent naming scheme ending inpart1,part2, etc. Defaults toInf.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 tojldopenwhen data is written.
Oceananigans.OutputWriters.FieldOutputs — MethodFieldOutputs(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
:uOceananigans.OutputWriters.NetCDFOutputWriter — TypeNetCDFOutputWriter <: AbstractOutputWriterAn output writer for writing to NetCDF files.
Oceananigans.OutputWriters.write_grid — Methodwrite_grid(model; filename="grid.nc", mode="c",
compression=0, attributes=Dict(), slice_kw...)Writes a grid.nc file that contains all the dimensions of the domain.
Keyword arguments
- `filename::String` : File name to be saved under
- `mode::String` : Netcdf file is opened in either clobber ("c") or append ("a") mode (Default: "c" )
- `compression::Int` : Defines the compression level of data (0-9, default 0)
- `attributes::Dict` : Attributes (default: Dict())Oceananigans.OutputWriters.write_output — Methodwrite_output(model, OutputWriter)For internal user only. Writes output to the netcdf file at specified intervals. Increments the Time dimension every time an output is written to the file.
Oceananigans.OutputWriters.Checkpointer — TypeCheckpointer{I, T, P, A} <: AbstractOutputWriterAn output writer for checkpointing models to a JLD2 file from which models can be restored.
Oceananigans.OutputWriters.Checkpointer — MethodCheckpointer(model; frequency=nothing, interval=nothing, dir=".", prefix="checkpoint",
force=false, verbose=false,
properties = [:architecture, :boundary_conditions, :grid, :clock,
:eos, :constants, :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 everynmodel iterations.interval::Int: Save output everytunits 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.
Oceananigans.OutputWriters.restore_from_checkpoint — Methodrestore_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.
Time steppers
Oceananigans.TimeSteppers.time_step! — Methodtime_step!(model; Nt, Δt, kwargs...)Step forward model Nt time steps with step size Δt.
The kwargs are passed to the time_step! function specific to model.timestepper.
Oceananigans.TimeSteppers.compute_w_from_continuity! — MethodCompute the vertical velocity w by integrating the continuity equation from the bottom upwards
`w^{n+1} = -∫ [∂/∂x (u^{n+1}) + ∂/∂y (v^{n+1})] dz`Oceananigans.TimeSteppers.AdamsBashforthTimeStepper — TypeAdamsBashforthTimeStepper(float_type, arch, grid, tracers, χ)Return an AdamsBashforthTimeStepper object with tendency fields on arch and grid and AB2 parameter χ.
Oceananigans.TimeSteppers.time_step! — Methodtime_step!(model{<:AdamsBashforthTimeStepper}, Nt, Δt; init_with_euler=true)Step forward model Nt time steps with step size Δt with an Adams-Bashforth timestepping method.
Oceananigans.TimeSteppers.time_step! — MethodStep forward one time step with a 2nd-order Adams-Bashforth method and pressure-correction substep.
Tubrulence closures
Oceananigans.TurbulenceClosures.AnisotropicMinimumDissipation — TypeAnisotropicMinimumDissipationAn alias for VerstappenAnisotropicMinimumDissipation.
Oceananigans.TurbulenceClosures.ConstantSmagorinsky — TypeConstantSmagorinskyAn alias for SmagorinskyLilly.
Oceananigans.TurbulenceClosures.IsotropicViscosity — TypeIsotropicViscosity{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.
Oceananigans.TurbulenceClosures.∂ⱼ_2ν_Σ₁ⱼ — Method∂ⱼ_2ν_Σ₁ⱼ(i, j, k, grid, closure, U, diffusivities)Return the $x$-component of the turbulent diffusive flux divergence:
∂x(2 ν Σ₁₁) + ∂y(2 ν Σ₁₁) + ∂z(2 ν Σ₁₁)
at the location fcc.
Oceananigans.TurbulenceClosures.∂ⱼ_2ν_Σ₂ⱼ — Method∂ⱼ_2ν_Σ₂ⱼ(i, j, k, grid, closure, U, diffusivities)Return the $y$-component of the turbulent diffusive flux divergence:
∂x(2 ν Σ₂₁) + ∂y(2 ν Σ₂₂) + ∂z(2 ν Σ₂₂)
at the location ccf.
Oceananigans.TurbulenceClosures.∂ⱼ_2ν_Σ₃ⱼ — Method∂ⱼ_2ν_Σ₃ⱼ(i, j, k, grid, closure, diffusivities)Return the $z$-component of the turbulent diffusive flux divergence:
∂x(2 ν Σ₃₁) + ∂y(2 ν Σ₃₂) + ∂z(2 ν Σ₃₃)
at the location ccf.
Oceananigans.TurbulenceClosures.cell_diffusion_timescale — MethodReturns the time-scale for diffusion on a regular grid across a single grid cell.
Oceananigans.TurbulenceClosures.AnisotropicBiharmonicDiffusivity — TypeAnisotropicBiharmonicDiffusivity{FT, KH, KV}Parameters for anisotropic biharmonic diffusivity models.
Oceananigans.TurbulenceClosures.AnisotropicBiharmonicDiffusivity — TypeAnisotropicBiharmonicDiffusivity(; ν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
Oceananigans.TurbulenceClosures.SmagorinskyLilly — TypeSmagorinskyLilly{FT} <: AbstractSmagorinsky{FT}Parameters for the Smagorinsky-Lilly turbulence closure.
Oceananigans.TurbulenceClosures.SmagorinskyLilly — TypeSmagorinskyLilly([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 N² 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)
Oceananigans.TurbulenceClosures.∇_κ_∇c — Method∇_κ_∇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.
Oceananigans.TurbulenceClosures.ConstantIsotropicDiffusivity — TypeConstantIsotropicDiffusivity{FT, K}Parameters for constant isotropic diffusivity models.
Oceananigans.TurbulenceClosures.ConstantIsotropicDiffusivity — TypeConstantIsotropicDiffusivity([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).
Oceananigans.TurbulenceClosures.VerstappenAnisotropicMinimumDissipation — TypeVerstappenAnisotropicMinimumDissipation{FT} <: AbstractAnisotropicMinimumDissipation{FT}Parameters for the anisotropic minimum dissipation large eddy simulation model proposed by Verstappen (2018) and described by Vreugdenhil & Taylor (2018).
Oceananigans.TurbulenceClosures.VerstappenAnisotropicMinimumDissipation — TypeVerstappenAnisotropicMinimumDissipation(FT=Float64; C=1/12, Cb=0.0, ν=ν₀, κ=κ₀)Returns parameters of type FT for the VerstappenAnisotropicMinimumDissipation turbulence closure.
Keyword arguments
- `C` : Poincaré constant
- `Cb` : Buoyancy modification multiplier (`Cb = 0` turns it off, `Cb = 1` turns it on)
- `ν` : 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.By default, 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, and molecular values are used for ν and κ.
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.
Oceananigans.TurbulenceClosures.∇_κ_∇c — Method∇_κ_∇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.
Oceananigans.TurbulenceClosures.BlasiusSmagorinsky — TypeBlasiusSmagorinsky{ML, FT}Parameters for the version of the Smagorinsky closure used in the UK Met Office code Blasius, according to Polton and Belcher (2007).
Oceananigans.TurbulenceClosures.BlasiusSmagorinsky — TypeBlasiusSmagorinsky(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.
Oceananigans.TurbulenceClosures.ConstantAnisotropicDiffusivity — TypeConstantAnisotropicDiffusivity{FT, KH, KV}Parameters for constant anisotropic diffusivity models.
Oceananigans.TurbulenceClosures.ConstantAnisotropicDiffusivity — TypeConstantAnisotropicDiffusivity(; ν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).
Oceananigans.TurbulenceClosures.RozemaAnisotropicMinimumDissipation — TypeRozemaAnisotropicMinimumDissipation(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 tracerSee Rozema et al., " (2015)
Oceananigans.TurbulenceClosures.TwoDimensionalLeith — TypeTwoDimensionalLeith{FT} <: AbstractLeith{FT}Parameters for the 2D Leith turbulence closure.
Oceananigans.TurbulenceClosures.TwoDimensionalLeith — TypeTwoDimensionalLeith([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
Oceananigans.TurbulenceClosures.∇_κ_∇c — Method∇_κ_∇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.
Utilities
Oceananigans.GiB — ConstantGiBA Float64 constant equal to 1024MiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 50GiB.
Oceananigans.KiB — ConstantKiBA Float64 constant equal to 1024.0. Useful for increasing the clarity of scripts, e.g. max_filesize = 250KiB.
Oceananigans.MiB — ConstantMiBA Float64 constant equal to 1024KiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 100MiB.
Oceananigans.TiB — ConstantTiBA Float64 constant equal to 1024GiB. Useful for increasing the clarity of scripts, e.g. max_filesize = 2TiB.
Oceananigans.day — ConstantdayA Float64 constant equal to 24hour. Useful for increasing the clarity of scripts, e.g. Δt = 0.5day.
Oceananigans.hour — ConstanthourA Float64 constant equal to 60minute. Useful for increasing the clarity of scripts, e.g. Δt = 3hour.
Oceananigans.minute — ConstantminuteA Float64 constant equal to 60second. Useful for increasing the clarity of scripts, e.g. Δt = 15minute.
Oceananigans.second — ConstantsecondA Float64 constant equal to 1.0. Useful for increasing the clarity of scripts, e.g. Δt = 1second.
Oceananigans.pretty_filesize — Functionpretty_filesize(s, suffix="B")Convert a floating point value s representing a file size to a more human-friendly formatted string with one decimal places with a suffix defaulting to "B". Depending on the value of s the string will be formatted to show s using an SI prefix from bytes, kiB (1024 bytes), MiB (1024² bytes), and so on up to YiB (1024⁸ bytes).
Oceananigans.prettytime — Methodprettytime(t)Convert a floating point value t representing an amount of time in seconds to a more human-friendly formatted string with three decimal places. Depending on the value of t the string will be formatted to show t in nanoseconds (ns), microseconds (μs), milliseconds (ms), seconds (s), minutes (min), hours (hr), or days (day).
Oceananigans.update_Δt! — Methodupdate_Δt!(wizard, model)Compute wizard.Δt given the velocities and diffusivities of model, and the parameters of wizard.
Abstract operations
Oceananigans.AbstractOperations.@unary — Macro@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}}
Oceananigans.AbstractOperations.@binary — Macro@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}}
Oceananigans.AbstractOperations.@multiary — Macro@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}}
Oceananigans.AbstractOperations.∂x — MethodReturn the x-derivative function acting at (X, Any, Any).
Oceananigans.AbstractOperations.∂x — Method∂x(a::Oceananigans.AbstractLocatedField)Return an abstract representation of a x-derivative acting on a.
Oceananigans.AbstractOperations.∂x — Method∂x(L::Tuple, a::Oceananigans.AbstractLocatedField)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.
Oceananigans.AbstractOperations.∂y — MethodReturn the y-derivative function acting at (Any, Y, Any).
Oceananigans.AbstractOperations.∂y — Method∂y(a::Oceananigans.AbstractLocatedField)Return an abstract representation of a y-derivative acting on a.
Oceananigans.AbstractOperations.∂y — Method∂y(L::Tuple, a::Oceananigans.AbstractLocatedField)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.
Oceananigans.AbstractOperations.∂z — MethodReturn the z-derivative function acting at (Any, Any, Z).
Oceananigans.AbstractOperations.∂z — Method∂z(a::Oceananigans.AbstractLocatedField)Return an abstract representation of a z-derivative acting on a.
Oceananigans.AbstractOperations.∂z — Method∂z(L::Tuple, a::Oceananigans.AbstractLocatedField)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.
Oceananigans.AbstractOperations.Computation — TypeComputation{R, T, O, G}Represents an operation performed over the elements of a field.
Oceananigans.AbstractOperations.Computation — Method(computation::Computation)(args...)Performs the compute(computation) and returns the result if isnothing(return_type), or the result after being converted to return_type.
Oceananigans.AbstractOperations.Computation — MethodComputation(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.
Oceananigans.AbstractOperations.compute! — Methodcompute!(computation::Computation)Perform a computation. The result is stored in computation.result.
Oceananigans.AbstractOperations.@at — Macro@at location abstract_operationModify the abstract_operation so that it returns values at location, where location is a 3-tuple of Faces and Cells.