MOM6
mom_eos_nemo Module Reference

Detailed Description

The equation of state using the expressions of Roquet et al. that are used in NEMO.

Data Types

interface  calculate_density_derivs_nemo
 For a given thermodynamic state, return the derivatives of density with conservative temperature and absolute salinity, the expressions derived for use with NEMO. More...
 
interface  calculate_density_nemo
 Compute the in situ density of sea water ([kg m-3]), or its anomaly with respect to a reference density, from absolute salinity (g/kg), conservative temperature (in deg C), and pressure [Pa], using the expressions derived for use with NEMO. More...
 

Functions/Subroutines

subroutine, public calculate_density_scalar_nemo (T, S, pressure, rho, rho_ref)
 This subroutine computes the in situ density of sea water (rho in [kg m-3]) from absolute salinity (S [g kg-1]), conservative temperature (T [degC]), and pressure [Pa]. It uses the expressions derived for use with NEMO. More...
 
subroutine, public calculate_density_array_nemo (T, S, pressure, rho, start, npts, rho_ref)
 This subroutine computes the in situ density of sea water (rho in [kg m-3]) from absolute salinity (S [g kg-1]), conservative temperature (T [degC]), and pressure [Pa]. It uses the expressions derived for use with NEMO. More...
 
subroutine calculate_density_derivs_array_nemo (T, S, pressure, drho_dT, drho_dS, start, npts)
 For a given thermodynamic state, calculate the derivatives of density with conservative temperature and absolute salinity, using the expressions derived for use with NEMO. More...
 
subroutine calculate_density_derivs_scalar_nemo (T, S, pressure, drho_dt, drho_ds)
 Wrapper to calculate_density_derivs_array for scalar inputs. More...
 
subroutine, public calculate_compress_nemo (T, S, pressure, rho, drho_dp, start, npts)
 Compute the in situ density of sea water (rho in [kg m-3]) and the compressibility (drho/dp = C_sound^-2, stored as drho_dp [s2 m-2]) from absolute salinity (sal in g/kg), conservative temperature (T [degC]), and pressure [Pa], using the expressions derived for use with NEMO. More...
 

Variables

real, parameter pa2db = 1.e-4
 Conversion factor between Pa and dbar.
 

Function/Subroutine Documentation

◆ calculate_compress_nemo()

subroutine, public mom_eos_nemo::calculate_compress_nemo ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(out)  rho,
real, dimension(:), intent(out)  drho_dp,
integer, intent(in)  start,
integer, intent(in)  npts 
)

Compute the in situ density of sea water (rho in [kg m-3]) and the compressibility (drho/dp = C_sound^-2, stored as drho_dp [s2 m-2]) from absolute salinity (sal in g/kg), conservative temperature (T [degC]), and pressure [Pa], using the expressions derived for use with NEMO.

Parameters
[in]tConservative temperature [degC].
[in]sAbsolute salinity [g/kg].
[in]pressurepressure [Pa].
[out]rhoIn situ density [kg m-3].
[out]drho_dpThe partial derivative of density with pressure (also the inverse of the square of sound speed) [s2 m-2].
[in]startThe starting point in the arrays.
[in]nptsThe number of values to calculate.

Definition at line 369 of file MOM_EOS_NEMO.F90.

369  real, intent(in), dimension(:) :: T !< Conservative temperature [degC].
370  real, intent(in), dimension(:) :: S !< Absolute salinity [g/kg].
371  real, intent(in), dimension(:) :: pressure !< pressure [Pa].
372  real, intent(out), dimension(:) :: rho !< In situ density [kg m-3].
373  real, intent(out), dimension(:) :: drho_dp !< The partial derivative of density with pressure
374  !! (also the inverse of the square of sound speed)
375  !! [s2 m-2].
376  integer, intent(in) :: start !< The starting point in the arrays.
377  integer, intent(in) :: npts !< The number of values to calculate.
378 
379  ! Local variables
380  real :: zs,zt,zp
381  integer :: j
382 
383  call calculate_density_array_nemo(t, s, pressure, rho, start, npts)
384  !
385  !NOTE: The following calculates the TEOS10 approximation to compressibility
386  ! since the corresponding NEMO approximation is not available yet.
387  !
388  do j=start,start+npts-1
389  !Conversions
390  zs = s(j) !gsw_sr_from_sp(S(j)) !Convert practical salinity to absolute salinity
391  zt = t(j) !gsw_ct_from_pt(S(j),T(j)) !Convert potantial temp to conservative temp
392  zp = pressure(j)* pa2db !Convert pressure from Pascal to decibar
393  call gsw_rho_first_derivatives(zs,zt,zp, drho_dp=drho_dp(j))
394  enddo

