MOM6
mom_eos_unesco Module Reference

Data Types

interface  calculate_density_unesco
 

Functions/Subroutines

subroutine, public calculate_density_scalar_unesco (T, S, pressure, rho)
 This subroutine computes the in situ density of sea water (rho in units of kg/m^3) from salinity (S in psu), potential temperature (T in deg C), and pressure in Pa. It uses the expression from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740. Coded by R. Hallberg, 7/00. More...
 
subroutine, public calculate_density_array_unesco (T, S, pressure, rho, start, npts)
 This subroutine computes the in situ density of sea water (rho in units of kg/m^3) from salinity (S in psu), potential temperature (T in deg C), and pressure in Pa. More...
 
subroutine, public calculate_density_derivs_unesco (T, S, pressure, drho_dT, drho_dS, start, npts)
 This subroutine calculates the partial derivatives of density with potential temperature and salinity. More...
 
subroutine, public calculate_compress_unesco (T, S, pressure, rho, drho_dp, start, npts)
 This subroutine computes the in situ density of sea water (rho) and the compressibility (drho/dp == C_sound^-2) at the given salinity, potential temperature, and pressure. More...
 

Variables

real, parameter r00 = 999.842594
 
real, parameter r10 = 6.793952e-2
 
real, parameter r20 = -9.095290e-3
 
real, parameter r30 = 1.001685e-4
 
real, parameter r40 = -1.120083e-6
 
real, parameter r50 = 6.536332e-9
 
real, parameter r01 = 0.824493
 
real, parameter r11 = -4.0899e-3
 
real, parameter r21 = 7.6438e-5
 
real, parameter r31 = -8.2467e-7
 
real, parameter r41 = 5.3875e-9
 
real, parameter r032 = -5.72466e-3
 
real, parameter r132 = 1.0227e-4
 
real, parameter r232 = -1.6546e-6
 
real, parameter r02 = 4.8314e-4
 
real, parameter s00 = 1.965933e4
 
real, parameter s10 = 1.444304e2
 
real, parameter s20 = -1.706103
 
real, parameter s30 = 9.648704e-3
 
real, parameter s40 = -4.190253e-5
 
real, parameter s01 = 52.84855
 
real, parameter s11 = -3.101089e-1
 
real, parameter s21 = 6.283263e-3
 
real, parameter s31 = -5.084188e-5
 
real, parameter s032 = 3.886640e-1
 
real, parameter s132 = 9.085835e-3
 
real, parameter s232 = -4.619924e-4
 
real, parameter sp00 = 3.186519
 
real, parameter sp10 = 2.212276e-2
 
real, parameter sp20 = -2.984642e-4
 
real, parameter sp30 = 1.956415e-6
 
real, parameter sp01 = 6.704388e-3
 
real, parameter sp11 = -1.847318e-4
 
real, parameter sp21 = 2.059331e-7
 
real, parameter sp032 = 1.480266e-4
 
real, parameter sp000 = 2.102898e-4
 
real, parameter sp010 = -1.202016e-5
 
real, parameter sp020 = 1.394680e-7
 
real, parameter sp001 = -2.040237e-6
 
real, parameter sp011 = 6.128773e-8
 
real, parameter sp021 = 6.207323e-10
 

Function/Subroutine Documentation

◆ calculate_compress_unesco()

subroutine, public mom_eos_unesco::calculate_compress_unesco ( 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 
)

This subroutine computes the in situ density of sea water (rho) and the compressibility (drho/dp == C_sound^-2) at the given salinity, potential temperature, and pressure.

Parameters
[in]tPotential temperature relative to the surface in C.
[in]sSalinity in PSU.
[in]pressurePressure in Pa.
[out]rhoIn situ density in kg m-3.
[out]drho_dpThe partial derivative of density with pressure (also the inverse of the square of sound speed) in s2 m-2.
[in]startThe starting point in the arrays.
[in]nptsThe number of values to calculate.

Definition at line 232 of file MOM_EOS_UNESCO.F90.

References r00, r01, r02, r032, r10, r11, r132, r20, r21, r232, r30, r31, r40, r41, r50, s00, s01, s032, s10, s11, s132, s20, s21, s232, s30, s31, s40, sp00, sp000, sp001, sp01, sp010, sp011, sp020, sp021, sp032, sp10, sp11, sp20, sp21, and sp30.

Referenced by mom_eos::calculate_compress().

232  real, intent(in), dimension(:) :: t !< Potential temperature relative to the surface
233  !! in C.
234  real, intent(in), dimension(:) :: s !< Salinity in PSU.
235  real, intent(in), dimension(:) :: pressure !< Pressure in Pa.
236  real, intent(out), dimension(:) :: rho !< In situ density in kg m-3.
237  real, intent(out), dimension(:) :: drho_dp !< The partial derivative of density with pressure
238  !! (also the inverse of the square of sound speed)
239  !! in s2 m-2.
240  integer, intent(in) :: start !< The starting point in the arrays.
241  integer, intent(in) :: npts !< The number of values to calculate.
242 
243 ! * This subroutine computes the in situ density of sea water (rho) *
244 ! * and the compressibility (drho/dp == C_sound^-2) at the given *
245 ! * salinity, potential temperature, and pressure. *
246 ! * *
247 ! * Arguments: T - potential temperature relative to the surface in C. *
248 ! * (in) S - salinity in PSU. *
249 ! * (in) pressure - pressure in Pa. *
250 ! * (out) rho - in situ density in kg m-3. *
251 ! * (out) drho_dp - the partial derivative of density with *
252 ! * pressure (also the inverse of the square of *
253 ! * sound speed) in s2 m-2. *
254 ! * (in) start - the starting point in the arrays. *
255 ! * (in) npts - the number of values to calculate. *
256 
257  real :: t_local, t2, t3, t4, t5; ! Temperature to the 1st - 5th power.
258  real :: s_local, s32, s2; ! Salinity to the 1st, 3/2, & 2nd power.
259  real :: p1, p2; ! Pressure (in bars) to the 1st and 2nd power.
260  real :: rho0; ! Density at 1 bar pressure, in kg m-3.
261  real :: ks; ! The secant bulk modulus in bar.
262  real :: ks_0, ks_1, ks_2;
263  real :: dks_dp; ! The derivative of the secant bulk modulus
264  ! with pressure, nondimensional.
265  integer :: j
266 
267  do j=start,start+npts-1
268  p1 = pressure(j)*1.0e-5; p2 = p1*p1;
269  t_local = t(j); t2 = t_local*t_local; t3 = t_local*t2; t4 = t2*t2; t5 = t3*t2;
270  s_local = s(j); s2 = s_local*s_local; s32 = s_local*sqrt(s_local);
271 
272 ! Compute rho(s,theta,p=0) - (same as rho(s,t_insitu,p=0) ).
273 
274  rho0 = r00 + r10*t_local + r20*t2 + r30*t3 + r40*t4 + r50*t5 + &
275  s_local*(r01 + r11*t_local + r21*t2 + r31*t3 + r41*t4) + &
276  s32*(r032 + r132*t_local + r232*t2) + r02*s2;
277 
278 ! Compute rho(s,theta,p), first calculating the secant bulk modulus.
279  ks_0 = s00 + s10*t_local + s20*t2 + s30*t3 + s40*t4 + &
280  s_local*(s01 + s11*t_local + s21*t2 + s31*t3) + s32*(s032 + s132*t_local + s232*t2);
281  ks_1 = sp00 + sp10*t_local + sp20*t2 + sp30*t3 + &
282  s_local*(sp01 + sp11*t_local + sp21*t2) + sp032*s32;
283  ks_2 = sp000 + sp010*t_local + sp020*t2 + s_local*(sp001 + sp011*t_local + sp021*t2);
284 
285  ks = ks_0 + p1*ks_1 + p2*ks_2;
286  dks_dp = ks_1 + 2.0*p1*ks_2;
287 
288  rho(j) = rho0*ks / (ks - p1);
289 ! The factor of 1.0e-5 is because pressure here is in bars, not Pa.
290  drho_dp(j) = 1.0e-5 * (rho(j) / (ks - p1)) * (1.0 - dks_dp*p1/ks);
291  enddo
Here is the caller graph for this function:

◆ calculate_density_array_unesco()

subroutine, public mom_eos_unesco::calculate_density_array_unesco ( 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 
)

This subroutine computes the in situ density of sea water (rho in units of kg/m^3) from salinity (S in psu), potential temperature (T in deg C), and pressure in Pa.

Parameters
[in]tPotential temperature relative to the surface in C.
[in]sSalinity in PSU.
[in]pressurePressure in Pa.
[out]rhoIn situ density in kg m-3.
[in]startThe starting point in the arrays.
[in]nptsThe number of values to calculate.

Definition at line 103 of file MOM_EOS_UNESCO.F90.

References r00, r01, r02, r032, r10, r11, r132, r20, r21, r232, r30, r31, r40, r41, r50, s00, s01, s032, s10, s11, s132, s20, s21, s232, s30, s31, s40, sp00, sp000, sp001, sp01, sp010, sp011, sp020, sp021, sp032, sp10, sp11, sp20, sp21, and sp30.

Referenced by calculate_density_scalar_unesco().

103  real, intent(in), dimension(:) :: t !< Potential temperature relative to the surface
104  !! in C.
105  real, intent(in), dimension(:) :: s !< Salinity in PSU.
106  real, intent(in), dimension(:) :: pressure !< Pressure in Pa.
107  real, intent(out), dimension(:) :: rho !< In situ density in kg m-3.
108  integer, intent(in) :: start !< The starting point in the arrays.
109  integer, intent(in) :: npts !< The number of values to calculate.
110 
111 ! * This subroutine computes the in situ density of sea water (rho in *
112 ! * units of kg/m^3) from salinity (S in psu), potential temperature *
113 ! * (T in deg C), and pressure in Pa. *
114 
115 ! * Arguments: T - potential temperature relative to the surface in C. *
116 ! * (in) S - salinity in PSU. *
117 ! * (in) pressure - pressure in Pa. *
118 ! * (out) rho - in situ density in kg m-3. *
119 ! * (in) start - the starting point in the arrays. *
120 ! * (in) npts - the number of values to calculate. *
121 
122  real :: t_local, t2, t3, t4, t5; ! Temperature to the 1st - 5th power.
123  real :: s_local, s32, s2; ! Salinity to the 1st, 3/2, & 2nd power.
124  real :: p1, p2; ! Pressure (in bars) to the 1st and 2nd power.
125  real :: rho0; ! Density at 1 bar pressure, in kg m-3.
126  real :: ks; ! The secant bulk modulus in bar.
127  integer :: j
128 
129  do j=start,start+npts-1
130  p1 = pressure(j)*1.0e-5; p2 = p1*p1;
131  t_local = t(j); t2 = t_local*t_local; t3 = t_local*t2; t4 = t2*t2; t5 = t3*t2;
132  s_local = s(j); s2 = s_local*s_local; s32 = s_local*sqrt(s_local);
133 
134 ! Compute rho(s,theta,p=0) - (same as rho(s,t_insitu,p=0) ).
135 
136  rho0 = r00 + r10*t_local + r20*t2 + r30*t3 + r40*t4 + r50*t5 + &
137  s_local*(r01 + r11*t_local + r21*t2 + r31*t3 + r41*t4) + &
138  s32*(r032 + r132*t_local + r232*t2) + r02*s2;
139 
140 ! Compute rho(s,theta,p), first calculating the secant bulk modulus.
141 
142  ks = s00 + s10*t_local + s20*t2 + s30*t3 + s40*t4 + s_local*(s01 + s11*t_local + s21*t2 + s31*t3) + &
143  s32*(s032 + s132*t_local + s232*t2) + &
144  p1*(sp00 + sp10*t_local + sp20*t2 + sp30*t3 + &
145  s_local*(sp01 + sp11*t_local + sp21*t2) + sp032*s32) + &
146  p2*(sp000 + sp010*t_local + sp020*t2 + s_local*(sp001 + sp011*t_local + sp021*t2));
147 
148  rho(j) = rho0*ks / (ks - p1);
149  enddo
Here is the caller graph for this function:

◆ calculate_density_derivs_unesco()

subroutine, public mom_eos_unesco::calculate_density_derivs_unesco ( 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 
)

This subroutine calculates the partial derivatives of density with potential temperature and salinity.

Parameters
[in]tPotential temperature relative to the surface in C.
[in]sSalinity in PSU.
[in]pressurePressure in Pa.
[out]drho_dtThe partial derivative of density with potential temperature, in kg m-3 K-1.
[out]drho_dsThe partial derivative of density with salinity, in kg m-3 psu-1.
[in]startThe starting point in the arrays.
[in]nptsThe number of values to calculate.

Definition at line 155 of file MOM_EOS_UNESCO.F90.

References r00, r01, r02, r032, r10, r11, r132, r20, r21, r232, r30, r31, r40, r41, r50, s00, s01, s032, s10, s11, s132, s20, s21, s232, s30, s31, s40, sp00, sp000, sp001, sp01, sp010, sp011, sp020, sp021, sp032, sp10, sp11, sp20, sp21, and sp30.

155  real, intent(in), dimension(:) :: t !< Potential temperature relative to the surface
156  !! in C.
157  real, intent(in), dimension(:) :: s !< Salinity in PSU.
158  real, intent(in), dimension(:) :: pressure !< Pressure in Pa.
159  real, intent(out), dimension(:) :: drho_dt !< The partial derivative of density with potential
160  !! temperature, in kg m-3 K-1.
161  real, intent(out), dimension(:) :: drho_ds !< The partial derivative of density with salinity,
162  !! in kg m-3 psu-1.
163  integer, intent(in) :: start !< The starting point in the arrays.
164  integer, intent(in) :: npts !< The number of values to calculate.
165 
166 ! * This subroutine calculates the partial derivatives of density *
167 ! * with potential temperature and salinity. *
168 ! * *
169 ! * Arguments: T - potential temperature relative to the surface in C. *
170 ! * (in) S - salinity in PSU. *
171 ! * (in) pressure - pressure in Pa. *
172 ! * (out) drho_dT - the partial derivative of density with *
173 ! * potential temperature, in kg m-3 K-1. *
174 ! * (out) drho_dS - the partial derivative of density with *
175 ! * salinity, in kg m-3 psu-1. *
176 ! * (in) start - the starting point in the arrays. *
177 ! * (in) npts - the number of values to calculate. *
178  real :: t_local, t2, t3, t4, t5; ! Temperature to the 1st - 5th power.
179  real :: s12, s_local, s32, s2; ! Salinity to the 1/2 - 2nd powers.
180  real :: p1, p2; ! Pressure (in bars) to the 1st & 2nd power.
181  real :: rho0; ! Density at 1 bar pressure, in kg m-3.
182  real :: ks; ! The secant bulk modulus, in bar.
183  real :: drho0_dt; ! Derivative of rho0 with T, in kg m-3 K-1.
184  real :: drho0_ds; ! Derivative of rho0 with S, kg m-3 psu-1.
185  real :: dks_dt; ! Derivative of ks with T, in bar K-1.
186  real :: dks_ds; ! Derivative of ks with S, in bar psu-1.
187  real :: denom; ! 1.0 / (ks - p1) in bar-1.
188  integer :: j
189 
190  do j=start,start+npts-1
191  p1 = pressure(j)*1.0e-5; p2 = p1*p1;
192  t_local = t(j); t2 = t_local*t_local; t3 = t_local*t2; t4 = t2*t2; t5 = t3*t2;
193  s_local = s(j); s2 = s_local*s_local; s12 = sqrt(s_local); s32 = s_local*s12;
194 
195 ! compute rho(s,theta,p=0) - (same as rho(s,t_insitu,p=0) )
196 
197  rho0 = r00 + r10*t_local + r20*t2 + r30*t3 + r40*t4 + r50*t5 + &
198  s_local*(r01 + r11*t_local + r21*t2 + r31*t3 + r41*t4) + &
199  s32*(r032 + r132*t_local + r232*t2) + r02*s2;
200  drho0_dt = r10 + 2.0*r20*t_local + 3.0*r30*t2 + 4.0*r40*t3 + 5.0*r50*t4 + &
201  s_local*(r11 + 2.0*r21*t_local + 3.0*r31*t2 + 4.0*r41*t3) + &
202  s32*(r132 + 2.0*r232*t_local);
203  drho0_ds = (r01 + r11*t_local + r21*t2 + r31*t3 + r41*t4) + &
204  1.5*s12*(r032 + r132*t_local + r232*t2) + 2.0*r02*s_local;
205 
206 ! compute rho(s,theta,p)
207 
208  ks = s00 + s10*t_local + s20*t2 + s30*t3 + s40*t4 + s_local*(s01 + s11*t_local + s21*t2 + s31*t3) + &
209  s32*(s032 + s132*t_local + s232*t2) + &
210  p1*(sp00 + sp10*t_local + sp20*t2 + sp30*t3 + &
211  s_local*(sp01 + sp11*t_local + sp21*t2) + sp032*s32) + &
212  p2*(sp000 + sp010*t_local + sp020*t2 + s_local*(sp001 + sp011*t_local + sp021*t2));
213  dks_dt = s10 + 2.0*s20*t_local + 3.0*s30*t2 + 4.0*s40*t3 + &
214  s_local*(s11 + 2.0*s21*t_local + 3.0*s31*t2) + s32*(s132 + 2.0*s232*t_local) + &
215  p1*(sp10 + 2.0*sp20*t_local + 3.0*sp30*t2 + s_local*(sp11 + 2.0*sp21*t_local)) + &
216  p2*(sp010 + 2.0*sp020*t_local + s_local*(sp011 + 2.0*sp021*t_local));
217  dks_ds = (s01 + s11*t_local + s21*t2 + s31*t3) + 1.5*s12*(s032 + s132*t_local + s232*t2) + &
218  p1*(sp01 + sp11*t_local + sp21*t2 + 1.5*sp032*s12) + &
219  p2*(sp001 + sp011*t_local + sp021*t2);
220 
221  denom = 1.0 / (ks - p1);
222  drho_dt(j) = denom*(ks*drho0_dt - rho0*p1*denom*dks_dt);
223  drho_ds(j) = denom*(ks*drho0_ds - rho0*p1*denom*dks_ds);
224  enddo
225 

