List of functions in WenoNeverworld
WenoNeverworld.checkpoint_outputs! — Methodfunction checkpoint_outputs!(simulation, output_prefix; overwrite_existing = true, checkpoint_time = 100days)attaches a Checkpointer to the simulation with prefix output_prefix that is saved every checkpoint_time	
WenoNeverworld.initialize_model! — Methodfunction initialize_model!(model, Val(interpolate), initial_buoyancy, grid, previous_grid, init_file, buoyancymodel)initializes the model according to interpolate or not on a finer/coarser grid Val(interpolate)
WenoNeverworld.reduced_outputs! — Methodfunction reduced_outputs!(simulation, output_prefix; overwrite_existing = true, 	
                                                     checkpoint_time    = 100days,	
                                                     snapshot_time      = 30days,	
                                                     surface_time       = 1days,	
                                                     bottom_time        = 1days)attaches four JLD2OutputWriters to simulation with prefix output_prefix	
Outputs attached
snapshots: snapshots ofu,v,wandbsaved everysnapshot_timesurface_fields: snapshots ofu,v,wandbat the surface saved everysurface_timebottom_fields: snapshots ofu,v,wandbat the bottom (k = 2) saved everybottom_timecheckpointer: checkpointer saved everycheckpoint_time
WenoNeverworld.standard_outputs! — Methodfunction standard_outputs!(simulation, output_prefix; overwrite_existing = true, 	
                                                      checkpoint_time    = 100days,	
                                                      snapshot_time      = 30days,	
                                                      surface_time       = 5days,	
                                                      average_time       = 30days,	
                                                      average_window     = average_time,	    
                                                      average_stride     = 10)attaches four JLD2OutputWriters to simulation with prefix output_prefix	
Outputs attached
snapshots: snapshots ofu,v,wandbsaved everysnapshot_timesurface_fields: snapshots ofu,v,wandbat the surface saved everysurface_timeaveraged_fields: averages ofu,v,w,b,ζ,ζ2,u2,v2,w2,b2,ub,vb, andwbsaved everyaverage_timewith a window ofaverage_windowand stride ofaverage_stridecheckpointer: checkpointer saved everycheckpoint_time
WenoNeverworld.vertically_averaged_outputs! — Methodfunction vertically_averaged_outputs!(simulation, output_prefix; overwrite_existing = true, 	
                                                                 average_time       = 30days,	
                                                                 average_window     = average_time,	    
                                                                 average_stride     = 10)attaches a JLD2OutputWriters to simulation with prefix output_prefix	
Outputs attached
vertically_averaged_outputs: average ofKEand heat content (integral of temperature in ᵒCm³)
WenoNeverworld.weno_neverworld_simulation — Methodfunction weno_neverworld_simulation(grid; 
                                    previous_grid = grid,
                                    μ_drag = 0.001,  
                                    convective_adjustment = default_convective_adjustment,
                                    vertical_diffusivity  = default_vertical_diffusivity,
                                    horizontal_closure    = nothing,
                                    coriolis = HydrostaticSphericalCoriolis(scheme = ActiveCellEnstrophyConserving()),
                                    free_surface = SplitExplicitFreeSurface(; grid, cfl = 0.75),
                                    momentum_advection = default_momentum_advection(grid.underlying_grid),
                                    tracer_advection   = WENO(grid.underlying_grid), 
                                    interp_init = false,
                                    init_file = nothing,
                                    Δt = 5minutes,
                                    stop_time = 10years,
                                    stop_iteration = Inf,
                                    initial_buoyancy = initial_buoyancy_parabola,
                                    wind_stress               = WindStressBoundaryCondition(),
                                    buoyancy_relaxation       = BuoyancyRelaxationBoundaryCondition(),
                                    tracer_boundary_condition = NamedTuple(),
                                    tracers = :b
                                    )returns a simulation object for the Neverworld simulation.