◆ calculate_density_array_nemo()

subroutine, public mom_eos_nemo::calculate_density_array_nemo ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(out)  rho,
integer, intent(in)  start,
integer, intent(in)  npts,
real, intent(in), optional  rho_ref 
)

This subroutine computes the in situ density of sea water (rho in [kg m-3]) from absolute salinity (S [g kg-1]), conservative temperature (T [degC]), and pressure [Pa]. It uses the expressions derived for use with NEMO.

Parameters
[in]tConservative temperature [degC].
[in]sAbsolute salinity [g kg-1].
[in]pressurepressure [Pa].
[out]rhoin situ density [kg m-3].
[in]startthe starting point in the arrays.
[in]nptsthe number of values to calculate.
[in]rho_refA reference density [kg m-3].

Definition at line 206 of file MOM_EOS_NEMO.F90.

206  real, dimension(:), intent(in) :: T !< Conservative temperature [degC].
207  real, dimension(:), intent(in) :: S !< Absolute salinity [g kg-1].
208  real, dimension(:), intent(in) :: pressure !< pressure [Pa].
209  real, dimension(:), intent(out) :: rho !< in situ density [kg m-3].
210  integer, intent(in) :: start !< the starting point in the arrays.
211  integer, intent(in) :: npts !< the number of values to calculate.
212  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
213 
214  ! Local variables
215  real :: zp, zt, zh, zs, zr0, zn, zn0, zn1, zn2, zn3, zs0
216  integer :: j
217 
218  do j=start,start+npts-1
219  !Conversions
220  zs = s(j) !gsw_sr_from_sp(S(j)) !Convert practical salinity to absolute salinity
221  zt = t(j) !gsw_ct_from_pt(S(j),T(j)) !Convert potential temp to conservative temp
222  zp = pressure(j)* pa2db !Convert pressure from Pascal to decibar
223 
224  !The following algorithm was provided by Roquet in a private communication.
225  !It is not necessarily the algorithm used in NEMO ocean!
226  zp = zp * r1_p0 !pressure
227  zt = zt * r1_t0 !temperature
228  zs = sqrt( abs( zs + rdeltas ) * r1_s0 ) ! square root salinity
229 
230  zn3 = eos013*zt &
231  & + eos103*zs+eos003
232 
233  zn2 = (eos022*zt &
234  & + eos112*zs+eos012)*zt &
235  & + (eos202*zs+eos102)*zs+eos002
236 
237  zn1 = (((eos041*zt &
238  & + eos131*zs+eos031)*zt &
239  & + (eos221*zs+eos121)*zs+eos021)*zt &
240  & + ((eos311*zs+eos211)*zs+eos111)*zs+eos011)*zt &
241  & + (((eos401*zs+eos301)*zs+eos201)*zs+eos101)*zs+eos001
242 
243  zn0 = (((((eos060*zt &
244  & + eos150*zs+eos050)*zt &
245  & + (eos240*zs+eos140)*zs+eos040)*zt &
246  & + ((eos330*zs+eos230)*zs+eos130)*zs+eos030)*zt &
247  & + (((eos420*zs+eos320)*zs+eos220)*zs+eos120)*zs+eos020)*zt &
248  & + ((((eos510*zs+eos410)*zs+eos310)*zs+eos210)*zs+eos110)*zs+eos010)*zt
249 
250  zs0 = (((((eos600*zs+eos500)*zs+eos400)*zs+eos300)*zs+eos200)*zs+eos100)*zs + eos000
251 
252  zr0 = (((((r05 * zp+r04) * zp+r03 ) * zp+r02 ) * zp+r01) * zp+r00) * zp
253 
254  if (present(rho_ref)) then
255  zn = ( ( zn3 * zp + zn2 ) * zp + zn1 ) * zp + (zn0 + (zs0 - rho_ref))
256  rho(j) = ( zn + zr0 ) ! density
257  else
258  zn = ( ( zn3 * zp + zn2 ) * zp + zn1 ) * zp + (zn0 + zs0)
259  rho(j) = ( zn + zr0 ) ! density
260  endif
261 
262  enddo

◆ calculate_density_derivs_array_nemo()

subroutine mom_eos_nemo::calculate_density_derivs_array_nemo ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(out)  drho_dT,
real, dimension(:), intent(out)  drho_dS,
integer, intent(in)  start,
integer, intent(in)  npts 
)
private

For a given thermodynamic state, calculate the derivatives of density with conservative temperature and absolute salinity, using the expressions derived for use with NEMO.

Parameters
[in]tConservative temperature [degC].
[in]sAbsolute salinity [g kg-1].
[in]pressurepressure [Pa].
[out]drho_dtThe partial derivative of density with potential temperature [kg m-3 degC-1].
[out]drho_dsThe partial derivative of density with salinity, in [kg m-3 ppt-1].
[in]startThe starting point in the arrays.
[in]nptsThe number of values to calculate.

Definition at line 268 of file MOM_EOS_NEMO.F90.

268  real, intent(in), dimension(:) :: T !< Conservative temperature [degC].
269  real, intent(in), dimension(:) :: S !< Absolute salinity [g kg-1].
270  real, intent(in), dimension(:) :: pressure !< pressure [Pa].
271  real, intent(out), dimension(:) :: drho_dT !< The partial derivative of density with potential
272  !! temperature [kg m-3 degC-1].
273  real, intent(out), dimension(:) :: drho_dS !< The partial derivative of density with salinity,
274  !! in [kg m-3 ppt-1].
275  integer, intent(in) :: start !< The starting point in the arrays.
276  integer, intent(in) :: npts !< The number of values to calculate.
277 
278  ! Local variables
279  real :: zp,zt , zh , zs , zr0, zn , zn0, zn1, zn2, zn3
280  integer :: j
281 
282  do j=start,start+npts-1
283  !Conversions
284  zs = s(j) !gsw_sr_from_sp(S(j)) !Convert practical salinity to absolute salinity
285  zt = t(j) !gsw_ct_from_pt(S(j),T(j)) !Convert potantial temp to conservative temp
286  zp = pressure(j)* pa2db !Convert pressure from Pascal to decibar
287 
288  !The following algorithm was provided by Roquet in a private communication.
289  !It is not necessarily the algorithm used in NEMO ocean!
290  zp = zp * r1_p0 ! pressure (first converted to decibar)
291  zt = zt * r1_t0 ! temperature
292  zs = sqrt( abs( zs + rdeltas ) * r1_s0 ) ! square root salinity
293  !
294  ! alpha
295  zn3 = alp003
296  !
297  zn2 = alp012*zt + alp102*zs+alp002
298  !
299  zn1 = ((alp031*zt &
300  & + alp121*zs+alp021)*zt &
301  & + (alp211*zs+alp111)*zs+alp011)*zt &
302  & + ((alp301*zs+alp201)*zs+alp101)*zs+alp001
303  !
304  zn0 = ((((alp050*zt &
305  & + alp140*zs+alp040)*zt &
306  & + (alp230*zs+alp130)*zs+alp030)*zt &
307  & + ((alp320*zs+alp220)*zs+alp120)*zs+alp020)*zt &
308  & + (((alp410*zs+alp310)*zs+alp210)*zs+alp110)*zs+alp010)*zt &
309  & + ((((alp500*zs+alp400)*zs+alp300)*zs+alp200)*zs+alp100)*zs+alp000
310  !
311  zn = ( ( zn3 * zp + zn2 ) * zp + zn1 ) * zp + zn0
312  !
313  drho_dt(j) = -zn
314  !
315  ! beta
316  !
317  zn3 = bet003
318  !
319  zn2 = bet012*zt + bet102*zs+bet002
320  !
321  zn1 = ((bet031*zt &
322  & + bet121*zs+bet021)*zt &
323  & + (bet211*zs+bet111)*zs+bet011)*zt &
324  & + ((bet301*zs+bet201)*zs+bet101)*zs+bet001
325  !
326  zn0 = ((((bet050*zt &
327  & + bet140*zs+bet040)*zt &
328  & + (bet230*zs+bet130)*zs+bet030)*zt &
329  & + ((bet320*zs+bet220)*zs+bet120)*zs+bet020)*zt &
330  & + (((bet410*zs+bet310)*zs+bet210)*zs+bet110)*zs+bet010)*zt &
331  & + ((((bet500*zs+bet400)*zs+bet300)*zs+bet200)*zs+bet100)*zs+bet000
332  !
333  zn = ( ( zn3 * zp + zn2 ) * zp + zn1 ) * zp + zn0
334  !
335  drho_ds(j) = zn / zs
336  enddo
337 

