namespace mom_eos_wright

Overview

The equation of state using the Wright 1997 expressions. More…

namespace mom_eos_wright {

// interfaces

interface calculate_density_derivs_wright;
interface calculate_density_second_derivs_wright;
interface calculate_density_wright;
interface calculate_spec_vol_wright;

// global functions

subroutine, public calculate_specvol_derivs_wright(
    T T,
    S S,
    pressure pressure,
    dSV_dT dSV_dT,
    dSV_dS dSV_dS,
    start start,
    npts npts
    );

subroutine, public calculate_compress_wright(
    T T,
    S S,
    pressure pressure,
    rho rho,
    drho_dp drho_dp,
    start start,
    npts npts
    );

subroutine, public int_density_dz_wright(
    T T,
    S S,
    z_t z_t,
    z_b z_b,
    rho_ref rho_ref,
    rho_0 rho_0,
    G_e G_e,
    HII HII,
    HIO HIO,
    dpa dpa,
    intz_dpa intz_dpa,
    intx_dpa intx_dpa,
    inty_dpa inty_dpa,
    bathyT bathyT,
    dz_neglect dz_neglect,
    useMassWghtInterp useMassWghtInterp
    );

subroutine, public int_spec_vol_dp_wright(
    T T,
    S S,
    p_t p_t,
    p_b p_b,
    spv_ref spv_ref,
    HI HI,
    dza dza,
    intp_dza intp_dza,
    intx_dza intx_dza,
    inty_dza inty_dza,
    halo_size halo_size,
    bathyP bathyP,
    dP_neglect dP_neglect,
    useMassWghtInterp useMassWghtInterp
    );

} // namespace mom_eos_wright

Detailed Documentation

The equation of state using the Wright 1997 expressions.

Global Functions

subroutine, public calculate_specvol_derivs_wright(
    T T,
    S S,
    pressure pressure,
    dSV_dT dSV_dT,
    dSV_dS dSV_dS,
    start start,
    npts npts
    )

For a given thermodynamic state, return the partial derivatives of specific volume with temperature and salinity.

Parameters:

t

Potential temperature relative to the surface [degC].

s

Salinity [PSU].

pressure

pressure [Pa].

dsv_dt

The partial derivative of specific volume with potential temperature [m3 kg-1 degC-1].

dsv_ds

The partial derivative of specific volume with salinity [m3 kg-1 / Pa].

start

The starting point in the arrays.

npts

The number of values to calculate.

subroutine, public calculate_compress_wright(
    T T,
    S S,
    pressure pressure,
    rho rho,
    drho_dp drho_dp,
    start start,
    npts npts
    )

This subroutine computes the in situ density of sea water (rho in [kg m-3]) and the compressibility (drho/dp = C_sound^-2) (drho_dp [s2 m-2]) from salinity (sal in psu), potential temperature (T [degC]), and pressure [Pa]. It uses the expressions from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740. Coded by R. Hallberg, 1/01.

Parameters:

t

Potential temperature relative to the surface [degC].

s

Salinity [PSU].

pressure

pressure [Pa].

rho

In situ density [kg m-3].

drho_dp

The partial derivative of density with pressure (also the inverse of the square of sound speed) [s2 m-2].

start

The starting point in the arrays.

npts

The number of values to calculate.

subroutine, public int_density_dz_wright(
    T T,
    S S,
    z_t z_t,
    z_b z_b,
    rho_ref rho_ref,
    rho_0 rho_0,
    G_e G_e,
    HII HII,
    HIO HIO,
    dpa dpa,
    intz_dpa intz_dpa,
    intx_dpa intx_dpa,
    inty_dpa inty_dpa,
    bathyT bathyT,
    dz_neglect dz_neglect,
    useMassWghtInterp useMassWghtInterp
    )

This subroutine calculates analytical and nearly-analytical integrals of pressure anomalies across layers, which are required for calculating the finite-volume form pressure accelerations in a Boussinesq model.

Parameters:

hii

The horizontal index type for the input arrays.

hio

The horizontal index type for the output arrays.

t

Potential temperature relative to the surface

s

Salinity [PSU].

z_t

Height at the top of the layer in depth units [Z ~> m].

z_b

Height at the top of the layer [Z ~> m].

rho_ref

A mean density [kg m-3], that is subtracted out to reduce the magnitude of each of the integrals. (The pressure is calucated as p~=-z*rho_0*G_e.)

rho_0

Density [kg m-3], that is used to calculate the pressure (as p~=-z*rho_0*G_e) used in the equation of state.

g_e

The Earth’s gravitational acceleration [m2 Z-1 s-2 ~> m s-2].

dpa

The change in the pressure anomaly across the

intz_dpa

The integral through the thickness of the layer

intx_dpa

The integral in x of the difference between the

inty_dpa

The integral in y of the difference between the

bathyt

The depth of the bathymetry [Z ~> m].

dz_neglect

A miniscule thickness change [Z ~> m].

usemasswghtinterp

If true, uses mass weighting to interpolate T/S for top and bottom integrals.

subroutine, public int_spec_vol_dp_wright(
    T T,
    S S,
    p_t p_t,
    p_b p_b,
    spv_ref spv_ref,
    HI HI,
    dza dza,
    intp_dza intp_dza,
    intx_dza intx_dza,
    inty_dza inty_dza,
    halo_size halo_size,
    bathyP bathyP,
    dP_neglect dP_neglect,
    useMassWghtInterp useMassWghtInterp
    )

This subroutine calculates analytical and nearly-analytical integrals in pressure across layers of geopotential anomalies, which are required for calculating the finite-volume form pressure accelerations in a non-Boussinesq model. There are essentially no free assumptions, apart from the use of Bode’s rule to do the horizontal integrals, and from a truncation in the series for log(1-eps/1+eps) that assumes that |eps| < 0.34.

Parameters:

hi

The ocean’s horizontal index type.

t

Potential temperature relative to the surface

s

Salinity [PSU].

p_t

Pressure at the top of the layer [Pa].

p_b

Pressure at the top of the layer [Pa].

spv_ref

A mean specific volume that is subtracted out to reduce the magnitude of each of the integrals [m3 kg-1]. The calculation is mathematically identical with different values of spv_ref, but this reduces the effects of roundoff.

dza

The change in the geopotential anomaly across

intp_dza

The integral in pressure through the layer of

intx_dza

The integral in x of the difference between the

inty_dza

The integral in y of the difference between the

halo_size

The width of halo points on which to calculate dza.

bathyp

The pressure at the bathymetry [Pa]

dp_neglect

A miniscule pressure change with the same units as p_t [Pa]

usemasswghtinterp

If true, uses mass weighting to interpolate T/S for top and bottom integrals.