◆ calculate_density_scalar_unesco()

subroutine, public mom_eos_unesco::calculate_density_scalar_unesco ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  rho 
)

This subroutine computes the in situ density of sea water (rho in units of kg/m^3) from salinity (S in psu), potential temperature (T in deg C), and pressure in Pa. It uses the expression from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740. Coded by R. Hallberg, 7/00.

Parameters
[in]tPotential temperature relative to the surface in C.
[in]sSalinity in PSU.
[in]pressurePressure in Pa.
[out]rhoIn situ density in kg m-3.

Definition at line 67 of file MOM_EOS_UNESCO.F90.

References calculate_density_array_unesco().

67 real, intent(in) :: t !< Potential temperature relative to the surface in C.
68 real, intent(in) :: s !< Salinity in PSU.
69 real, intent(in) :: pressure !< Pressure in Pa.
70 real, intent(out) :: rho !< In situ density in kg m-3.
71 
72 ! * Arguments: T - potential temperature relative to the surface in C. *
73 ! * (in) S - salinity in PSU. *
74 ! * (in) pressure - pressure in Pa. *
75 ! * (out) rho - in situ density in kg m-3. *
76 ! * (in) start - the starting point in the arrays. *
77 ! * (in) npts - the number of values to calculate. *
78 
79 ! *====================================================================*
80 ! * This subroutine computes the in situ density of sea water (rho in *
81 ! * units of kg/m^3) from salinity (S in psu), potential temperature *
82 ! * (T in deg C), and pressure in Pa. It uses the expression from *
83 ! * Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740. *
84 ! * Coded by R. Hallberg, 7/00 *
85 ! *====================================================================*
86 
87  real, dimension(1) :: t0, s0, pressure0
88  real, dimension(1) :: rho0
89 
90  t0(1) = t
91  s0(1) = s
92  pressure0(1) = pressure
93 
94  call calculate_density_array_unesco(t0, s0, pressure0, rho0, 1, 1)
95  rho = rho0(1)
96 
Here is the call graph for this function:

Variable Documentation

◆ r00

real, parameter mom_eos_unesco::r00 = 999.842594
private

Definition at line 41 of file MOM_EOS_UNESCO.F90.

Referenced by calculate_compress_unesco(), calculate_density_array_unesco(), and calculate_density_derivs_unesco().

41 real, parameter :: r00 = 999.842594, r10 = 6.793952e-2, r20 = -9.095290e-3, &
42  r30 = 1.001685e-4, r40 = -1.120083e-6, r50 = 6.536332e-9, r01 = 0.824493, &
43  r11 = -4.0899e-3, r21 = 7.6438e-5, r31 = -8.2467e-7, r41 = 5.3875e-9, &
44  r032 = -5.72466e-3, r132 = 1.0227e-4, r232 = -1.6546e-6, r02 = 4.8314e-4;

◆ r01

real, parameter mom_eos_unesco::r01 = 0.824493
private

◆ r02

real, parameter mom_eos_unesco::r02 = 4.8314e-4
private

◆ r032

real, parameter mom_eos_unesco::r032 = -5.72466e-3
private

◆ r10

real, parameter mom_eos_unesco::r10 = 6.793952e-2
private

◆ r11

