pymoc.modules.Psi_SO

class pymoc.modules.Psi_SO(z=None, y=None, b=None, bs=None, tau=None, f=0.00012, rho=1030, L=10000000.0, KGM=1000.0, c=None, bvp_with_Ek=False, Hsill=None, HEk=None, Htapertop=None, Htaperbot=None, smax=0.01)

Southern Ocean Overturning Transport Model

Instances of this class represent a 1D model of the overturning transport in the Southern Ocean interior, calculated based on a density profile in the adjoining basin, and local surface buoyancy and surface wind stress in the SO.

Parameters
  • z (ndarray) – Vertical depth levels of overturning grid. Units: m

  • y (ndarray) – Meridional overturning grid. Units: m

  • b (float, function, or ndarray) – Vertical buoyancy profile from the adjoining basin, on the north side of the ACC. Units: m/s2

  • bs (float, function, or ndarray) – Surface level buoyancy boundary condition. Can be a constant, or an array or function in y. Units: m/s2

  • tau (float, function, or ndarray) – Surface wind stress. Can be a constant, or an array or function in y. Units: N/m2

  • f (float) – Coriolis parameter. Units s-1

  • rho (float) – Density of sea water for Boussinesq approximation. Units: kg/m3

  • L (float) – Zonal length of the modeled ACC. Units: m

  • KGM (float) – Gent & McWilliams (GM) eddy diffusivity coefficient. Units:

  • c (float) – Phase speed cutoff for smoothing when solving the GM boundary value problem. Units: m/s

  • bvp_with_Ek (logical) – Whether to enforce the boundary condition that Psi_GM=-Psi_Ek at the ocean surface and bottom when solving the boundary value problem for the GM streamfunction.

  • Hsill (float) – Height above the bottom at which the Ekman streamfunction is tapered. Units: m

  • Hek (float) – Depth of the surface Ekman layer. Units: m

  • Htapertop (float) – Height of the quadratic surface tapering layer for the GM streamfunction. Units: m

  • Htaperbot (float) – Height of the quadratic bottom tapering layer for the GM streamfunction. Units: m

  • smax (float) – Maximum slope of the GM streamfunction, above which Psi_GM is clipped. Units: m-1

__init__(z=None, y=None, b=None, bs=None, tau=None, f=0.00012, rho=1030, L=10000000.0, KGM=1000.0, c=None, bvp_with_Ek=False, Hsill=None, HEk=None, Htapertop=None, Htaperbot=None, smax=0.01)
Parameters
  • z (ndarray) – Vertical depth levels of overturning grid. Units: m

  • y (ndarray) – Meridional overturning grid. Units: m

  • b (float, function, or ndarray) – Vertical buoyancy profile from the adjoining basin, on the north side of the ACC. Units: m/s2

  • bs (float, function, or ndarray) – Surface level buoyancy boundary condition. Can be a constant, or an array or function in y. Units: m/s2

  • tau (float, function, or ndarray) – Surface wind stress. Can be a constant, or an array or function in y. Units: N/m2

  • f (float) – Coriolis parameter. Units s-1

  • rho (float) – Density of sea water for Boussinesq approximation. Units: kg/m3

  • L (float) – Zonal length of the modeled ACC. Units: m

  • KGM (float) – Gent & McWilliams (GM) eddy diffusivity coefficient. Units:

  • c (float) – Phase speed cutoff for smoothing when solving the GM boundary value problem. Units: m/s

  • bvp_with_Ek (logical) – Whether to enforce the boundary condition that Psi_GM=-Psi_Ek at the ocean surface and bottom when solving the boundary value problem for the GM streamfunction.

  • Hsill (float) – Height above the bottom at which the Ekman streamfunction is tapered. Units: m

  • Hek (float) – Depth of the surface Ekman layer. Units: m

  • Htapertop (float) – Height of the quadratic surface tapering layer for the GM streamfunction. Units: m

  • Htaperbot (float) – Height of the quadratic bottom tapering layer for the GM streamfunction. Units: m

  • smax (float) – Maximum slope of the GM streamfunction, above which Psi_GM is clipped. Units: m-1

Methods

__init__([z, y, b, bs, tau, f, rho, L, KGM, …])

param z

Vertical depth levels of overturning grid. Units: m

bc_GM(ya, yb)

Calculate the residuals of boundary conditions for the eddy driven transport boundary value problem.

calc_Ekman()

Compute the Ekman transport from the wind stress averaged from the northern boundary of the domain to the latitude of the northernmost outcropped isopycnal.

calc_GM()

Compute the eddy (Gent & Mcwilliams) transport based on the meridionally averaged isopycnal slope.

calc_N2()

Calculate the buouyancy (Brunt-Väisällä) frequency profile for the Southern Ocean

calc_bottom_taper(H, z)

Calculate the quadratic tapering profile relative to the ocean floor.

calc_top_taper(H, z[, scalar])

Calculate the quadratic tapering profile relative to the ocean surface.

solve()

Compute the residual overturning transport in the Southern Ocean.

update([b, bs])

Update the vertical buoyancy profile and surface buoyancy, based on changes in the adjoining basin and/or in the surface boundary conditions.

ys(b)

Inversion function of \(bs\left(y\right)\).