Additional Functionalities

Format Conversions

MITgcm.findtilesFunction
findtiles(ni::Int,nj::Int,mygrid::gcmgrid)
findtiles(ni::Int,nj::Int,grid::String="LatLonCap",GridParentDir="../inputs/GRID_LLC90/")

Return a MeshArray map of tile indices, mytiles["tileNo"], for tile size ni,nj and extract grid variables accordingly.

source
MITgcm.cube2compactFunction
cube2compact(x::Array)

Reshape from e.g. size (192, 32, 5) in cube format to (32, 192, 5) in compact format.

source
MITgcm.compact2cubeFunction
compact2cube(x::Array)

Reshape from e.g. size (32, 192, 5) in cube format to (192, 32, 5) in compact format.

source

Formulas, Parameters

MITgcm.SeaWaterDensityFunction

SeaWaterDensity(Θ,Σ,Π,Π0)

Compute potential density (ρP), in situ density (ρI), and density referenced to PREF (Π0 in decibars) from potential temperature (Θ in °C), salinity (Σ in psu) and pressure (Π in decibars) according to the UNESCO / Jackett & McDougall 1994 equation of state.

Credits: code based on a Matlab implementation by B. Ferron Reference: https://www.jodc.go.jp/info/iocdoc/UNESCOtech/059832eb.pdf Check value: ρI = 1041.83267kg/m^3 for Θ=3°Celcius, Σ=35psu, Π=3000dbar

(ρP,ρI,ρR) = SeaWaterDensity(3.,35.5,3000.)
isapprox(ρI,1041.83267, rtol=1e-6)
source
MITgcm.MixedLayerDepthFunction

MixedLayerDepth(Θ,Σ,Δ,mthd)

Compute mixed layer depth from potential temperature (Θ in °C), salinity (Σ in psu) and depth (Δ in method) according to various formulas (mthd == "BM", "Suga", "Kara"). Inputs must be dense vectors without any missing value (or NaN, etc).

D=collect(0.0:1.0:500.0); tmp=(1.0.-tanh.(5*(-1 .+ 2/D[end]*D)));
T=2.0 .+ 8.0*tmp; S=34.0 .+ 0.5*tmp;
(ρP,ρI,ρR) = SeaWaterDensity(T,S,D);

mld=MixedLayerDepth(T,S,D,"BM"); isapprox(mld,134.0)

using Plots
plot(ρP,-D,w=2,label="Potential Density",ylabel="Depth")
plot!(vec([ρP[1] ρP[end]]),-fill(mld,2),label="Mixed Layer Depth",w=2,c="black",s=:dash)
source