Function index

OceanRasterConversions exported functions

OceanRasterConversions.OceanRasterConversionsModule

Module to convert variables depth, practical salinity and potential temperature to the TEOS-10 standard variables pressure, absolute salinity and conservative temperature (respectively) from a Raster, RasterStack or RasterSeries. A few chosen seawater variables can then be computed from these state variables.

source
OceanRasterConversions.Sₚ_to_SₐMethod
function Sₚ_to_Sₐ(Sₚ::Raster)
function Sₚ_to_Sₐ(stack::RasterStack, Sₚ::Symbol)
function Sₚ_to_Sₐ(series::RasterSeries, Sₚ::Symbol)

Convert a Raster of practical salinity (Sₚ) to absolute salinity (Sₐ) using gsw_sa_from_sp from GibbsSeaWater.jl. This conversion depends on pressure. If converting from a RasterStack or RasterSeries, the symbol for the practical salinity in the RasterStack/Series must be passed in as a Symbol –- that is if the variable name is SALT the RasterStack/Series, the Symbol :SALT must be passed in.

source
OceanRasterConversions.convert_ocean_varsMethod
function convert_ocean_vars(raster::RasterStack, var_names::NamedTuple;
                            ref_pressure = nothing,
                            with_α = false,
                            with_β = false)
function convert_ocean_vars(raster::Rasterseries, var_names::NamedTuple;
                            ref_pressure = nothing,
                            with_α = false,
                            with_β = false)

Convert ocean variables depth, practical salinity and potential temperature to pressure, absolute salinity and conservative temperature. All conversions are done using the julia implementation of TEOS-10 GibbsSeaWater.jl. A new Raster is returned that contains the variables pressure, absolute salinity, conservative temperature and density (either in-situ or referenced to a user defined reference pressure). As pressure depends on latitude and depth, it is added as a new variable –- that is, each longitude, latitude, depth and time have a variable for pressure. A density variable is also computed which, by default, is in-situ density. Potential density at a reference pressure can be computed instead by passing a the keyword argument ref_pressure. Optional keyword arguments with_α and with_β allow the thermal expansion and haline contraction coefficients (respectively) to be computed and added to the returned RasterStack/Series.

The name of the variables for potential temperature and practical salinity must be passed in as a NamedTuple of the form (Sₚ = :salt_name, θ = :potential_temp_name) where :potential_temp_name and :salt_name are the name of the potential temperature and practical salinity in the Raster.

source
OceanRasterConversions.depth_to_pressureMethod
function depth_to_pressure(raster::Raster)
function depth_to_pressure(stack::RasterStack)

Convert the depth dimension (Z) to pressure using gsw_p_from_z from GibbsSeaWater.jl. Note that pressure depends on depth and latitude so the returned pressure is stored as a variable in the resulting Raster rather than replacing the vertical depth dimension.

source
OceanRasterConversions.get_αMethod
function get_α(Sₐ::Raster, Θ::Raster, p::Raster)
function get_α(stack::RasterStack, var_names::NamedTuple)
function get_α(series::RasterSeries, var_names::NamedTuple)

Compute the thermal exapnsion coefficient, α, using gsw_alpha from GibbsSeaWater.jl. To compute α from a RasterStack or RasterSeries the variable names must be passed into the function as a NamedTuple in the form (Sₐ = :salt_var, Θ = :temp_var, p = :pressure_var). The returned Raster will have the same dimensions as Rasterstack that is passed in.

source
OceanRasterConversions.get_βMethod
function get_β(Sₐ::Raster, Θ::Raster, p::Raster)
function get_β(stack::RasterStack, var_names::NamedTuple)
function get_β(series::RasterSeries, var_names::NamedTuple)

Compute the haline contraction coefficient, β, using gsw_beta from GibbsSeaWater.jl. To compute β from a RasterStack or RasterSeries the variable names must be passed into the function as a NamedTuple in the form (Sₐ = :salt_var, Θ = :temp_var, p = :pressure_var). The returned Raster will have the same dimensions as Rasterstack that is passed in.

source
OceanRasterConversions.get_ρMethod
function get_ρ(Sₐ::Raster, Θ::Raster, p::Raster)
function get_ρ(stack::RasterStack, var_names::NamedTuple)
function get_ρ(series::RasterStack, var_names::NamedTuple)

Compute in-situ density, ρ, using gsw_rho from GibbsSeaWater.jl. This computation depends on absolute salinity (Sₐ), conservative temperature (Θ) and pressure (p). To compute ρ from a RasterStack or RasterSeries the variable names must be passed into the function as a NamedTuple in the form (Sₐ = :salt_var, Θ = :temp_var, p = :pressure_var). The returned Raster will have the same dimensions as Rasterstack that is passed in.

source
OceanRasterConversions.get_σₚMethod
function get_σₚ(Sₐ::Raster, Θ::Raster, p::Number)
function get_σₚ(stack::RasterStack, var_names::NamedTuple)
function get_σₚ(series::RasterStack, var_names::NamedTuple)

Compute potential density at reference pressure p, σₚ, using gsw_rho from GibbsSeaWater.jl. This computation depends on absolute salinity (Sₐ), conservative temperature (Θ) and a user entered reference pressure (p). Compute and return the potential density σₚ at reference pressure p from a RasterStack or RasterSeries. This computation depends on absolute salinity Sₐ, conservative temperature Θ and a reference pressure p. The variable names must be passed into the function as a NamedTuple in the form (Sₐ = :salt_var, Θ = :temp_var, p = ref_pressure). Note p in this case is a number. The returned Raster will have the same dimensions as Rasterstack that is passed in.

source
OceanRasterConversions.θ_to_ΘMethod
function θ_to_Θ(θ::Raster, Sₐ::Raster)
function θ_to_Θ(stack::RasterStack, var_names::NamedTuple)
function θ_to_Θ(series::RasterSeries, var_names::NamedTuple)

Convert a Raster of potential temperature (θ) to conservative temperature (Θ) using gsw_ct_from_pt from GibbsSeaWater.jl. This conversion depends on absolute salinity. If converting from a from a RasterStack or a RasterSeries, the var_names must be passed in as for convert_ocean_vars –- that is, as a named tuple in the form (Sₚ = :salt_name, θ = :potential_temp_name) where :potential_temp_name and :salt_name are the name of the potential temperature and salinity in the RasterStack.

source

OceanRasterConversions private functions