Library
Documenting the public user interface.
Boundary conditions
Oceananigans.BoundaryCondition — Type.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 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.CoordinateBoundaryConditions — Type.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.
Oceananigans.Dirchlet — Type.DirchletAn alias for the Value boundary condition type.
Oceananigans.FieldBoundaryConditions — Type.FieldBoundaryConditionsAn alias for NamedTuple{(:x, :y, :z)} that represents a set of three CoordinateBoundaryConditions applied to a field along x, y, and z.
Oceananigans.FieldBoundaryConditions — Method.FieldBoundaryConditions(x, y, z)Construct a FieldBoundaryConditions using a CoordinateBoundaryCondition for each of the x, y, and z coordinates.
Oceananigans.Flux — Type.FluxA type specifying a boundary condition on the flux of a field.
Oceananigans.Gradient — Type.GradientA type specifying a boundary condition on the derivative or gradient of a field. Also called a Neumann boundary condition.
Oceananigans.Neumann — Type.NeumannAn alias for the Gradient boundary condition type.
Oceananigans.Periodic — Type.PeriodicA type specifying a periodic boundary condition.
A condition may not be specified with a Periodic boundary condition.
Oceananigans.Value — Type.ValueA type specifying a boundary condition on the value of a field. Also called a Dirchlet boundary condition.
Oceananigans.ChannelBCs — Method.ChannelBCs(; 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 — Method.ChannelSolutionBCs(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 — Method.HorizontallyPeriodicBCs(; 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 — Method.HorizontallyPeriodicSolutionBCs(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 — Method.SolutionBoundaryConditions(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 — Type.BuoyancyTracer <: AbstractBuoyancy{Nothing}Type indicating that the tracer b represents buoyancy.
Oceananigans.LinearEquationOfState — Type.LinearEquationOfState{FT} <: AbstractEquationOfStateLinear equation of state for seawater.
Oceananigans.LinearEquationOfState — Type.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
Default constants are taken from Table 1.2 (page 33) of Vallis, "Atmospheric and Oceanic Fluid Dynamics: Fundamentals and Large-Scale Circulation" (2ed, 2017).
Oceananigans.SeawaterBuoyancy — Type.SeawaterBuoyancy{G, EOS} <: AbstractBuoyancy{EOS}Buoyancy model for temperature- and salt-stratified seawater.
Oceananigans.SeawaterBuoyancy — Type.SeawaterBuoyancy([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 — Type.Clock{T<:Number}
Clock{T}(time, iteration)Keeps track of the current time and iteration number.
Coriolis
Oceananigans.FPlane — Type.FPlane{T} <: AbstractRotationA parameter object for constant rotation around a vertical axis.
Oceananigans.FPlane — Type.FPlane([T=Float64;] f)Returns a parameter object for constant rotation at the angular frequency 2f, and therefore with background vorticity f, around a vertical axis.
Also called FPlane, after the "f-plane" approximation for the local effect of Earth's rotation in a planar coordinate system tangent to the Earth's surface.
Diagnostics
Oceananigans.CFL — Type.CFL{D, S}An object for computing the Courant-Freidrichs-Lewy (CFL) number.
Oceananigans.CFL — Method.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
Oceananigans.FieldMaximum — Type.FieldMaximum(mapping, field)An object for calculating the maximum of a mapping function applied element-wise to field.
Examples
julia> model = BasicModel(N=(16, 16, 16), L=(1, 1, 1));
julia> max_abs_u = FieldMaximum(abs, model.velocities.u);
julia> max_w² = FieldMaximum(x->x^2, model.velocities.w);Oceananigans.HorizontalAverage — Type.HorizontalAverage{F, R, P, I, Ω} <: AbstractDiagnosticA diagnostic for computing horizontal average of a field or the product of multiple fields.
Oceananigans.HorizontalAverage — Method.HorizontalAverage(model, fields; frequency=nothing, interval=nothing, return_type=Array)Construct a HorizontalAverage diagnostic for model.
After the horizontal average is computed it will be stored in the profile property.
The HorizontalAverage can be used as a callable object that computes and returns the horizontal average.
If a single field is passed to fields the the horizontal average of that single field will be computed. If multiple fields are passed to fields, then the horizontal average of their product will be computed.
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.
Warning
Right now taking products of multiple fields does not take into account their locations on the staggered grid and no attempt is made to interpolate all the different fields onto a common location before calculating the product.
Oceananigans.NaNChecker — Type.NaNChecker{F} <: AbstractDiagnosticA diagnostic that checks for NaN values and aborts the simulation if any are found.
Oceananigans.NaNChecker — Method.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).
Oceananigans.Timeseries — Type.Timeseries{D, Ω, I, T, TT} <: AbstractDiagnosticA diagnostic for collecting and storing timeseries.
Oceananigans.Timeseries — Method.Timeseries(diagnostic, model; frequency=nothing, interval=nothing)A Timeseries Diagnostic that records a time series of diagnostic(model).
Example
julia> model = BasicModel(N=(16, 16, 16), L=(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.Timeseries — Method.Timeseries(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 = BasicModel(N=(16, 16, 16), L=(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.AdvectiveCFL — Method.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 = BasicModel(N=(16, 16, 16), L=(8, 8, 8));
julia> cfl = AdvectiveCFL(1.0);
julia> data(model.velocities.u) .= π;
julia> cfl(model)
6.283185307179586Oceananigans.DiffusiveCFL — Method.DiffusiveCFL(Δt)Returns an object for computing the 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.
Example
julia> model = BasicModel(N=(16, 16, 16), L=(1, 1, 1));
julia> cfl = DiffusiveCFL(0.1);
julia> cfl(model)
2.688e-5Fields
Oceananigans.CellField — Type.CellFieldA field defined at the cell centers. Used for pressure and tracers.
Oceananigans.CellField — Method.CellField([T=eltype(grid)], arch, grid)Return a CellField on architecture arch and grid.
Oceananigans.FaceFieldX — Type.FaceFieldXA field defined at the faces along the x-direction. Used for horizontal velocity u.
Oceananigans.FaceFieldX — Method.FaceFieldX([T=eltype(grid)], arch, grid)Return a FaceFieldX on architecture arch and grid.
Oceananigans.FaceFieldY — Type.FaceFieldYA field defined at the faces along the y-direction. Used for horizontal velocity v.
Oceananigans.FaceFieldY — Method.FaceFieldY([T=eltype(grid)], arch, grid)Return a FaceFieldY on architecture arch and grid.
Oceananigans.FaceFieldZ — Type.FaceFieldYA field defined at the faces along the z-direction. Used for vertical velocity w.
Oceananigans.FaceFieldZ — Method.FaceFieldZ([T=eltype(grid)], arch, grid)Return a FaceFieldZ on architecture arch and grid.
Oceananigans.Field — Type.Field{LX, LY, LZ, A, G} <: AbstractField{A, G}A field defined at the location (LX, LY, LZ) which can be either Cell or Face.
Oceananigans.Field — Method.Field(Lx, Ly, Lz, data, grid)Construct a Field on grid using the array data with location defined by Lx, Ly, and Lz which are Cell or Face.
Oceananigans.Field — Method.Field(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 — Method.Field(L::Tuple, arch::AbstractArchitecture, grid)Construct a Field on architecture arch and grid with location defined by the tuple L of length 3 whose elements are Cell or Face.
Oceananigans.data — Method.Returns a view over the interior points of the field.data.
Oceananigans.set! — Method.Set the CPU field u data to the function f(x, y, z).
Oceananigans.set! — Method.Set the CPU field u to the array v.
Oceananigans.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 = Model(grid=RegularCartesianGrid(N=(32, 32, 32), L=(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 — Type.SimpleForcing{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 — Method.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, ω=π))Oceananigans.ModelForcing — Method.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 = Model(forcing=ModelForcing(u=u_forcing))
Grids
Oceananigans.RegularCartesianGrid — Type.RegularCartesianGrid{T<:AbstractFloat, R<:AbstractRange} <: AbstractGrid{T}A Cartesian grid with with constant grid spacings Δx, Δy, and Δz between cell centers and cell faces.
Oceananigans.RegularCartesianGrid — Method.RegularCartesianGrid([T=Float64]; N, L)Creates a RegularCartesianGrid with N = (Nx, Ny, Nz) grid points and domain size L = (Lx, Ly, Lz), where constants are stored using floating point values of type T.
Additional 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(N=(32, 32, 32), L=(1, 1, 1))
RegularCartesianGrid{Float64}
resolution (Nx, Ny, Nz) = (32, 32, 32)
halo size (Hx, Hy, Hz) = (1, 1, 1)
domain (Lx, Ly, Lz) = (1.0, 1.0, 1.0)
grid spacing (Δx, Δy, Δz) = (0.03125, 0.03125, 0.03125)julia> grid = RegularCartesianGrid(Float32; N=(32, 32, 16), L=(8, 8, 2))
RegularCartesianGrid{Float32}
resolution (Nx, Ny, Nz) = (32, 32, 16)
halo size (Hx, Hy, Hz) = (1, 1, 1)
domain (Lx, Ly, Lz) = (8.0f0, 8.0f0, 2.0f0)
grid spacing (Δx, Δy, Δz) = (0.25f0, 0.25f0, 0.125f0)Models
Oceananigans.Model — Method.Model(; grid, kwargs...)Construct an Oceananigans.jl model on grid.
Keyword arguments
grid: (required) The resolution and discrete geometry on whichmodelis solved. Currently the only option isRegularCartesianGrid.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.BasicModel — Method.BasicModel(; N, L, ν=ν₀, κ=κ₀, float_type=Float64, kwargs...)Construct a "Basic" Model with resolution N, domain extent L, precision float_type, and constant isotropic viscosity and diffusivity ν, and κ.
Additional kwargs are passed to the regular Model constructor.
Oceananigans.ChannelModel — Method.ChannelModel(; 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.
Output writers
Oceananigans.Checkpointer — Type.Checkpointer{I, T, P, A} <: AbstractOutputWriterAn output writer for checkpointing models to a JLD2 file from which models can be restored.
Oceananigans.Checkpointer — Method.Checkpointer(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.FieldOutput — Type.FieldOutput([return_type=Array], field)Returns a FieldOutput type intended for use with the JLD2OutputWriter. Calling FieldOutput(model) returns return_type(field.data.parent).
Oceananigans.JLD2OutputWriter — Type.JLD2OutputWriter{F, I, O, IF, IN, KW} <: AbstractOutputWriterAn output writer for writing to JLD2 files.
Oceananigans.JLD2OutputWriter — Method.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)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.NetCDFOutputWriter — Type.NetCDFOutputWriter <: AbstractOutputWriterAn output writer for writing to NetCDF files.
Oceananigans.FieldOutputs — Method.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
:uOceananigans.restore_from_checkpoint — Method.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.
Oceananigans.write_grid — Method.write_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())Time steppers
Oceananigans.time_step! — Method.time_step!(model, Nt, Δt; init_with_euler=true)Step forward model Nt time steps with step size Δt.
If init_with_euler is set to true, then the first step will be taken using a first-order forward Euler method.
Utilities
Oceananigans.pretty_filesize — Function.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).
Oceananigans.prettytime — Method.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).
Oceananigans.update_Δt! — Method.update_Δt!(wizard, model)Compute wizard.Δt given the velocities and diffusivities of model, and the parameters of wizard.