MOM6
mom_eos::calculate_density Interface Reference

Detailed Description

Calculates density of sea water from T, S and P.

Definition at line 68 of file MOM_EOS.F90.

Private functions

subroutine calculate_density_scalar (T, S, pressure, rho, EOS, rho_ref, scale)
 Calls the appropriate subroutine to calculate density of sea water for scalar inputs. If rho_ref is present, the anomaly with respect to rho_ref is returned. The pressure and density can be rescaled with the US. If both the US and scale arguments are present the density scaling uses the product of the two scaling factors. More...
 
subroutine calculate_density_array (T, S, pressure, rho, start, npts, EOS, rho_ref, scale)
 Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs. If rho_ref is present, the anomaly with respect to rho_ref is returned. More...
 
subroutine calculate_density_1d (T, S, pressure, rho, EOS, dom, rho_ref, scale)
 Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs, potentially limiting the domain of indices that are worked on. If rho_ref is present, the anomaly with respect to rho_ref is returned. More...
 
subroutine calculate_stanley_density_scalar (T, S, pressure, Tvar, TScov, Svar, rho, EOS, rho_ref, scale)
 Calls the appropriate subroutine to calculate density of sea water for scalar inputs including the variance of T, S and covariance of T-S. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned. The density can be rescaled using rho_ref. More...
 
subroutine calculate_stanley_density_array (T, S, pressure, Tvar, TScov, Svar, rho, start, npts, EOS, rho_ref, scale)
 Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including the variance of T, S and covariance of T-S. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned. More...
 
subroutine calculate_stanley_density_1d (T, S, pressure, Tvar, TScov, Svar, rho, EOS, dom, rho_ref, scale)
 Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including the variance of T, S and covariance of T-S, potentially limiting the domain of indices that are worked on. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned. More...
 

Detailed Description

Calculates density of sea water from T, S and P.

Definition at line 68 of file MOM_EOS.F90.

Functions and subroutines

◆ calculate_density_1d()

subroutine mom_eos::calculate_density::calculate_density_1d ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(inout)  rho,
type(eos_type), pointer  EOS,
integer, dimension(2), intent(in), optional  dom,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs, potentially limiting the domain of indices that are worked on. If rho_ref is present, the anomaly with respect to rho_ref is returned.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [R L2 T-2 ~> Pa]
[in,out]rhoDensity (in-situ if pressure is local) [R ~> kg m-3]
eosEquation of state structure
[in]domThe domain of indices to work on, taking into account that arrays start at 1.
[in]rho_refA reference density [kg m-3]
[in]scaleA multiplicative factor by which to scale density in combination with scaling given by US [various]

Definition at line 359 of file MOM_EOS.F90.

360  real, dimension(:), intent(in) :: T !< Potential temperature referenced to the surface [degC]
361  real, dimension(:), intent(in) :: S !< Salinity [ppt]
362  real, dimension(:), intent(in) :: pressure !< Pressure [R L2 T-2 ~> Pa]
363  real, dimension(:), intent(inout) :: rho !< Density (in-situ if pressure is local) [R ~> kg m-3]
364  type(EOS_type), pointer :: EOS !< Equation of state structure
365  integer, dimension(2), optional, intent(in) :: dom !< The domain of indices to work on, taking
366  !! into account that arrays start at 1.
367  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3]
368  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density
369  !! in combination with scaling given by US [various]
370  ! Local variables
371  real :: p_scale ! A factor to convert pressure to units of Pa [Pa T2 R-1 L-2 ~> 1]
372  real :: rho_scale ! A factor to convert density from kg m-3 to the desired units [R m3 kg-1 ~> 1]
373  real :: rho_unscale ! A factor to convert density from R to kg m-3 [kg m-3 R-1 ~> 1]
374  real :: rho_reference ! rho_ref converted to [kg m-3]
375  real, dimension(size(rho)) :: pres ! Pressure converted to [Pa]
376  integer :: i, is, ie, npts
377 
378  if (.not.associated(eos)) call mom_error(fatal, &
379  "calculate_density_1d called with an unassociated EOS_type EOS.")
380 
381  if (present(dom)) then
382  is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
383  else
384  is = 1 ; ie = size(rho) ; npts = 1 + ie - is
385  endif
386 
387  p_scale = eos%RL2_T2_to_Pa
388  rho_unscale = eos%R_to_kg_m3
389 
390  if ((p_scale == 1.0) .and. (rho_unscale == 1.0)) then
391  call calculate_density_array(t, s, pressure, rho, is, npts, eos, rho_ref=rho_ref)
392  elseif (present(rho_ref)) then ! This is the same as above, but with some extra work to rescale variables.
393  do i=is,ie ; pres(i) = p_scale * pressure(i) ; enddo
394  rho_reference = rho_unscale*rho_ref
395  call calculate_density_array(t, s, pres, rho, is, npts, eos, rho_ref=rho_reference)
396  else ! There is rescaling of variables, but rho_ref is not present. Passing a 0 value of rho_ref
397  ! changes answers at roundoff for some equations of state, like Wright and UNESCO.
398  do i=is,ie ; pres(i) = p_scale * pressure(i) ; enddo
399  call calculate_density_array(t, s, pres, rho, is, npts, eos)
400  endif
401 
402  rho_scale = eos%kg_m3_to_R
403  if (present(scale)) rho_scale = rho_scale * scale
404  if (rho_scale /= 1.0) then ; do i=is,ie
405  rho(i) = rho_scale * rho(i)
406  enddo ; endif
407 