Arguments:
grid: the grid on which the simulation is to be run
Keyword arguments:
previous_grid: the grid on whichinit_filehas been generated, if we restart frominit_fileμ_drag: the drag coefficient for the quadratic bottom drag, default:0.001convective_adjustment: the convective adjustment scheme, default:RiBasedVerticalDiffusivity()vertical_diffusivity: the vertical diffusivity scheme, default:VerticalScalarDiffusivity(ν=1e-4, κ=3e-5)horizontal_closure: the horizontal closure scheme, default:nothingcoriolis: the coriolis scheme, default:HydrostaticSphericalCoriolis(scheme = ActiveCellEnstrophyConserving())free_surface: the free surface scheme, default: SplitExplicitFreeSurface(; grid, cfl = 0.75)momentum_advection: the momentum advection scheme, default:VectorInvariant(vorticity_scheme = WENO(order = 9), vertical_scheme = WENO(grid))tracer_advection: the tracer advection scheme, default:WENO(grid)interp_init: whether to interpolate the initial conditions frominit_filetogrid, default: falseinit_file: the file from which to read the initial conditions, default:nothingΔt: the time step, default:5minutesstop_time: the time at which to stop the simulation, default:10yearsstop_iteration: the iteration at which to stop the simulation, default: Infinitial_buoyancy: the initial buoyancy field in case ofinit_file = nothing, function of(x, y, z)default:initial_buoyancy_parabolawind_stress: the wind stress boundary condition, default:WindStressBoundaryCondition()(seesrc/neverworld_initial_and_boundary_conditions.jl)buoyancy_relaxation: the buoyancy relaxation boundary condition, default:BuoyancyRelaxationBoundaryCondition()(seesrc/neverworld_initial_and_boundary_conditions.jl)tracer_boundary_condition: boundary conditions for tracers outside:b, default: nothingtracers: the tracers to be advected, default::b
WenoNeverworld.Diagnostics.AreaField — FunctionAreaField(grid, loc=(Center, Center, Nothing); indices = default_indices(3))Returns a two-dimensional field containing the cell horizontal areas at location loc with indices indices.
WenoNeverworld.Diagnostics.DeformationRadius — MethodDeformationRadius(f::Dict, i)Returns the two-dimensional deformation vorticity at time index i.
WenoNeverworld.Diagnostics.DensityField — MethodDensityField(b::Field; ρ₀ = 1000.0, g = 9.80655)Returns the three-dimensional density given a buoyancy field b.
WenoNeverworld.Diagnostics.HeightField — FunctionHeightField(grid, loc = (Center, Center, Center))Returns a three-dimensional field containing the cell vertical spacing at location loc.
WenoNeverworld.Diagnostics.KineticEnergy — MethodKineticEnergy(f::Dict, i)Returns the three-dimensional kinetic energy at time index i.
WenoNeverworld.Diagnostics.PotentialVorticity — MethodPotentialVorticity(f::Dict, i)Returns the three-dimensional potential vorticity at time index i.
WenoNeverworld.Diagnostics.Stratification — MethodStratification(f::Dict, i)Returns the three-dimensional stratification at time index i.
WenoNeverworld.Diagnostics.VerticalVorticity — MethodVerticalVorticity(f::Dict, i)Returns the three-dimensional vertical vorticity at time index i.
WenoNeverworld.Diagnostics.VolumeField — FunctionVolumeField(grid, loc=(Center, Center, Center);  indices = default_indices(3))Returns a three-dimensional field containing the cell volumes at location loc with indices indices.
WenoNeverworld.Diagnostics.add_kinetic_energy_and_vorticity_to_timeseries! — Methodadds the kinetic energy and vertical vorticity to an instantaneous timeseries
WenoNeverworld.Diagnostics.add_kinetic_energy_to_averaged_timeseries! — Methodadds the kinetic energy to a timeseries of values averaged in time. The latter must contains u2 and v2
WenoNeverworld.Diagnostics.all_fieldtimeseries — Functionall_fieldtimeseries(filename, dir = nothing; variables = ("u", "v", "w", "b"), checkpointer = false, number_files = nothing)Load and return a dictionary of field time series data.
Arguments
filename: The name of the file containing the field data.dir: The directory where the field data files are located. Defaults to "./".variables: A tuple of variable names to load. Defaults to("u", "v", "w", "b").checkpointer: A boolean indicating whether to read checkpointers or time series. Defaults tofalse.number_files: The number of files to load. Defaults tonothing.
Returns
A dictionary where the keys are symbols representing the variable names and the values are FieldTimeSeries objects.
WenoNeverworld.Diagnostics.checkpoint_fields — Methodreturns a nametuple of (u, v, w, b) from the data in file
WenoNeverworld.Diagnostics.compress_all_restarts — Methodcompress_all_restarts(resolution, ranks, dir; output_prefix = "weno_thirtytwo", remove_restart = false, leave_last_file = true)Compresses all restart files in the specified directory.
Arguments
resolution: The resolution of the restart files.ranks: The number of ranks used for the simulations.dir: The directory containing the restart files.
Keyword Arguments
output_prefix: The prefix for the compressed files. Default is "weno_thirtytwo".remove_restart: Whether to remove the original restart files after compression. Default isfalse.leave_last_file: Whether to leave the last file uncompressed. Default istrue.
WenoNeverworld.Diagnostics.compress_restart_file — Functioncompress_restart_file(resolution, ranks, iteration, folder = "../")Compresses the restart files for a given simulation.
Arguments
resolution: The resolution of the simulation.ranks: The number of ranks used in the simulation.iteration: The iteration of the checkpoint.folder: The folder where the restart files are located. Default is"../".
Examples
julia> compress_restart_file(1/32, 8, 0)WenoNeverworld.Diagnostics.integral_available_potential_energy — Methodintegral_available_potential_energy(b::FieldTimeSeries; stride = 1, start_time = 1, end_time = length(u.times))Compute the integral available potential energy (APE) over time for a given FieldTimeSeries b.
Arguments
b::FieldTimeSeries: The field time series containing bouyancy data.stride::Int: The stride value for iterating over the time steps. Default is 1.start_time::Int: The starting time step for integration. Default is 1.end_time::Int: The ending time step for integration. Default is the length ofu.times.
Returns
energy::Vector{Float64}: The vector of integrated APE values over time.
WenoNeverworld.Diagnostics.integral_kinetic_energy — Methodintegral_kinetic_energy(u::FieldTimeSeries, v::FieldTimeSeries; stride = 1, start_time = 1, end_time = length(u.times))Compute the integral kinetic energy over time for the given field time series u and v.
Arguments
u::FieldTimeSeries: The field time series for the u-component of the velocity.v::FieldTimeSeries: The field time series for the v-component of the velocity.stride::Int: The stride between time steps to consider. Default is 1.start_time::Int: The starting time step to consider. Default is 1.end_time::Int: The ending time step to consider. Default is the length ofu.times.
Returns
energy::Vector{Float64}: The computed integral of kinetic energy over time.
WenoNeverworld.Diagnostics.limit_timeseries! — Methodlimit the timeseries to times
WenoNeverworld.Diagnostics.limit_timeseries! — Methodlimit the timeseries to times
WenoNeverworld.Diagnostics.propagate — Methodpropagate(fields...; func)Propagates the function func with inputs fields... through time.
Arguments
fields: The input fieldsfunc: The function to apply to the fields at each time step.
Returns
field_output: The output of functionfuncas aFieldTimeSeriesobject.
WenoNeverworld.Diagnostics.reduce_output_size! — Methodsaves a new file with name new_file_name with the last limit_to timeseries
WenoNeverworld.Auxiliaries.continue_downwards! — Methodcontinue_downwards!(field)continue downwards a field with missing values at field[i, j, k] == 0
the continue_downwards! implementation is inspired by https://github.com/CliMA/ClimaOcean.jl/pull/60
WenoNeverworld.Auxiliaries.cubic_interpolate — Methodfunction cubic_interpolate(x, x1, x2, y1, y2, d1, d2)returns a cubic function between points (x1, y1) and (x2, y2) with derivative d1 and d2
WenoNeverworld.Auxiliaries.exponential_profile — Methodutility profiles (atan, exponential, and parabolic)
WenoNeverworld.Auxiliaries.increase_simulation_Δt! — Methodfunction increase_simulation_Δt!(simulation; cutoff_time = 20days, new_Δt = 2minutes)utility to update the Δt of a simulation after a certain cutoff_time with new_Δt.	 Note: this function adds a callback to simulation, so the order of increase_simulation_Δt! 	 matters (i.e. the Δt will be updated based on the order of increase_simulation_Δt! specified)	
WenoNeverworld.Auxiliaries.propagate_horizontally! — Methodpropagate_horizontally!(field; max_iter = Inf)propagate horizontally a field with missing values at field[i, j, k] == 0
disclaimer: the propagate_horizontally! implementation is inspired by https://github.com/CliMA/ClimaOcean.jl/pull/60
WenoNeverworld.Auxiliaries.regrid_field — Methodfunction regridded_field(old_vector, old_grid, new_grid, loc)interpolate old_vector (living on loc) from old_grid to new_grid 	
WenoNeverworld.Auxiliaries.update_simulation_clock! — Methodfunction update_simulation_clock!(simulation, init_file)updates the clock of simulation with the time in init_file	
WenoNeverworld.NeverworldGrids.NeverworldGrid — Functionfunction NeverworldGrid(arch, degree, FT::DataType = Float64; H = 7, longitude = (-2, 62), latitude = (-70, 0), bathymetry_params = NeverWorldBathymetryParameters(), longitudinal_extent = 60)builds a LatitudeLongitudeGrid with a specified bathymetry
Arguments
arch: architecture of the grid, can beCPU()orGPU()orDistributedresolution: resolution in degrees.FT: (optional) floating point precision (default =Float64)
Keyword Arguments
H: halo size,Intlongitudinal_extent: size of the actual domain in longitudinal direction,Numberlongitude: longitudinal extremes of the domain,Tuple. Note: this keyword must be at leastlongitude_extent + resolution * 2Hto allow for correct advection stencilslatitude: latitudinal extremes of the domainbathymetry_params: parameters for the neverworld bathymetry, seeneverworld_bathymetry.jlz_faces: array containing the z faces
WenoNeverworld.NeverworldGrids.bottom_ridge_xy — Methodsmoothed coasts for the inlet and outlet of the channel
WenoNeverworld.NeverworldGrids.exponential_z_faces — Methodfunction exponential_z_faces(; Nz = 69, Lz = 4000.0, e_folding = 0.06704463421863584)generates an array of exponential z faces
WenoNeverworld.NeverworldBoundaries.BuoyancyRelaxationBoundaryCondition — TypeBuoyancyRelaxationBoundaryCondition(func = (y, t) -> parabolic_scaling(y); ΔB = ΔB, λ = 7days)Buoyancy relaxation profile which implements a latitude-time dependent boundary condition following:
b = Δz_surface / λ * (b_surf - ΔB * func(φ, t))
Arguments:
- func: function which takes the latitude φ and time t and returns a scalar
 