real, parameter mom_eos_unesco::r11 = -4.0899e-3
private

◆ r132

real, parameter mom_eos_unesco::r132 = 1.0227e-4
private

◆ r20

real, parameter mom_eos_unesco::r20 = -9.095290e-3
private

◆ r21

real, parameter mom_eos_unesco::r21 = 7.6438e-5
private

◆ r232

real, parameter mom_eos_unesco::r232 = -1.6546e-6
private

◆ r30

real, parameter mom_eos_unesco::r30 = 1.001685e-4
private

◆ r31

real, parameter mom_eos_unesco::r31 = -8.2467e-7
private

◆ r40

real, parameter mom_eos_unesco::r40 = -1.120083e-6
private

◆ r41

real, parameter mom_eos_unesco::r41 = 5.3875e-9
private

◆ r50

real, parameter mom_eos_unesco::r50 = 6.536332e-9
private

◆ s00

real, parameter mom_eos_unesco::s00 = 1.965933e4
private

Definition at line 50 of file MOM_EOS_UNESCO.F90.

Referenced by calculate_compress_unesco(), calculate_density_array_unesco(), and calculate_density_derivs_unesco().

50 real, parameter :: s00 = 1.965933e4, s10 = 1.444304e2, s20 = -1.706103, &
51  s30 = 9.648704e-3, s40 = -4.190253e-5, s01 = 52.84855, s11 = -3.101089e-1, &
52  s21 = 6.283263e-3, s31 = -5.084188e-5, s032 = 3.886640e-1, s132 = 9.085835e-3, &
53  s232 = -4.619924e-4, sp00 = 3.186519, sp10 = 2.212276e-2, sp20 = -2.984642e-4, &
54  sp30 = 1.956415e-6, sp01 = 6.704388e-3, sp11 = -1.847318e-4, sp21 = 2.059331e-7, &
55  sp032 = 1.480266e-4, sp000 = 2.102898e-4, sp010 = -1.202016e-5, sp020 = 1.394680e-7, &
56  sp001 = -2.040237e-6, sp011 = 6.128773e-8, sp021 = 6.207323e-10;

◆ s01

real, parameter mom_eos_unesco::s01 = 52.84855
private

◆ s032

real, parameter mom_eos_unesco::s032 = 3.886640e-1
private

◆ s10

real, parameter mom_eos_unesco::s10 = 1.444304e2
private

◆ s11

real, parameter mom_eos_unesco::s11 = -3.101089e-1
private

◆ s132

real, parameter mom_eos_unesco::s132 = 9.085835e-3
private

◆ s20

real, parameter mom_eos_unesco::s20 = -1.706103
private

◆ s21

real, parameter mom_eos_unesco::s21 = 6.283263e-3
private

◆ s232

real, parameter mom_eos_unesco::s232 = -4.619924e-4
private

◆ s30

real, parameter mom_eos_unesco::s30 = 9.648704e-3
private

◆ s31

real, parameter mom_eos_unesco::s31 = -5.084188e-5
private

◆ s40

real, parameter mom_eos_unesco::s40 = -4.190253e-5
private

◆ sp00

real, parameter mom_eos_unesco::sp00 = 3.186519
private

◆ sp000

real, parameter mom_eos_unesco::sp000 = 2.102898e-4
private

◆ sp001

real, parameter mom_eos_unesco::sp001 = -2.040237e-6
private

◆ sp01

real, parameter mom_eos_unesco::sp01 = 6.704388e-3
private

◆ sp010

real, parameter mom_eos_unesco::sp010 = -1.202016e-5
private

◆ sp011

real, parameter mom_eos_unesco::sp011 = 6.128773e-8
private

◆ sp020

real, parameter mom_eos_unesco::sp020 = 1.394680e-7
private

◆ sp021

real, parameter mom_eos_unesco::sp021 = 6.207323e-10
private

◆ sp032

real, parameter mom_eos_unesco::sp032 = 1.480266e-4
private

◆ sp10

real, parameter mom_eos_unesco::sp10 = 2.212276e-2
private

◆ sp11

real, parameter mom_eos_unesco::sp11 = -1.847318e-4
private

◆ sp20

real, parameter mom_eos_unesco::sp20 = -2.984642e-4
private

◆ sp21

real, parameter mom_eos_unesco::sp21 = 2.059331e-7
private

◆ sp30

real, parameter mom_eos_unesco::sp30 = 1.956415e-6
private