◆ calculate_density_derivs_scalar_nemo()

subroutine mom_eos_nemo::calculate_density_derivs_scalar_nemo ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  drho_dt,
real, intent(out)  drho_ds 
)
private

Wrapper to calculate_density_derivs_array for scalar inputs.

Parameters
[in]tPotential temperature relative to the surface [degC].
[in]sSalinity [g kg-1].
[in]pressurePressure [Pa].
[out]drho_dtThe partial derivative of density with potential temperature [kg m-3 degC-1].
[out]drho_dsThe partial derivative of density with salinity, in [kg m-3 ppt-1].

Definition at line 342 of file MOM_EOS_NEMO.F90.

342  real, intent(in) :: T !< Potential temperature relative to the surface [degC].
343  real, intent(in) :: S !< Salinity [g kg-1].
344  real, intent(in) :: pressure !< Pressure [Pa].
345  real, intent(out) :: drho_dT !< The partial derivative of density with potential
346  !! temperature [kg m-3 degC-1].
347  real, intent(out) :: drho_dS !< The partial derivative of density with salinity,
348  !! in [kg m-3 ppt-1].
349  ! Local variables
350  real :: al0, p0, lambda
351  integer :: j
352  real, dimension(1) :: T0, S0, pressure0
353  real, dimension(1) :: drdt0, drds0
354 
355  t0(1) = t
356  s0(1) = s
357  pressure0(1) = pressure
358 
359  call calculate_density_derivs_array_nemo(t0, s0, pressure0, drdt0, drds0, 1, 1)
360  drho_dt = drdt0(1)
361  drho_ds = drds0(1)

◆ calculate_density_scalar_nemo()

subroutine, public mom_eos_nemo::calculate_density_scalar_nemo ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  rho,
real, intent(in), optional  rho_ref 
)

This subroutine computes the in situ density of sea water (rho in [kg m-3]) from absolute salinity (S [g kg-1]), conservative temperature (T [degC]), and pressure [Pa]. It uses the expressions derived for use with NEMO.

Parameters
[in]tConservative temperature [degC].
[in]sAbsolute salinity [g kg-1].
[in]pressurepressure [Pa].
[out]rhoIn situ density [kg m-3].
[in]rho_refA reference density [kg m-3].

Definition at line 181 of file MOM_EOS_NEMO.F90.

181  real, intent(in) :: T !< Conservative temperature [degC].
182  real, intent(in) :: S !< Absolute salinity [g kg-1].
183  real, intent(in) :: pressure !< pressure [Pa].
184  real, intent(out) :: rho !< In situ density [kg m-3].
185  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
186 
187  real :: al0, p0, lambda
188  integer :: j
189  real, dimension(1) :: T0, S0, pressure0
190  real, dimension(1) :: rho0
191 
192  t0(1) = t
193  s0(1) = s
194  pressure0(1) = pressure
195 
196  call calculate_density_array_nemo(t0, s0, pressure0, rho0, 1, 1, rho_ref)
197  rho = rho0(1)
198