MOM6
mom_eos_linear::calculate_spec_vol_linear Interface Reference

Detailed Description

Compute the specific volume of sea water (in m^3/kg), or its anomaly from a reference value, using a simple linear equation of state from salinity (in psu), potential temperature (in deg C) and pressure [Pa].

Definition at line 34 of file MOM_EOS_linear.F90.

Private functions

subroutine calculate_spec_vol_scalar_linear (T, S, pressure, specvol, Rho_T0_S0, dRho_dT, dRho_dS, spv_ref)
 This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC]) and pressure [Pa], using a trivial linear equation of state for density. If spv_ref is present, specvol is an anomaly from spv_ref. More...
 
subroutine calculate_spec_vol_array_linear (T, S, pressure, specvol, start, npts, Rho_T0_S0, dRho_dT, dRho_dS, spv_ref)
 This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC]) and pressure [Pa], using a trivial linear equation of state for density. If spv_ref is present, specvol is an anomaly from spv_ref. More...
 

Functions and subroutines

◆ calculate_spec_vol_array_linear()

subroutine mom_eos_linear::calculate_spec_vol_linear::calculate_spec_vol_array_linear ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(out)  specvol,
integer, intent(in)  start,
integer, intent(in)  npts,
real, intent(in)  Rho_T0_S0,
real, intent(in)  dRho_dT,
real, intent(in)  dRho_dS,
real, intent(in), optional  spv_ref 
)
private

This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC]) and pressure [Pa], using a trivial linear equation of state for density. If spv_ref is present, specvol is an anomaly from spv_ref.

Parameters
[in]tpotential temperature relative to the surface [degC].
[in]sSalinity [PSU].
[in]pressurePressure [Pa].
[out]specvolin situ specific volume [m3 kg-1].
[in]startthe starting point in the arrays.
[in]nptsthe number of values to calculate.
[in]rho_t0_s0The density at T=0, S=0 [kg m-3].
[in]drho_dtThe derivatives of density with temperature [kg m-3 degC-1].
[in]drho_dsThe derivatives of density with salinity [kg m-3 ppt-1].
[in]spv_refA reference specific volume [m3 kg-1].

Definition at line 138 of file MOM_EOS_linear.F90.

138  real, dimension(:), intent(in) :: T !< potential temperature relative to the surface
139  !! [degC].
140  real, dimension(:), intent(in) :: S !< Salinity [PSU].
141  real, dimension(:), intent(in) :: pressure !< Pressure [Pa].
142  real, dimension(:), intent(out) :: specvol !< in situ specific volume [m3 kg-1].
143  integer, intent(in) :: start !< the starting point in the arrays.
144  integer, intent(in) :: npts !< the number of values to calculate.
145  real, intent(in) :: Rho_T0_S0 !< The density at T=0, S=0 [kg m-3].
146  real, intent(in) :: dRho_dT !< The derivatives of density with temperature [kg m-3 degC-1].
147  real, intent(in) :: dRho_dS !< The derivatives of density with salinity [kg m-3 ppt-1].
148  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1].
149  ! Local variables
150  integer :: j
151 
152  if (present(spv_ref)) then ; do j=start,start+npts-1
153  specvol(j) = ((1.0 - rho_t0_s0*spv_ref) + spv_ref*(drho_dt*t(j) + drho_ds*s(j))) / &
154  ( rho_t0_s0 + (drho_dt*t(j) + drho_ds*s(j)))
155  enddo ; else ; do j=start,start+npts-1
156  specvol(j) = 1.0 / ( rho_t0_s0 + (drho_dt*t(j) + drho_ds*s(j)))
157  enddo ; endif
158 

◆ calculate_spec_vol_scalar_linear()

subroutine mom_eos_linear::calculate_spec_vol_linear::calculate_spec_vol_scalar_linear ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  specvol,
real, intent(in)  Rho_T0_S0,
real, intent(in)  dRho_dT,
real, intent(in)  dRho_dS,
real, intent(in), optional  spv_ref 
)
private

This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC]) and pressure [Pa], using a trivial linear equation of state for density. If spv_ref is present, specvol is an anomaly from spv_ref.

Parameters
[in]tpotential temperature relative to the surface [degC].
[in]sSalinity [PSU].
[in]pressurePressure [Pa].
[out]specvolIn situ specific volume [m3 kg-1].
[in]rho_t0_s0The density at T=0, S=0 [kg m-3].
[in]drho_dtThe derivatives of density with temperature [kg m-3 degC-1].
[in]drho_dsThe derivatives of density with salinity [kg m-3 ppt-1].
[in]spv_refA reference specific volume [m3 kg-1].

Definition at line 111 of file MOM_EOS_linear.F90.

111  real, intent(in) :: T !< potential temperature relative to the surface
112  !! [degC].
113  real, intent(in) :: S !< Salinity [PSU].
114  real, intent(in) :: pressure !< Pressure [Pa].
115  real, intent(out) :: specvol !< In situ specific volume [m3 kg-1].
116  real, intent(in) :: Rho_T0_S0 !< The density at T=0, S=0 [kg m-3].
117  real, intent(in) :: dRho_dT !< The derivatives of density with temperature [kg m-3 degC-1].
118  real, intent(in) :: dRho_dS !< The derivatives of density with salinity [kg m-3 ppt-1].
119  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1].
120  ! Local variables
121  integer :: j
122 
123  if (present(spv_ref)) then
124  specvol = ((1.0 - rho_t0_s0*spv_ref) + spv_ref*(drho_dt*t + drho_ds*s)) / &
125  ( rho_t0_s0 + (drho_dt*t + drho_ds*s))
126  else
127  specvol = 1.0 / ( rho_t0_s0 + (drho_dt*t + drho_ds*s))
128  endif
129 

The documentation for this interface was generated from the following file: