MOM6
MOM_EOS_Wright.F90
Go to the documentation of this file.
1 !> The equation of state using the Wright 1997 expressions
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 !***********************************************************************
7 !* The subroutines in this file implement the equation of state for *
8 !* sea water using the formulae given by Wright, 1997, J. Atmos. *
9 !* Ocean. Tech., 14, 735-740. Coded by R. Hallberg, 7/00. *
10 !***********************************************************************
11 
12 use mom_hor_index, only : hor_index_type
13 
14 implicit none ; private
15 
16 #include <MOM_memory.h>
17 
22 
23 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
24 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
25 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
26 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
27 
28 
29 !> Compute the in situ density of sea water (in [kg m-3]), or its anomaly with respect to
30 !! a reference density, from salinity (in psu), potential temperature (in deg C), and pressure [Pa],
31 !! using the expressions from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
34 end interface calculate_density_wright
35 
36 !> Compute the in situ specific volume of sea water (in [m3 kg-1]), or an anomaly with respect
37 !! to a reference specific volume, from salinity (in psu), potential temperature (in deg C), and
38 !! pressure [Pa], using the expressions from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
41 end interface calculate_spec_vol_wright
42 
43 !> For a given thermodynamic state, return the derivatives of density with temperature and salinity
46 end interface
47 
48 !> For a given thermodynamic state, return the second derivatives of density with various combinations
49 !! of temperature, salinity, and pressure
52 end interface
53 
54 !>@{ Parameters in the Wright equation of state
55 !real :: a0, a1, a2, b0, b1, b2, b3, b4, b5, c0, c1, c2, c3, c4, c5
56 ! One of the two following blocks of values should be commented out.
57 ! Following are the values for the full range formula.
58 !
59 !real, parameter :: a0 = 7.133718e-4, a1 = 2.724670e-7, a2 = -1.646582e-7
60 !real, parameter :: b0 = 5.613770e8, b1 = 3.600337e6, b2 = -3.727194e4
61 !real, parameter :: b3 = 1.660557e2, b4 = 6.844158e5, b5 = -8.389457e3
62 !real, parameter :: c0 = 1.609893e5, c1 = 8.427815e2, c2 = -6.931554
63 !real, parameter :: c3 = 3.869318e-2, c4 = -1.664201e2, c5 = -2.765195
64 
65 
66 ! Following are the values for the reduced range formula.
67 real, parameter :: a0 = 7.057924e-4, a1 = 3.480336e-7, a2 = -1.112733e-7 ! a0/a1 ~= 2028 ; a0/a2 ~= -6343
68 real, parameter :: b0 = 5.790749e8, b1 = 3.516535e6, b2 = -4.002714e4 ! b0/b1 ~= 165 ; b0/b4 ~= 974
69 real, parameter :: b3 = 2.084372e2, b4 = 5.944068e5, b5 = -9.643486e3
70 real, parameter :: c0 = 1.704853e5, c1 = 7.904722e2, c2 = -7.984422 ! c0/c1 ~= 216 ; c0/c4 ~= -740
71 real, parameter :: c3 = 5.140652e-2, c4 = -2.302158e2, c5 = -3.079464
72 !!@}
73 
74 contains
75 
76 !> This subroutine computes the in situ density of sea water (rho in
77 !! [kg m-3]) from salinity (S [PSU]), potential temperature
78 !! (T [degC]), and pressure [Pa]. It uses the expression from
79 !! Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
80 subroutine calculate_density_scalar_wright(T, S, pressure, rho, rho_ref)
81  real, intent(in) :: T !< Potential temperature relative to the surface [degC].
82  real, intent(in) :: S !< Salinity [PSU].
83  real, intent(in) :: pressure !< pressure [Pa].
84  real, intent(out) :: rho !< In situ density [kg m-3].
85  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
86 
87 ! *====================================================================*
88 ! * This subroutine computes the in situ density of sea water (rho in *
89 ! * [kg m-3]) from salinity (S [PSU]), potential temperature *
90 ! * (T [degC]), and pressure [Pa]. It uses the expression from *
91 ! * Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740. *
92 ! * Coded by R. Hallberg, 7/00 *
93 ! *====================================================================*
94 
95  real, dimension(1) :: T0, S0, pressure0, rho0
96 
97  t0(1) = t
98  s0(1) = s
99  pressure0(1) = pressure
100 
101  call calculate_density_array_wright(t0, s0, pressure0, rho0, 1, 1, rho_ref)
102  rho = rho0(1)
103 
104 end subroutine calculate_density_scalar_wright
105 
106 !> This subroutine computes the in situ density of sea water (rho in
107 !! [kg m-3]) from salinity (S [PSU]), potential temperature
108 !! (T [degC]), and pressure [Pa]. It uses the expression from
109 !! Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
110 subroutine calculate_density_array_wright(T, S, pressure, rho, start, npts, rho_ref)
111  real, dimension(:), intent(in) :: T !< potential temperature relative to the surface [degC].
112  real, dimension(:), intent(in) :: S !< salinity [PSU].
113  real, dimension(:), intent(in) :: pressure !< pressure [Pa].
114  real, dimension(:), intent(out) :: rho !< in situ density [kg m-3].
115  integer, intent(in) :: start !< the starting point in the arrays.
116  integer, intent(in) :: npts !< the number of values to calculate.
117  real, optional, intent(in) :: rho_ref !< A reference density [kg m-3].
118 
119  ! Original coded by R. Hallberg, 7/00, anomaly coded in 3/18.
120  ! Local variables
121  real :: al0, p0, lambda
122  real :: al_TS, p_TSp, lam_TS, pa_000
123  integer :: j
124 
125  if (present(rho_ref)) pa_000 = (b0*(1.0 - a0*rho_ref) - rho_ref*c0)
126  if (present(rho_ref)) then ; do j=start,start+npts-1
127  al_ts = a1*t(j) +a2*s(j)
128  al0 = a0 + al_ts
129  p_tsp = pressure(j) + (b4*s(j) + t(j) * (b1 + (t(j)*(b2 + b3*t(j)) + b5*s(j))))
130  lam_ts = c4*s(j) + t(j) * (c1 + (t(j)*(c2 + c3*t(j)) + c5*s(j)))
131 
132  ! The following two expressions are mathematically equivalent.
133  ! rho(j) = (b0 + p0_TSp) / ((c0 + lam_TS) + al0*(b0 + p0_TSp)) - rho_ref
134  rho(j) = (pa_000 + (p_tsp - rho_ref*(p_tsp*al0 + (b0*al_ts + lam_ts)))) / &
135  ( (c0 + lam_ts) + al0*(b0 + p_tsp) )
136  enddo ; else ; do j=start,start+npts-1
137  al0 = (a0 + a1*t(j)) +a2*s(j)
138  p0 = (b0 + b4*s(j)) + t(j) * (b1 + t(j)*(b2 + b3*t(j)) + b5*s(j))
139  lambda = (c0 +c4*s(j)) + t(j) * (c1 + t(j)*(c2 + c3*t(j)) + c5*s(j))
140  rho(j) = (pressure(j) + p0) / (lambda + al0*(pressure(j) + p0))
141  enddo ; endif
142 
143 end subroutine calculate_density_array_wright
144 
145 !> This subroutine computes the in situ specific volume of sea water (specvol in
146 !! [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC])
147 !! and pressure [Pa]. It uses the expression from
148 !! Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
149 !! If spv_ref is present, specvol is an anomaly from spv_ref.
150 subroutine calculate_spec_vol_scalar_wright(T, S, pressure, specvol, spv_ref)
151  real, intent(in) :: T !< potential temperature relative to the surface [degC].
152  real, intent(in) :: S !< salinity [PSU].
153  real, intent(in) :: pressure !< pressure [Pa].
154  real, intent(out) :: specvol !< in situ specific volume [m3 kg-1].
155  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1].
156 
157  ! Local variables
158  real, dimension(1) :: T0, S0, pressure0, spv0
159 
160  t0(1) = t ; s0(1) = s ; pressure0(1) = pressure
161 
162  call calculate_spec_vol_array_wright(t0, s0, pressure0, spv0, 1, 1, spv_ref)
163  specvol = spv0(1)
165 
166 !> This subroutine computes the in situ specific volume of sea water (specvol in
167 !! [m3 kg-1]) from salinity (S [PSU]), potential temperature (T [degC])
168 !! and pressure [Pa]. It uses the expression from
169 !! Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
170 !! If spv_ref is present, specvol is an anomaly from spv_ref.
171 subroutine calculate_spec_vol_array_wright(T, S, pressure, specvol, start, npts, spv_ref)
172  real, dimension(:), intent(in) :: T !< potential temperature relative to the
173  !! surface [degC].
174  real, dimension(:), intent(in) :: S !< salinity [PSU].
175  real, dimension(:), intent(in) :: pressure !< pressure [Pa].
176  real, dimension(:), intent(out) :: specvol !< in situ specific volume [m3 kg-1].
177  integer, intent(in) :: start !< the starting point in the arrays.
178  integer, intent(in) :: npts !< the number of values to calculate.
179  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1].
180 
181  ! Local variables
182  real :: al0, p0, lambda
183  integer :: j
184 
185  do j=start,start+npts-1
186  al0 = (a0 + a1*t(j)) +a2*s(j)
187  p0 = (b0 + b4*s(j)) + t(j) * (b1 + t(j)*((b2 + b3*t(j))) + b5*s(j))
188  lambda = (c0 +c4*s(j)) + t(j) * (c1 + t(j)*((c2 + c3*t(j))) + c5*s(j))
189 
190  if (present(spv_ref)) then
191  specvol(j) = (lambda + (al0 - spv_ref)*(pressure(j) + p0)) / (pressure(j) + p0)
192  else
193  specvol(j) = (lambda + al0*(pressure(j) + p0)) / (pressure(j) + p0)
194  endif
195  enddo
196 end subroutine calculate_spec_vol_array_wright
197 
198 !> For a given thermodynamic state, return the thermal/haline expansion coefficients
199 subroutine calculate_density_derivs_array_wright(T, S, pressure, drho_dT, drho_dS, start, npts)
200  real, intent(in), dimension(:) :: T !< Potential temperature relative to the
201  !! surface [degC].
202  real, intent(in), dimension(:) :: S !< Salinity [PSU].
203  real, intent(in), dimension(:) :: pressure !< pressure [Pa].
204  real, intent(out), dimension(:) :: drho_dT !< The partial derivative of density with potential
205  !! temperature [kg m-3 degC-1].
206  real, intent(out), dimension(:) :: drho_dS !< The partial derivative of density with salinity,
207  !! in [kg m-3 PSU-1].
208  integer, intent(in) :: start !< The starting point in the arrays.
209  integer, intent(in) :: npts !< The number of values to calculate.
210 
211  ! Local variables
212  real :: al0, p0, lambda, I_denom2
213  integer :: j
214 
215  do j=start,start+npts-1
216  al0 = (a0 + a1*t(j)) + a2*s(j)
217  p0 = (b0 + b4*s(j)) + t(j) * (b1 + t(j)*((b2 + b3*t(j))) + b5*s(j))
218  lambda = (c0 +c4*s(j)) + t(j) * (c1 + t(j)*((c2 + c3*t(j))) + c5*s(j))
219 
220  i_denom2 = 1.0 / (lambda + al0*(pressure(j) + p0))
221  i_denom2 = i_denom2 *i_denom2
222  drho_dt(j) = i_denom2 * &
223  (lambda* (b1 + t(j)*(2.0*b2 + 3.0*b3*t(j)) + b5*s(j)) - &
224  (pressure(j)+p0) * ( (pressure(j)+p0)*a1 + &
225  (c1 + t(j)*(c2*2.0 + c3*3.0*t(j)) + c5*s(j)) ))
226  drho_ds(j) = i_denom2 * (lambda* (b4 + b5*t(j)) - &
227  (pressure(j)+p0) * ( (pressure(j)+p0)*a2 + (c4 + c5*t(j)) ))
228  enddo
229 
231 
232 !> The scalar version of calculate_density_derivs which promotes scalar inputs to a 1-element array and then
233 !! demotes the output back to a scalar
234 subroutine calculate_density_derivs_scalar_wright(T, S, pressure, drho_dT, drho_dS)
235  real, intent(in) :: T !< Potential temperature relative to the surface [degC].
236  real, intent(in) :: S !< Salinity [PSU].
237  real, intent(in) :: pressure !< pressure [Pa].
238  real, intent(out) :: drho_dT !< The partial derivative of density with potential
239  !! temperature [kg m-3 degC-1].
240  real, intent(out) :: drho_dS !< The partial derivative of density with salinity,
241  !! in [kg m-3 PSU-1].
242 
243  ! Local variables needed to promote the input/output scalars to 1-element arrays
244  real, dimension(1) :: T0, S0, P0
245  real, dimension(1) :: drdt0, drds0
246 
247  t0(1) = t
248  s0(1) = s
249  p0(1) = pressure
250  call calculate_density_derivs_array_wright(t0, s0, p0, drdt0, drds0, 1, 1)
251  drho_dt = drdt0(1)
252  drho_ds = drds0(1)
253 
255 
256 !> Second derivatives of density with respect to temperature, salinity, and pressure
257 subroutine calculate_density_second_derivs_array_wright(T, S, P, drho_ds_ds, drho_ds_dt, drho_dt_dt, &
258  drho_ds_dp, drho_dt_dp, start, npts)
259  real, dimension(:), intent(in ) :: T !< Potential temperature referenced to 0 dbar [degC]
260  real, dimension(:), intent(in ) :: S !< Salinity [PSU]
261  real, dimension(:), intent(in ) :: P !< Pressure [Pa]
262  real, dimension(:), intent( out) :: drho_ds_ds !< Partial derivative of beta with respect
263  !! to S [kg m-3 PSU-2]
264  real, dimension(:), intent( out) :: drho_ds_dt !< Partial derivative of beta with respcct
265  !! to T [kg m-3 PSU-1 degC-1]
266  real, dimension(:), intent( out) :: drho_dt_dt !< Partial derivative of alpha with respect
267  !! to T [kg m-3 degC-2]
268  real, dimension(:), intent( out) :: drho_ds_dp !< Partial derivative of beta with respect
269  !! to pressure [kg m-3 PSU-1 Pa-1]
270  real, dimension(:), intent( out) :: drho_dt_dp !< Partial derivative of alpha with respect
271  !! to pressure [kg m-3 degC-1 Pa-1]
272  integer, intent(in ) :: start !< Starting index in T,S,P
273  integer, intent(in ) :: npts !< Number of points to loop over
274 
275  ! Local variables
276  real :: z0, z1, z2, z3, z4, z5, z6 ,z7, z8, z9, z10, z11, z2_2, z2_3
277  integer :: j
278  ! Based on the above expression with common terms factored, there probably exists a more numerically stable
279  ! and/or efficient expression
280 
281  do j = start,start+npts-1
282  z0 = t(j)*(b1 + b5*s(j) + t(j)*(b2 + b3*t(j)))
283  z1 = (b0 + p(j) + b4*s(j) + z0)
284  z3 = (b1 + b5*s(j) + t(j)*(2.*b2 + 2.*b3*t(j)))
285  z4 = (c0 + c4*s(j) + t(j)*(c1 + c5*s(j) + t(j)*(c2 + c3*t(j))))
286  z5 = (b1 + b5*s(j) + t(j)*(b2 + b3*t(j)) + t(j)*(b2 + 2.*b3*t(j)))
287  z6 = c1 + c5*s(j) + t(j)*(c2 + c3*t(j)) + t(j)*(c2 + 2.*c3*t(j))
288  z7 = (c4 + c5*t(j) + a2*z1)
289  z8 = (c1 + c5*s(j) + t(j)*(2.*c2 + 3.*c3*t(j)) + a1*z1)
290  z9 = (a0 + a2*s(j) + a1*t(j))
291  z10 = (b4 + b5*t(j))
292  z11 = (z10*z4 - z1*z7)
293  z2 = (c0 + c4*s(j) + t(j)*(c1 + c5*s(j) + t(j)*(c2 + c3*t(j))) + z9*z1)
294  z2_2 = z2*z2
295  z2_3 = z2_2*z2
296 
297  drho_ds_ds(j) = (z10*(c4 + c5*t(j)) - a2*z10*z1 - z10*z7)/z2_2 - (2.*(c4 + c5*t(j) + z9*z10 + a2*z1)*z11)/z2_3
298  drho_ds_dt(j) = (z10*z6 - z1*(c5 + a2*z5) + b5*z4 - z5*z7)/z2_2 - (2.*(z6 + z9*z5 + a1*z1)*z11)/z2_3
299  drho_dt_dt(j) = (z3*z6 - z1*(2.*c2 + 6.*c3*t(j) + a1*z5) + (2.*b2 + 4.*b3*t(j))*z4 - z5*z8)/z2_2 - &
300  (2.*(z6 + z9*z5 + a1*z1)*(z3*z4 - z1*z8))/z2_3
301  drho_ds_dp(j) = (-c4 - c5*t(j) - 2.*a2*z1)/z2_2 - (2.*z9*z11)/z2_3
302  drho_dt_dp(j) = (-c1 - c5*s(j) - t(j)*(2.*c2 + 3.*c3*t(j)) - 2.*a1*z1)/z2_2 - (2.*z9*(z3*z4 - z1*z8))/z2_3
303  enddo
304 
306 
307 !> Second derivatives of density with respect to temperature, salinity, and pressure for scalar inputs. Inputs
308 !! promoted to 1-element array and output demoted to scalar
309 subroutine calculate_density_second_derivs_scalar_wright(T, S, P, drho_ds_ds, drho_ds_dt, drho_dt_dt, &
310  drho_ds_dp, drho_dt_dp)
311  real, intent(in ) :: T !< Potential temperature referenced to 0 dbar
312  real, intent(in ) :: S !< Salinity [PSU]
313  real, intent(in ) :: P !< pressure [Pa]
314  real, intent( out) :: drho_ds_ds !< Partial derivative of beta with respect
315  !! to S [kg m-3 PSU-2]
316  real, intent( out) :: drho_ds_dt !< Partial derivative of beta with respcct
317  !! to T [kg m-3 PSU-1 degC-1]
318  real, intent( out) :: drho_dt_dt !< Partial derivative of alpha with respect
319  !! to T [kg m-3 degC-2]
320  real, intent( out) :: drho_ds_dp !< Partial derivative of beta with respect
321  !! to pressure [kg m-3 PSU-1 Pa-1]
322  real, intent( out) :: drho_dt_dp !< Partial derivative of alpha with respect
323  !! to pressure [kg m-3 degC-1 Pa-1]
324  ! Local variables
325  real, dimension(1) :: T0, S0, P0
326  real, dimension(1) :: drdsds, drdsdt, drdtdt, drdsdp, drdtdp
327 
328  t0(1) = t
329  s0(1) = s
330  p0(1) = p
331  call calculate_density_second_derivs_array_wright(t0, s0, p0, drdsds, drdsdt, drdtdt, drdsdp, drdtdp, 1, 1)
332  drho_ds_ds = drdsds(1)
333  drho_ds_dt = drdsdt(1)
334  drho_dt_dt = drdtdt(1)
335  drho_ds_dp = drdsdp(1)
336  drho_dt_dp = drdtdp(1)
337 
339 
340 !> For a given thermodynamic state, return the partial derivatives of specific volume
341 !! with temperature and salinity
342 subroutine calculate_specvol_derivs_wright(T, S, pressure, dSV_dT, dSV_dS, start, npts)
343  real, intent(in), dimension(:) :: t !< Potential temperature relative to the surface [degC].
344  real, intent(in), dimension(:) :: s !< Salinity [PSU].
345  real, intent(in), dimension(:) :: pressure !< pressure [Pa].
346  real, intent(out), dimension(:) :: dsv_dt !< The partial derivative of specific volume with
347  !! potential temperature [m3 kg-1 degC-1].
348  real, intent(out), dimension(:) :: dsv_ds !< The partial derivative of specific volume with
349  !! salinity [m3 kg-1 / Pa].
350  integer, intent(in) :: start !< The starting point in the arrays.
351  integer, intent(in) :: npts !< The number of values to calculate.
352 
353  ! Local variables
354  real :: al0, p0, lambda, i_denom
355  integer :: j
356 
357  do j=start,start+npts-1
358 ! al0 = (a0 + a1*T(j)) + a2*S(j)
359  p0 = (b0 + b4*s(j)) + t(j) * (b1 + t(j)*((b2 + b3*t(j))) + b5*s(j))
360  lambda = (c0 +c4*s(j)) + t(j) * (c1 + t(j)*((c2 + c3*t(j))) + c5*s(j))
361 
362  ! SV = al0 + lambda / (pressure(j) + p0)
363 
364  i_denom = 1.0 / (pressure(j) + p0)
365  dsv_dt(j) = (a1 + i_denom * (c1 + t(j)*((2.0*c2 + 3.0*c3*t(j))) + c5*s(j))) - &
366  (i_denom**2 * lambda) * (b1 + t(j)*((2.0*b2 + 3.0*b3*t(j))) + b5*s(j))
367  dsv_ds(j) = (a2 + i_denom * (c4 + c5*t(j))) - &
368  (i_denom**2 * lambda) * (b4 + b5*t(j))
369  enddo
370 
371 end subroutine calculate_specvol_derivs_wright
372 
373 !> This subroutine computes the in situ density of sea water (rho in
374 !! [kg m-3]) and the compressibility (drho/dp = C_sound^-2)
375 !! (drho_dp [s2 m-2]) from salinity (sal in psu), potential
376 !! temperature (T [degC]), and pressure [Pa]. It uses the expressions
377 !! from Wright, 1997, J. Atmos. Ocean. Tech., 14, 735-740.
378 !! Coded by R. Hallberg, 1/01
379 subroutine calculate_compress_wright(T, S, pressure, rho, drho_dp, start, npts)
380  real, intent(in), dimension(:) :: t !< Potential temperature relative to the surface [degC].
381  real, intent(in), dimension(:) :: s !< Salinity [PSU].
382  real, intent(in), dimension(:) :: pressure !< pressure [Pa].
383  real, intent(out), dimension(:) :: rho !< In situ density [kg m-3].
384  real, intent(out), dimension(:) :: drho_dp !< The partial derivative of density with pressure
385  !! (also the inverse of the square of sound speed)
386  !! [s2 m-2].
387  integer, intent(in) :: start !< The starting point in the arrays.
388  integer, intent(in) :: npts !< The number of values to calculate.
389 
390  ! Coded by R. Hallberg, 1/01
391  ! Local variables
392  real :: al0, p0, lambda, i_denom
393  integer :: j
394 
395  do j=start,start+npts-1
396  al0 = (a0 + a1*t(j)) +a2*s(j)
397  p0 = (b0 + b4*s(j)) + t(j) * (b1 + t(j)*((b2 + b3*t(j))) + b5*s(j))
398  lambda = (c0 +c4*s(j)) + t(j) * (c1 + t(j)*((c2 + c3*t(j))) + c5*s(j))
399 
400  i_denom = 1.0 / (lambda + al0*(pressure(j) + p0))
401  rho(j) = (pressure(j) + p0) * i_denom
402  drho_dp(j) = lambda * i_denom * i_denom
403  enddo
404 end subroutine calculate_compress_wright
405 
406 !> This subroutine calculates analytical and nearly-analytical integrals of
407 !! pressure anomalies across layers, which are required for calculating the
408 !! finite-volume form pressure accelerations in a Boussinesq model.
409 subroutine int_density_dz_wright(T, S, z_t, z_b, rho_ref, rho_0, G_e, HII, HIO, &
410  dpa, intz_dpa, intx_dpa, inty_dpa, &
411  bathyT, dz_neglect, useMassWghtInterp)
412  type(hor_index_type), intent(in) :: hii !< The horizontal index type for the input arrays.
413  type(hor_index_type), intent(in) :: hio !< The horizontal index type for the output arrays.
414  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed), &
415  intent(in) :: t !< Potential temperature relative to the surface
416  !! [degC].
417  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed), &
418  intent(in) :: s !< Salinity [PSU].
419  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed), &
420  intent(in) :: z_t !< Height at the top of the layer in depth units [Z ~> m].
421  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed), &
422  intent(in) :: z_b !< Height at the top of the layer [Z ~> m].
423  real, intent(in) :: rho_ref !< A mean density [kg m-3], that is subtracted out
424  !! to reduce the magnitude of each of the integrals.
425  !! (The pressure is calucated as p~=-z*rho_0*G_e.)
426  real, intent(in) :: rho_0 !< Density [kg m-3], that is used to calculate the
427  !! pressure (as p~=-z*rho_0*G_e) used in the
428  !! equation of state.
429  real, intent(in) :: g_e !< The Earth's gravitational acceleration [m2 Z-1 s-2 ~> m s-2].
430  real, dimension(HIO%isd:HIO%ied,HIO%jsd:HIO%jed), &
431  intent(out) :: dpa !< The change in the pressure anomaly across the
432  !! layer [Pa].
433  real, dimension(HIO%isd:HIO%ied,HIO%jsd:HIO%jed), &
434  optional, intent(out) :: intz_dpa !< The integral through the thickness of the layer
435  !! of the pressure anomaly relative to the anomaly
436  !! at the top of the layer [Pa Z ~> Pa m].
437  real, dimension(HIO%IsdB:HIO%IedB,HIO%jsd:HIO%jed), &
438  optional, intent(out) :: intx_dpa !< The integral in x of the difference between the
439  !! pressure anomaly at the top and bottom of the
440  !! layer divided by the x grid spacing [Pa].
441  real, dimension(HIO%isd:HIO%ied,HIO%JsdB:HIO%JedB), &
442  optional, intent(out) :: inty_dpa !< The integral in y of the difference between the
443  !! pressure anomaly at the top and bottom of the
444  !! layer divided by the y grid spacing [Pa].
445  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed), &
446  optional, intent(in) :: bathyt !< The depth of the bathymetry [Z ~> m].
447  real, optional, intent(in) :: dz_neglect !< A miniscule thickness change [Z ~> m].
448  logical, optional, intent(in) :: usemasswghtinterp !< If true, uses mass weighting to
449  !! interpolate T/S for top and bottom integrals.
450 
451  ! Local variables
452  real, dimension(HII%isd:HII%ied,HII%jsd:HII%jed) :: al0_2d, p0_2d, lambda_2d
453  real :: al0, p0, lambda
454  real :: rho_anom ! The density anomaly from rho_ref [kg m-3].
455  real :: eps, eps2, rem
456  real :: gxrho, i_rho
457  real :: p_ave, i_al0, i_lzz
458  real :: dz ! The layer thickness [Z ~> m].
459  real :: hwght ! A pressure-thickness below topography [Z ~> m].
460  real :: hl, hr ! Pressure-thicknesses of the columns to the left and right [Z ~> m].
461  real :: idenom ! The inverse of the denominator in the weights [Z-2 ~> m-2].
462  real :: hwt_ll, hwt_lr ! hWt_LA is the weighted influence of A on the left column [nondim].
463  real :: hwt_rl, hwt_rr ! hWt_RA is the weighted influence of A on the right column [nondim].
464  real :: wt_l, wt_r ! The linear weights of the left and right columns [nondim].
465  real :: wtt_l, wtt_r ! The weights for tracers from the left and right columns [nondim].
466  real :: intz(5) ! The integrals of density with height at the
467  ! 5 sub-column locations [Pa].
468  logical :: do_massweight ! Indicates whether to do mass weighting.
469  real, parameter :: c1_3 = 1.0/3.0, c1_7 = 1.0/7.0 ! Rational constants.
470  real, parameter :: c1_9 = 1.0/9.0, c1_90 = 1.0/90.0 ! Rational constants.
471  integer :: is, ie, js, je, isq, ieq, jsq, jeq, i, j, ioff, joff, m
472 
473  ioff = hio%idg_offset - hii%idg_offset
474  joff = hio%jdg_offset - hii%jdg_offset
475 
476  ! These array bounds work for the indexing convention of the input arrays, but
477  ! on the computational domain defined for the output arrays.
478  isq = hio%IscB + ioff ; ieq = hio%IecB + ioff
479  jsq = hio%JscB + joff ; jeq = hio%JecB + joff
480  is = hio%isc + ioff ; ie = hio%iec + ioff
481  js = hio%jsc + joff ; je = hio%jec + joff
482 
483  gxrho = g_e * rho_0
484  i_rho = 1.0 / rho_0
485 
486  do_massweight = .false.
487  if (present(usemasswghtinterp)) then ; if (usemasswghtinterp) then
488  do_massweight = .true.
489  ! if (.not.present(bathyT)) call MOM_error(FATAL, "int_density_dz_generic: "//&
490  ! "bathyT must be present if useMassWghtInterp is present and true.")
491  ! if (.not.present(dz_neglect)) call MOM_error(FATAL, "int_density_dz_generic: "//&
492  ! "dz_neglect must be present if useMassWghtInterp is present and true.")
493  endif ; endif
494 
495  do j=jsq,jeq+1 ; do i=isq,ieq+1
496  al0_2d(i,j) = (a0 + a1*t(i,j)) + a2*s(i,j)
497  p0_2d(i,j) = (b0 + b4*s(i,j)) + t(i,j) * (b1 + t(i,j)*((b2 + b3*t(i,j))) + b5*s(i,j))
498  lambda_2d(i,j) = (c0 +c4*s(i,j)) + t(i,j) * (c1 + t(i,j)*((c2 + c3*t(i,j))) + c5*s(i,j))
499 
500  al0 = al0_2d(i,j) ; p0 = p0_2d(i,j) ; lambda = lambda_2d(i,j)
501 
502  dz = z_t(i,j) - z_b(i,j)
503  p_ave = -0.5*gxrho*(z_t(i,j)+z_b(i,j))
504 
505  i_al0 = 1.0 / al0
506  i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
507  eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
508 
509 ! rho(j) = (pressure(j) + p0) / (lambda + al0*(pressure(j) + p0))
510 
511  rho_anom = (p0 + p_ave)*(i_lzz*i_al0) - rho_ref
512  rem = i_rho * (lambda * i_al0**2) * eps2 * &
513  (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
514  dpa(i-ioff,j-joff) = g_e*rho_anom*dz - 2.0*eps*rem
515  if (present(intz_dpa)) &
516  intz_dpa(i-ioff,j-joff) = 0.5*g_e*rho_anom*dz**2 - dz*(1.0+eps)*rem
517  enddo ; enddo
518 
519  if (present(intx_dpa)) then ; do j=js,je ; do i=isq,ieq
520  ! hWght is the distance measure by which the cell is violation of
521  ! hydrostatic consistency. For large hWght we bias the interpolation of
522  ! T & S along the top and bottom integrals, akin to thickness weighting.
523  hwght = 0.0
524  if (do_massweight) &
525  hwght = max(0., -bathyt(i,j)-z_t(i+1,j), -bathyt(i+1,j)-z_t(i,j))
526  if (hwght > 0.) then
527  hl = (z_t(i,j) - z_b(i,j)) + dz_neglect
528  hr = (z_t(i+1,j) - z_b(i+1,j)) + dz_neglect
529  hwght = hwght * ( (hl-hr)/(hl+hr) )**2
530  idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
531  hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
532  hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
533  else
534  hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
535  endif
536 
537  intz(1) = dpa(i-ioff,j-joff) ; intz(5) = dpa(i+1-ioff,j-joff)
538  do m=2,4
539  wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
540  wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
541 
542  al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i+1,j)
543  p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i+1,j)
544  lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i+1,j)
545 
546  dz = wt_l*(z_t(i,j) - z_b(i,j)) + wt_r*(z_t(i+1,j) - z_b(i+1,j))
547  p_ave = -0.5*gxrho*(wt_l*(z_t(i,j)+z_b(i,j)) + &
548  wt_r*(z_t(i+1,j)+z_b(i+1,j)))
549 
550  i_al0 = 1.0 / al0
551  i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
552  eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
553 
554  intz(m) = g_e*dz*((p0 + p_ave)*(i_lzz*i_al0) - rho_ref) - 2.0*eps * &
555  i_rho * (lambda * i_al0**2) * eps2 * &
556  (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
557  enddo
558  ! Use Bode's rule to integrate the values.
559  intx_dpa(i-ioff,j-joff) = c1_90*(7.0*(intz(1)+intz(5)) + 32.0*(intz(2)+intz(4)) + &
560  12.0*intz(3))
561  enddo ; enddo ; endif
562 
563  if (present(inty_dpa)) then ; do j=jsq,jeq ; do i=is,ie
564  ! hWght is the distance measure by which the cell is violation of
565  ! hydrostatic consistency. For large hWght we bias the interpolation of
566  ! T & S along the top and bottom integrals, akin to thickness weighting.
567  hwght = 0.0
568  if (do_massweight) &
569  hwght = max(0., -bathyt(i,j)-z_t(i,j+1), -bathyt(i,j+1)-z_t(i,j))
570  if (hwght > 0.) then
571  hl = (z_t(i,j) - z_b(i,j)) + dz_neglect
572  hr = (z_t(i,j+1) - z_b(i,j+1)) + dz_neglect
573  hwght = hwght * ( (hl-hr)/(hl+hr) )**2
574  idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
575  hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
576  hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
577  else
578  hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
579  endif
580 
581  intz(1) = dpa(i-ioff,j-joff) ; intz(5) = dpa(i-ioff,j+1-joff)
582  do m=2,4
583  wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
584  wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
585 
586  al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i,j+1)
587  p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i,j+1)
588  lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i,j+1)
589 
590  dz = wt_l*(z_t(i,j) - z_b(i,j)) + wt_r*(z_t(i,j+1) - z_b(i,j+1))
591  p_ave = -0.5*gxrho*(wt_l*(z_t(i,j)+z_b(i,j)) + &
592  wt_r*(z_t(i,j+1)+z_b(i,j+1)))
593 
594  i_al0 = 1.0 / al0
595  i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
596  eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
597 
598  intz(m) = g_e*dz*((p0 + p_ave)*(i_lzz*i_al0) - rho_ref) - 2.0*eps * &
599  i_rho * (lambda * i_al0**2) * eps2 * &
600  (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
601  enddo
602  ! Use Bode's rule to integrate the values.
603  inty_dpa(i-ioff,j-joff) = c1_90*(7.0*(intz(1)+intz(5)) + 32.0*(intz(2)+intz(4)) + &
604  12.0*intz(3))
605  enddo ; enddo ; endif
606 end subroutine int_density_dz_wright
607 
608 !> This subroutine calculates analytical and nearly-analytical integrals in
609 !! pressure across layers of geopotential anomalies, which are required for
610 !! calculating the finite-volume form pressure accelerations in a non-Boussinesq
611 !! model. There are essentially no free assumptions, apart from the use of
612 !! Bode's rule to do the horizontal integrals, and from a truncation in the
613 !! series for log(1-eps/1+eps) that assumes that |eps| < 0.34.
614 subroutine int_spec_vol_dp_wright(T, S, p_t, p_b, spv_ref, HI, dza, &
615  intp_dza, intx_dza, inty_dza, halo_size, &
616  bathyP, dP_neglect, useMassWghtInterp)
617  type(hor_index_type), intent(in) :: hi !< The ocean's horizontal index type.
618  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
619  intent(in) :: t !< Potential temperature relative to the surface
620  !! [degC].
621  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
622  intent(in) :: s !< Salinity [PSU].
623  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
624  intent(in) :: p_t !< Pressure at the top of the layer [Pa].
625  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
626  intent(in) :: p_b !< Pressure at the top of the layer [Pa].
627  real, intent(in) :: spv_ref !< A mean specific volume that is subtracted out
628  !! to reduce the magnitude of each of the integrals [m3 kg-1]. The calculation is
629  !! mathematically identical with different values of spv_ref, but this reduces the
630  !! effects of roundoff.
631  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
632  intent(out) :: dza !< The change in the geopotential anomaly across
633  !! the layer [m2 s-2].
634  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
635  optional, intent(out) :: intp_dza !< The integral in pressure through the layer of
636  !! the geopotential anomaly relative to the anomaly
637  !! at the bottom of the layer [Pa m2 s-2].
638  real, dimension(HI%IsdB:HI%IedB,HI%jsd:HI%jed), &
639  optional, intent(out) :: intx_dza !< The integral in x of the difference between the
640  !! geopotential anomaly at the top and bottom of
641  !! the layer divided by the x grid spacing [m2 s-2].
642  real, dimension(HI%isd:HI%ied,HI%JsdB:HI%JedB), &
643  optional, intent(out) :: inty_dza !< The integral in y of the difference between the
644  !! geopotential anomaly at the top and bottom of
645  !! the layer divided by the y grid spacing [m2 s-2].
646  integer, optional, intent(in) :: halo_size !< The width of halo points on which to calculate
647  !! dza.
648  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
649  optional, intent(in) :: bathyp !< The pressure at the bathymetry [Pa]
650  real, optional, intent(in) :: dp_neglect !< A miniscule pressure change with
651  !! the same units as p_t [Pa]
652  logical, optional, intent(in) :: usemasswghtinterp !< If true, uses mass weighting
653  !! to interpolate T/S for top and bottom integrals.
654 
655  ! Local variables
656  real, dimension(HI%isd:HI%ied,HI%jsd:HI%jed) :: al0_2d, p0_2d, lambda_2d
657  real :: al0, p0, lambda
658  real :: p_ave
659  real :: rem, eps, eps2
660  real :: alpha_anom ! The depth averaged specific volume anomaly [m3 kg-1].
661  real :: dp ! The pressure change through a layer [Pa].
662  real :: hwght ! A pressure-thickness below topography [Pa].
663  real :: hl, hr ! Pressure-thicknesses of the columns to the left and right [Pa].
664  real :: idenom ! The inverse of the denominator in the weights [Pa-2].
665  real :: hwt_ll, hwt_lr ! hWt_LA is the weighted influence of A on the left column [nondim].
666  real :: hwt_rl, hwt_rr ! hWt_RA is the weighted influence of A on the right column [nondim].
667  real :: wt_l, wt_r ! The linear weights of the left and right columns [nondim].
668  real :: wtt_l, wtt_r ! The weights for tracers from the left and right columns [nondim].
669  real :: intp(5) ! The integrals of specific volume with pressure at the
670  ! 5 sub-column locations [m2 s-2].
671  logical :: do_massweight ! Indicates whether to do mass weighting.
672  real, parameter :: c1_3 = 1.0/3.0, c1_7 = 1.0/7.0 ! Rational constants.
673  real, parameter :: c1_9 = 1.0/9.0, c1_90 = 1.0/90.0 ! Rational constants.
674  integer :: isq, ieq, jsq, jeq, ish, ieh, jsh, jeh, i, j, m, halo
675 
676  isq = hi%IscB ; ieq = hi%IecB ; jsq = hi%JscB ; jeq = hi%JecB
677  halo = 0 ; if (present(halo_size)) halo = max(halo_size,0)
678  ish = hi%isc-halo ; ieh = hi%iec+halo ; jsh = hi%jsc-halo ; jeh = hi%jec+halo
679  if (present(intx_dza)) then ; ish = min(isq,ish) ; ieh = max(ieq+1,ieh); endif
680  if (present(inty_dza)) then ; jsh = min(jsq,jsh) ; jeh = max(jeq+1,jeh); endif
681 
682  do_massweight = .false.
683  if (present(usemasswghtinterp)) then ; if (usemasswghtinterp) then
684  do_massweight = .true.
685 ! if (.not.present(bathyP)) call MOM_error(FATAL, "int_spec_vol_dp_generic: "//&
686 ! "bathyP must be present if useMassWghtInterp is present and true.")
687 ! if (.not.present(dP_neglect)) call MOM_error(FATAL, "int_spec_vol_dp_generic: "//&
688 ! "dP_neglect must be present if useMassWghtInterp is present and true.")
689  endif ; endif
690 
691  do j=jsh,jeh ; do i=ish,ieh
692  al0_2d(i,j) = (a0 + a1*t(i,j)) + a2*s(i,j)
693  p0_2d(i,j) = (b0 + b4*s(i,j)) + t(i,j) * (b1 + t(i,j)*((b2 + b3*t(i,j))) + b5*s(i,j))
694  lambda_2d(i,j) = (c0 +c4*s(i,j)) + t(i,j) * (c1 + t(i,j)*((c2 + c3*t(i,j))) + c5*s(i,j))
695 
696  al0 = al0_2d(i,j) ; p0 = p0_2d(i,j) ; lambda = lambda_2d(i,j)
697  dp = p_b(i,j) - p_t(i,j)
698  p_ave = 0.5*(p_t(i,j)+p_b(i,j))
699 
700  eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
701  alpha_anom = al0 + lambda / (p0 + p_ave) - spv_ref
702  rem = lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
703  dza(i,j) = alpha_anom*dp + 2.0*eps*rem
704  if (present(intp_dza)) &
705  intp_dza(i,j) = 0.5*alpha_anom*dp**2 - dp*(1.0-eps)*rem
706  enddo ; enddo
707 
708  if (present(intx_dza)) then ; do j=hi%jsc,hi%jec ; do i=isq,ieq
709  ! hWght is the distance measure by which the cell is violation of
710  ! hydrostatic consistency. For large hWght we bias the interpolation of
711  ! T & S along the top and bottom integrals, akin to thickness weighting.
712  hwght = 0.0
713  if (do_massweight) &
714  hwght = max(0., bathyp(i,j)-p_t(i+1,j), bathyp(i+1,j)-p_t(i,j))
715  if (hwght > 0.) then
716  hl = (p_b(i,j) - p_t(i,j)) + dp_neglect
717  hr = (p_b(i+1,j) - p_t(i+1,j)) + dp_neglect
718  hwght = hwght * ( (hl-hr)/(hl+hr) )**2
719  idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
720  hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
721  hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
722  else
723  hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
724  endif
725 
726  intp(1) = dza(i,j) ; intp(5) = dza(i+1,j)
727  do m=2,4
728  wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
729  wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
730 
731  ! T, S, and p are interpolated in the horizontal. The p interpolation
732  ! is linear, but for T and S it may be thickness wekghted.
733  al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i+1,j)
734  p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i+1,j)
735  lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i+1,j)
736 
737  dp = wt_l*(p_b(i,j) - p_t(i,j)) + wt_r*(p_b(i+1,j) - p_t(i+1,j))
738  p_ave = 0.5*(wt_l*(p_t(i,j)+p_b(i,j)) + wt_r*(p_t(i+1,j)+p_b(i+1,j)))
739 
740  eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
741  intp(m) = (al0 + lambda / (p0 + p_ave) - spv_ref)*dp + 2.0*eps* &
742  lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
743  enddo
744  ! Use Bode's rule to integrate the values.
745  intx_dza(i,j) = c1_90*(7.0*(intp(1)+intp(5)) + 32.0*(intp(2)+intp(4)) + &
746  12.0*intp(3))
747  enddo ; enddo ; endif
748 
749  if (present(inty_dza)) then ; do j=jsq,jeq ; do i=hi%isc,hi%iec
750  ! hWght is the distance measure by which the cell is violation of
751  ! hydrostatic consistency. For large hWght we bias the interpolation of
752  ! T & S along the top and bottom integrals, akin to thickness weighting.
753  hwght = 0.0
754  if (do_massweight) &
755  hwght = max(0., bathyp(i,j)-p_t(i,j+1), bathyp(i,j+1)-p_t(i,j))
756  if (hwght > 0.) then
757  hl = (p_b(i,j) - p_t(i,j)) + dp_neglect
758  hr = (p_b(i,j+1) - p_t(i,j+1)) + dp_neglect
759  hwght = hwght * ( (hl-hr)/(hl+hr) )**2
760  idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
761  hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
762  hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
763  else
764  hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
765  endif
766 
767  intp(1) = dza(i,j) ; intp(5) = dza(i,j+1)
768  do m=2,4
769  wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
770  wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
771 
772  ! T, S, and p are interpolated in the horizontal. The p interpolation
773  ! is linear, but for T and S it may be thickness wekghted.
774  al0 = wt_l*al0_2d(i,j) + wt_r*al0_2d(i,j+1)
775  p0 = wt_l*p0_2d(i,j) + wt_r*p0_2d(i,j+1)
776  lambda = wt_l*lambda_2d(i,j) + wt_r*lambda_2d(i,j+1)
777 
778  dp = wt_l*(p_b(i,j) - p_t(i,j)) + wt_r*(p_b(i,j+1) - p_t(i,j+1))
779  p_ave = 0.5*(wt_l*(p_t(i,j)+p_b(i,j)) + wt_r*(p_t(i,j+1)+p_b(i,j+1)))
780 
781  eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
782  intp(m) = (al0 + lambda / (p0 + p_ave) - spv_ref)*dp + 2.0*eps* &
783  lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
784  enddo
785  ! Use Bode's rule to integrate the values.
786  inty_dza(i,j) = c1_90*(7.0*(intp(1)+intp(5)) + 32.0*(intp(2)+intp(4)) + &
787  12.0*intp(3))
788  enddo ; enddo ; endif
789 end subroutine int_spec_vol_dp_wright
790 
791 end module mom_eos_wright
mom_eos_wright::a2
real, parameter a2
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:67
mom_eos_wright::calculate_specvol_derivs_wright
subroutine, public calculate_specvol_derivs_wright(T, S, pressure, dSV_dT, dSV_dS, start, npts)
For a given thermodynamic state, return the partial derivatives of specific volume with temperature a...
Definition: MOM_EOS_Wright.F90:343
mom_eos_wright::b0
real, parameter b0
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:68
mom_eos_wright::calculate_density_derivs_array_wright
subroutine calculate_density_derivs_array_wright(T, S, pressure, drho_dT, drho_dS, start, npts)
For a given thermodynamic state, return the thermal/haline expansion coefficients.
Definition: MOM_EOS_Wright.F90:200
mom_eos_wright::calculate_density_derivs_scalar_wright
subroutine calculate_density_derivs_scalar_wright(T, S, pressure, drho_dT, drho_dS)
The scalar version of calculate_density_derivs which promotes scalar inputs to a 1-element array and ...
Definition: MOM_EOS_Wright.F90:235
mom_eos_wright::calculate_density_second_derivs_array_wright
subroutine calculate_density_second_derivs_array_wright(T, S, P, drho_ds_ds, drho_ds_dt, drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
Second derivatives of density with respect to temperature, salinity, and pressure.
Definition: MOM_EOS_Wright.F90:259
mom_eos_wright::int_density_dz_wright
subroutine, public int_density_dz_wright(T, S, z_t, z_b, rho_ref, rho_0, G_e, HII, HIO, dpa, intz_dpa, intx_dpa, inty_dpa, bathyT, dz_neglect, useMassWghtInterp)
This subroutine calculates analytical and nearly-analytical integrals of pressure anomalies across la...
Definition: MOM_EOS_Wright.F90:412
mom_eos_wright::a1
real, parameter a1
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:67
mom_eos_wright::c2
real, parameter c2
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:70
mom_eos_wright::b1
real, parameter b1
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:68
mom_eos_wright::b5
real, parameter b5
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:69
mom_hor_index
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Definition: MOM_hor_index.F90:2
mom_eos_wright::c4
real, parameter c4
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:71
mom_eos_wright::calculate_density_second_derivs_scalar_wright
subroutine calculate_density_second_derivs_scalar_wright(T, S, P, drho_ds_ds, drho_ds_dt, drho_dt_dt, drho_ds_dp, drho_dt_dp)
Second derivatives of density with respect to temperature, salinity, and pressure for scalar inputs....
Definition: MOM_EOS_Wright.F90:311
mom_eos_wright
The equation of state using the Wright 1997 expressions.
Definition: MOM_EOS_Wright.F90:2
mom_eos_wright::calculate_compress_wright
subroutine, public calculate_compress_wright(T, S, pressure, rho, drho_dp, start, npts)
This subroutine computes the in situ density of sea water (rho in [kg m-3]) and the compressibility (...
Definition: MOM_EOS_Wright.F90:380
mom_eos_wright::c1
real, parameter c1
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:70
mom_eos_wright::calculate_density_array_wright
subroutine calculate_density_array_wright(T, S, pressure, rho, start, npts, rho_ref)
This subroutine computes the in situ density of sea water (rho in [kg m-3]) from salinity (S [PSU]),...
Definition: MOM_EOS_Wright.F90:111
mom_eos_wright::calculate_density_wright
Compute the in situ density of sea water (in [kg m-3]), or its anomaly with respect to a reference de...
Definition: MOM_EOS_Wright.F90:32
mom_eos_wright::a0
real, parameter a0
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:67
mom_eos_wright::b3
real, parameter b3
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:69
mom_eos_wright::c0
real, parameter c0
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:70
mom_eos_wright::calculate_spec_vol_array_wright
subroutine calculate_spec_vol_array_wright(T, S, pressure, specvol, start, npts, spv_ref)
This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinit...
Definition: MOM_EOS_Wright.F90:172
mom_hor_index::hor_index_type
Container for horizontal index ranges for data, computational and global domains.
Definition: MOM_hor_index.F90:15
mom_eos_wright::calculate_density_derivs_wright
For a given thermodynamic state, return the derivatives of density with temperature and salinity.
Definition: MOM_EOS_Wright.F90:44
mom_eos_wright::calculate_spec_vol_scalar_wright
subroutine calculate_spec_vol_scalar_wright(T, S, pressure, specvol, spv_ref)
This subroutine computes the in situ specific volume of sea water (specvol in [m3 kg-1]) from salinit...
Definition: MOM_EOS_Wright.F90:151
mom_eos_wright::calculate_density_scalar_wright
subroutine calculate_density_scalar_wright(T, S, pressure, rho, rho_ref)
This subroutine computes the in situ density of sea water (rho in [kg m-3]) from salinity (S [PSU]),...
Definition: MOM_EOS_Wright.F90:81
mom_eos_wright::b2
real, parameter b2
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:68
mom_eos_wright::calculate_spec_vol_wright
Compute the in situ specific volume of sea water (in [m3 kg-1]), or an anomaly with respect to a refe...
Definition: MOM_EOS_Wright.F90:39
mom_eos_wright::calculate_density_second_derivs_wright
For a given thermodynamic state, return the second derivatives of density with various combinations o...
Definition: MOM_EOS_Wright.F90:50
mom_eos_wright::b4
real, parameter b4
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:69
mom_eos_wright::int_spec_vol_dp_wright
subroutine, public int_spec_vol_dp_wright(T, S, p_t, p_b, spv_ref, HI, dza, intp_dza, intx_dza, inty_dza, halo_size, bathyP, dP_neglect, useMassWghtInterp)
This subroutine calculates analytical and nearly-analytical integrals in pressure across layers of ge...
Definition: MOM_EOS_Wright.F90:617
mom_eos_wright::c3
real, parameter c3
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:71
mom_eos_wright::c5
real, parameter c5
Parameters in the Wright equation of state.
Definition: MOM_EOS_Wright.F90:71