Keyword arguments:
- ΔB: buoyancy difference between the equator and the poles, default: 6.0e-2
 - λ: restoring time-scale, default: 7days
 
WenoNeverworld.NeverworldBoundaries.WindStressBoundaryCondition — MethodWindStressBoundaryCondition(; φs = default_φs, τs = default_τs)Wind stess boundary condition which implements a piecewise cubic interpolation between points φs (Tuple) and τs (Tuple).
WenoNeverworld.Parameterizations.EnergyBackScattering — Typestruct EnergyBackScattering{FT} <: AbstractTurbulenceClosure{ExplicitTimeDiscretization, 3}Energy backscattering turbulence closure model. This struct represents a turbulence closure model based on the energy backscattering principle.  It is a parameterization of the turbulent momentum flux in a fluid flow.  The model is implemented as a struct with a type parameter FT representing the floating-point type used for calculations.
Arguments
ν::FT: The kinematic anti-viscosity of the fluid.
reference: Zanna, L., Bolton, T. (2020). Data-driven equation discovery of ocean mesoscale closures. Geophysical Research Letters, 47, e2020GL088376. https://doi.org/10.1029/2020GL088376
WenoNeverworld.Parameterizations.QGLeith — Typestruct QGLeith{FT, M, S} <: AbstractScalarDiffusivity{ExplicitTimeDiscretization, HorizontalFormulation, 2}QGLeith is a struct representing the Leith scalar diffusivity parameterization for quasi-geostrophic models.
Fields
C: The coefficient for the diffusivity parameterization.min_N²: The minimum value for the squared buoyancy frequency.isopycnal_tensor: The isopycnal tensor model used for the diffusivity calculation.slope_limiter: The slope limiter used for the diffusivity calculation.
Constructors
QGLeith(FT::DataType = Float64; C=FT(1.0), min_N² = FT(1e-20), isopycnal_model=SmallSlopeIsopycnalTensor(), slope_limiter=FluxTapering(1e-2)): Construct a QGLeith object with optional parameters.
WenoNeverworld.Parameterizations.Δ²ᶜᶜᶜ — MethodReturn the filter width for an Horizontal closure on a general grid.