◆ calculate_density_array()

subroutine mom_eos::calculate_density::calculate_density_array ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(inout)  rho,
integer, intent(in)  start,
integer, intent(in)  npts,
type(eos_type), pointer  EOS,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs. If rho_ref is present, the anomaly with respect to rho_ref is returned.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [Pa] or [R L2 T-2 ~> Pa]
[in,out]rhoDensity (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
[in]startStart index for computation
[in]nptsNumber of point to compute
eosEquation of state structure
[in]rho_refA reference density [kg m-3]
[in]scaleA multiplicative factor by which to scale density in combination with scaling given by US [various]

Definition at line 262 of file MOM_EOS.F90.

263  real, dimension(:), intent(in) :: T !< Potential temperature referenced to the surface [degC]
264  real, dimension(:), intent(in) :: S !< Salinity [ppt]
265  real, dimension(:), intent(in) :: pressure !< Pressure [Pa] or [R L2 T-2 ~> Pa]
266  real, dimension(:), intent(inout) :: rho !< Density (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
267  integer, intent(in) :: start !< Start index for computation
268  integer, intent(in) :: npts !< Number of point to compute
269  type(EOS_type), pointer :: EOS !< Equation of state structure
270  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3]
271  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density
272  !! in combination with scaling given by US [various]
273  integer :: j
274 
275  if (.not.associated(eos)) call mom_error(fatal, &
276  "calculate_density_array called with an unassociated EOS_type EOS.")
277 
278  select case (eos%form_of_EOS)
279  case (eos_linear)
280  call calculate_density_linear(t, s, pressure, rho, start, npts, &
281  eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
282  case (eos_unesco)
283  call calculate_density_unesco(t, s, pressure, rho, start, npts, rho_ref)
284  case (eos_wright)
285  call calculate_density_wright(t, s, pressure, rho, start, npts, rho_ref)
286  case (eos_teos10)
287  call calculate_density_teos10(t, s, pressure, rho, start, npts, rho_ref)
288  case (eos_nemo)
289  call calculate_density_nemo(t, s, pressure, rho, start, npts, rho_ref)
290  case default
291  call mom_error(fatal, "calculate_density_array: EOS%form_of_EOS is not valid.")
292  end select
293 
294  if (present(scale)) then ; if (scale /= 1.0) then ; do j=start,start+npts-1
295  rho(j) = scale * rho(j)
296  enddo ; endif ; endif
297 

◆ calculate_density_scalar()

subroutine mom_eos::calculate_density::calculate_density_scalar ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  rho,
type(eos_type), pointer  EOS,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate density of sea water for scalar inputs. If rho_ref is present, the anomaly with respect to rho_ref is returned. The pressure and density can be rescaled with the US. If both the US and scale arguments are present the density scaling uses the product of the two scaling factors.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [Pa] or [R L2 T-2 ~> Pa]
[out]rhoDensity (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
eosEquation of state structure
[in]rho_refA reference density [kg m-3]
[in]scaleA multiplicative factor by which to scale density in combination with scaling given by US [various]

Definition at line 165 of file MOM_EOS.F90.

166  real, intent(in) :: T !< Potential temperature referenced to the surface [degC]
167  real, intent(in) :: S !< Salinity [ppt]
168  real, intent(in) :: pressure !< Pressure [Pa] or [R L2 T-2 ~> Pa]
169  real, intent(out) :: rho !< Density (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
170  type(EOS_type), pointer :: EOS !< Equation of state structure
171  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3]
172  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density in
173  !! combination with scaling given by US [various]
174 
175  real :: rho_scale ! A factor to convert density from kg m-3 to the desired units [R m3 kg-1 ~> 1]
176  real :: p_scale ! A factor to convert pressure to units of Pa [Pa T2 R-1 L-2 ~> 1]
177 
178  if (.not.associated(eos)) call mom_error(fatal, &
179  "calculate_density_scalar called with an unassociated EOS_type EOS.")
180 
181  p_scale = eos%RL2_T2_to_Pa
182 
183  select case (eos%form_of_EOS)
184  case (eos_linear)
185  call calculate_density_linear(t, s, p_scale*pressure, rho, &
186  eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
187  case (eos_unesco)
188  call calculate_density_unesco(t, s, p_scale*pressure, rho, rho_ref)
189  case (eos_wright)
190  call calculate_density_wright(t, s, p_scale*pressure, rho, rho_ref)
191  case (eos_teos10)
192  call calculate_density_teos10(t, s, p_scale*pressure, rho, rho_ref)
193  case (eos_nemo)
194  call calculate_density_nemo(t, s, p_scale*pressure, rho, rho_ref)
195  case default
196  call mom_error(fatal, "calculate_density_scalar: EOS is not valid.")
197  end select
198 
199  rho_scale = eos%kg_m3_to_R
200  if (present(scale)) rho_scale = rho_scale * scale
201  rho = rho_scale * rho
202 

◆ calculate_stanley_density_1d()

subroutine mom_eos::calculate_density::calculate_stanley_density_1d ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(in)  Tvar,
real, dimension(:), intent(in)  TScov,
real, dimension(:), intent(in)  Svar,
real, dimension(:), intent(inout)  rho,
type(eos_type), pointer  EOS,
integer, dimension(2), intent(in), optional  dom,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including the variance of T, S and covariance of T-S, potentially limiting the domain of indices that are worked on. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [R L2 T-2 ~> Pa]
[in]tvarVariance of potential temperature [degC2]
[in]tscovCovariance of potential temperature and salinity [degC ppt]
[in]svarVariance of salinity [ppt2]
[in,out]rhoDensity (in-situ if pressure is local) [R ~> kg m-3]
eosEquation of state structure
[in]domThe domain of indices to work on, taking into account that arrays start at 1.
[in]rho_refA reference density [kg m-3]
[in]scaleA multiplicative factor by which to scale density in combination with scaling given by US [various]

Definition at line 416 of file MOM_EOS.F90.

417  real, dimension(:), intent(in) :: T !< Potential temperature referenced to the surface [degC]
418  real, dimension(:), intent(in) :: S !< Salinity [ppt]
419  real, dimension(:), intent(in) :: pressure !< Pressure [R L2 T-2 ~> Pa]
420  real, dimension(:), intent(in) :: Tvar !< Variance of potential temperature [degC2]
421  real, dimension(:), intent(in) :: TScov !< Covariance of potential temperature and salinity [degC ppt]
422  real, dimension(:), intent(in) :: Svar !< Variance of salinity [ppt2]
423  real, dimension(:), intent(inout) :: rho !< Density (in-situ if pressure is local) [R ~> kg m-3]
424  type(EOS_type), pointer :: EOS !< Equation of state structure
425  integer, dimension(2), optional, intent(in) :: dom !< The domain of indices to work on, taking
426  !! into account that arrays start at 1.
427  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3]
428  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density
429  !! in combination with scaling given by US [various]
430  ! Local variables
431  real :: p_scale ! A factor to convert pressure to units of Pa [Pa T2 R-1 L-2 ~> 1]
432  real :: rho_scale ! A factor to convert density from kg m-3 to the desired units [R m3 kg-1 ~> 1]
433  real, dimension(size(rho)) :: pres ! Pressure converted to [Pa]
434  real, dimension(size(T)) :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp ! Second derivatives of density wrt T,S,p
435  integer :: i, is, ie, npts
436 
437  if (.not.associated(eos)) call mom_error(fatal, &
438  "calculate_density_1d called with an unassociated EOS_type EOS.")
439 
440  if (present(dom)) then
441  is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
442  else
443  is = 1 ; ie = size(rho) ; npts = 1 + ie - is
444  endif
445 
446  p_scale = eos%RL2_T2_to_Pa
447  do i=is,ie
448  pres(i) = p_scale * pressure(i)
449  enddo
450 
451  select case (eos%form_of_EOS)
452  case (eos_linear)
453  call calculate_density_linear(t, s, pres, rho, 1, npts, &
454  eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
455  call calculate_density_second_derivs_linear(t, s, pres, d2rdss, d2rdst, &
456  d2rdtt, d2rdsp, d2rdtp, 1, npts)
457  case (eos_wright)
458  call calculate_density_wright(t, s, pres, rho, 1, npts, rho_ref)
459  call calculate_density_second_derivs_wright(t, s, pres, d2rdss, d2rdst, &
460  d2rdtt, d2rdsp, d2rdtp, 1, npts)
461  case (eos_teos10)
462  call calculate_density_teos10(t, s, pres, rho, 1, npts, rho_ref)
463  call calculate_density_second_derivs_teos10(t, s, pres, d2rdss, d2rdst, &
464  d2rdtt, d2rdsp, d2rdtp, 1, npts)
465  case default
466  call mom_error(fatal, "calculate_stanley_density_scalar: EOS is not valid.")
467  end select
468 
469  ! Equation 25 of Stanley et al., 2020.
470  do i=is,ie
471  rho(i) = rho(i) &
472  + ( 0.5 * d2rdtt(i) * tvar(i) + ( d2rdst(i) * tscov(i) + 0.5 * d2rdss(i) * svar(i) ) )
473  enddo
474 
475  rho_scale = eos%kg_m3_to_R
476  if (present(scale)) rho_scale = rho_scale * scale
477  if (rho_scale /= 1.0) then ; do i=is,ie
478  rho(i) = rho_scale * rho(i)
479  enddo ; endif
480 

◆ calculate_stanley_density_array()

subroutine mom_eos::calculate_density::calculate_stanley_density_array ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(in)  Tvar,
real, dimension(:), intent(in)  TScov,
real, dimension(:), intent(in)  Svar,
real, dimension(:), intent(inout)  rho,
integer, intent(in)  start,
integer, intent(in)  npts,
type(eos_type), pointer  EOS,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including the variance of T, S and covariance of T-S. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [Pa]
[in]tvarVariance of potential temperature referenced to the surface [degC2]
[in]tscovCovariance of potential temperature and salinity [degC ppt]
[in]svarVariance of salinity [ppt2]
[in,out]rhoDensity (in-situ if pressure is local) [kg m-3]
[in]startStart index for computation
[in]nptsNumber of point to compute
eosEquation of state structure
[in]rho_refA reference density [kg m-3].
[in]scaleA multiplicative factor by which to scale density from kg m-3 to the desired units [R m3 kg-1]

Definition at line 305 of file MOM_EOS.F90.

306  real, dimension(:), intent(in) :: T !< Potential temperature referenced to the surface [degC]
307  real, dimension(:), intent(in) :: S !< Salinity [ppt]
308  real, dimension(:), intent(in) :: pressure !< Pressure [Pa]
309  real, dimension(:), intent(in) :: Tvar !< Variance of potential temperature referenced to the surface [degC2]
310  real, dimension(:), intent(in) :: TScov !< Covariance of potential temperature and salinity [degC ppt]
311  real, dimension(:), intent(in) :: Svar !< Variance of salinity [ppt2]
312  real, dimension(:), intent(inout) :: rho !< Density (in-situ if pressure is local) [kg m-3]
313  integer, intent(in) :: start !< Start index for computation
314  integer, intent(in) :: npts !< Number of point to compute
315  type(EOS_type), pointer :: EOS !< Equation of state structure
316  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
317  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density
318  !! from kg m-3 to the desired units [R m3 kg-1]
319  ! Local variables
320  real, dimension(size(T)) :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp ! Second derivatives of density wrt T,S,p
321  integer :: j
322 
323  if (.not.associated(eos)) call mom_error(fatal, &
324  "calculate_density_array called with an unassociated EOS_type EOS.")
325 
326  select case (eos%form_of_EOS)
327  case (eos_linear)
328  call calculate_density_linear(t, s, pressure, rho, start, npts, &
329  eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
330  call calculate_density_second_derivs_linear(t, s, pressure, d2rdss, d2rdst, &
331  d2rdtt, d2rdsp, d2rdtp, start, npts)
332  case (eos_wright)
333  call calculate_density_wright(t, s, pressure, rho, start, npts, rho_ref)
334  call calculate_density_second_derivs_wright(t, s, pressure, d2rdss, d2rdst, &
335  d2rdtt, d2rdsp, d2rdtp, start, npts)
336  case (eos_teos10)
337  call calculate_density_teos10(t, s, pressure, rho, start, npts, rho_ref)
338  call calculate_density_second_derivs_teos10(t, s, pressure, d2rdss, d2rdst, &
339  d2rdtt, d2rdsp, d2rdtp, start, npts)
340  case default
341  call mom_error(fatal, "calculate_stanley_density_array: EOS%form_of_EOS is not valid.")
342  end select
343 
344  ! Equation 25 of Stanley et al., 2020.
345  do j=start,start+npts-1
346  rho(j) = rho(j) &
347  + ( 0.5 * d2rdtt(j) * tvar(j) + ( d2rdst(j) * tscov(j) + 0.5 * d2rdss(j) * svar(j) ) )
348  enddo
349 
350  if (present(scale)) then ; if (scale /= 1.0) then ; do j=start,start+npts-1
351  rho(j) = scale * rho(j)
352  enddo ; endif ; endif
353 

◆ calculate_stanley_density_scalar()

subroutine mom_eos::calculate_density::calculate_stanley_density_scalar ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(in)  Tvar,
real, intent(in)  TScov,
real, intent(in)  Svar,
real, intent(out)  rho,
type(eos_type), pointer  EOS,
real, intent(in), optional  rho_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate density of sea water for scalar inputs including the variance of T, S and covariance of T-S. The calculation uses only the second order correction in a series as discussed in Stanley et al., 2020. If rho_ref is present, the anomaly with respect to rho_ref is returned. The density can be rescaled using rho_ref.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]tvarVariance of potential temperature referenced to the surface [degC2]
[in]tscovCovariance of potential temperature and salinity [degC ppt]
[in]svarVariance of salinity [ppt2]
[in]pressurePressure [Pa]
[out]rhoDensity (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
eosEquation of state structure
[in]rho_refA reference density [kg m-3].
[in]scaleA multiplicative factor by which to scale density from kg m-3 to the desired units [R m3 kg-1]

Definition at line 211 of file MOM_EOS.F90.

212  real, intent(in) :: T !< Potential temperature referenced to the surface [degC]
213  real, intent(in) :: S !< Salinity [ppt]
214  real, intent(in) :: Tvar !< Variance of potential temperature referenced to the surface [degC2]
215  real, intent(in) :: TScov !< Covariance of potential temperature and salinity [degC ppt]
216  real, intent(in) :: Svar !< Variance of salinity [ppt2]
217  real, intent(in) :: pressure !< Pressure [Pa]
218  real, intent(out) :: rho !< Density (in-situ if pressure is local) [kg m-3] or [R ~> kg m-3]
219  type(EOS_type), pointer :: EOS !< Equation of state structure
220  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
221  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale density
222  !! from kg m-3 to the desired units [R m3 kg-1]
223  ! Local variables
224  real :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp ! Second derivatives of density wrt T,S,p
225  real :: rho_scale ! A factor to convert density from kg m-3 to the desired units [R m3 kg-1 ~> 1]
226  real :: p_scale ! A factor to convert pressure to units of Pa [Pa T2 R-1 L-2 ~> 1]
227 
228  if (.not.associated(eos)) call mom_error(fatal, &
229  "calculate_stanley_density_scalar called with an unassociated EOS_type EOS.")
230 
231  p_scale = eos%RL2_T2_to_Pa
232 
233  select case (eos%form_of_EOS)
234  case (eos_linear)
235  call calculate_density_linear(t, s, p_scale*pressure, rho, &
236  eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
237  call calculate_density_second_derivs_linear(t, s, pressure, d2rdss, d2rdst, &
238  d2rdtt, d2rdsp, d2rdtp)
239  case (eos_wright)
240  call calculate_density_wright(t, s, p_scale*pressure, rho, rho_ref)
241  call calculate_density_second_derivs_wright(t, s, pressure, d2rdss, d2rdst, &
242  d2rdtt, d2rdsp, d2rdtp)
243  case (eos_teos10)
244  call calculate_density_teos10(t, s, p_scale*pressure, rho, rho_ref)
245  call calculate_density_second_derivs_teos10(t, s, pressure, d2rdss, d2rdst, &
246  d2rdtt, d2rdsp, d2rdtp)
247  case default
248  call mom_error(fatal, "calculate_stanley_density_scalar: EOS is not valid.")
249  end select
250 
251  ! Equation 25 of Stanley et al., 2020.
252  rho = rho + ( 0.5 * d2rdtt * tvar + ( d2rdst * tscov + 0.5 * d2rdss * svar ) )
253 
254  rho_scale = eos%kg_m3_to_R
255  if (present(scale)) rho_scale = rho_scale * scale
256  rho = rho_scale * rho
257 

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