14 implicit none ;
private
16 #include <MOM_memory.h>
21 public int_density_dz_wright, int_spec_vol_dp_wright
33 module procedure calculate_density_scalar_wright, calculate_density_array_wright
40 module procedure calculate_spec_vol_scalar_wright, calculate_spec_vol_array_wright
45 module procedure calculate_density_derivs_scalar_wright, calculate_density_derivs_array_wright
51 module procedure calculate_density_second_derivs_scalar_wright, calculate_density_second_derivs_array_wright
67 real,
parameter :: a0 = 7.057924e-4, a1 = 3.480336e-7, a2 = -1.112733e-7
68 real,
parameter :: b0 = 5.790749e8, b1 = 3.516535e6, b2 = -4.002714e4
69 real,
parameter :: b3 = 2.084372e2, b4 = 5.944068e5, b5 = -9.643486e3
70 real,
parameter :: c0 = 1.704853e5, c1 = 7.904722e2, c2 = -7.984422
71 real,
parameter :: c3 = 5.140652e-2, c4 = -2.302158e2, c5 = -3.079464
80 subroutine calculate_density_scalar_wright(T, S, pressure, rho, rho_ref)
83 real,
intent(in) :: pressure
84 real,
intent(out) :: rho
85 real,
optional,
intent(in) :: rho_ref
95 real,
dimension(1) :: T0, S0, pressure0, rho0
99 pressure0(1) = pressure
101 call calculate_density_array_wright(t0, s0, pressure0, rho0, 1, 1, rho_ref)
104 end subroutine calculate_density_scalar_wright
110 subroutine calculate_density_array_wright(T, S, pressure, rho, start, npts, rho_ref)
111 real,
dimension(:),
intent(in) :: T
112 real,
dimension(:),
intent(in) :: S
113 real,
dimension(:),
intent(in) :: pressure
114 real,
dimension(:),
intent(inout) :: rho
115 integer,
intent(in) :: start
116 integer,
intent(in) :: npts
117 real,
optional,
intent(in) :: rho_ref
121 real :: al0, p0, lambda
122 real :: al_TS, p_TSp, lam_TS, pa_000
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)
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)))
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))
143 end subroutine calculate_density_array_wright
150 subroutine calculate_spec_vol_scalar_wright(T, S, pressure, specvol, spv_ref)
151 real,
intent(in) :: T
152 real,
intent(in) :: S
153 real,
intent(in) :: pressure
154 real,
intent(out) :: specvol
155 real,
optional,
intent(in) :: spv_ref
158 real,
dimension(1) :: T0, S0, pressure0, spv0
160 t0(1) = t ; s0(1) = s ; pressure0(1) = pressure
162 call calculate_spec_vol_array_wright(t0, s0, pressure0, spv0, 1, 1, spv_ref)
164 end subroutine calculate_spec_vol_scalar_wright
171 subroutine calculate_spec_vol_array_wright(T, S, pressure, specvol, start, npts, spv_ref)
172 real,
dimension(:),
intent(in) :: T
174 real,
dimension(:),
intent(in) :: S
175 real,
dimension(:),
intent(in) :: pressure
176 real,
dimension(:),
intent(inout) :: specvol
177 integer,
intent(in) :: start
178 integer,
intent(in) :: npts
179 real,
optional,
intent(in) :: spv_ref
182 real :: al0, p0, lambda
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))
190 if (
present(spv_ref))
then
191 specvol(j) = (lambda + (al0 - spv_ref)*(pressure(j) + p0)) / (pressure(j) + p0)
193 specvol(j) = (lambda + al0*(pressure(j) + p0)) / (pressure(j) + p0)
196 end subroutine calculate_spec_vol_array_wright
199 subroutine calculate_density_derivs_array_wright(T, S, pressure, drho_dT, drho_dS, start, npts)
200 real,
intent(in),
dimension(:) :: T
202 real,
intent(in),
dimension(:) :: S
203 real,
intent(in),
dimension(:) :: pressure
204 real,
intent(inout),
dimension(:) :: drho_dT
206 real,
intent(inout),
dimension(:) :: drho_dS
208 integer,
intent(in) :: start
209 integer,
intent(in) :: npts
212 real :: al0, p0, lambda, I_denom2
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))
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)) ))
230 end subroutine calculate_density_derivs_array_wright
234 subroutine calculate_density_derivs_scalar_wright(T, S, pressure, drho_dT, drho_dS)
235 real,
intent(in) :: T
236 real,
intent(in) :: S
237 real,
intent(in) :: pressure
238 real,
intent(out) :: drho_dT
240 real,
intent(out) :: drho_dS
244 real,
dimension(1) :: T0, S0, P0
245 real,
dimension(1) :: drdt0, drds0
250 call calculate_density_derivs_array_wright(t0, s0, p0, drdt0, drds0, 1, 1)
254 end subroutine calculate_density_derivs_scalar_wright
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
260 real,
dimension(:),
intent(in ) :: S
261 real,
dimension(:),
intent(in ) :: P
262 real,
dimension(:),
intent(inout) :: drho_ds_ds
264 real,
dimension(:),
intent(inout) :: drho_ds_dt
266 real,
dimension(:),
intent(inout) :: drho_dt_dt
268 real,
dimension(:),
intent(inout) :: drho_ds_dp
270 real,
dimension(:),
intent(inout) :: drho_dt_dp
272 integer,
intent(in ) :: start
273 integer,
intent(in ) :: npts
276 real :: z0, z1, z2, z3, z4, z5, z6 ,z7, z8, z9, z10, z11, z2_2, z2_3
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))
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)
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
305 end subroutine calculate_density_second_derivs_array_wright
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
312 real,
intent(in ) :: S
313 real,
intent(in ) :: P
314 real,
intent( out) :: drho_ds_ds
316 real,
intent( out) :: drho_ds_dt
318 real,
intent( out) :: drho_dt_dt
320 real,
intent( out) :: drho_ds_dp
322 real,
intent( out) :: drho_dt_dp
325 real,
dimension(1) :: T0, S0, P0
326 real,
dimension(1) :: drdsds, drdsdt, drdtdt, drdsdp, drdtdp
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)
338 end subroutine calculate_density_second_derivs_scalar_wright
342 subroutine calculate_specvol_derivs_wright(T, S, pressure, dSV_dT, dSV_dS, start, npts)
343 real,
intent(in),
dimension(:) :: t
344 real,
intent(in),
dimension(:) :: s
345 real,
intent(in),
dimension(:) :: pressure
346 real,
intent(inout),
dimension(:) :: dsv_dt
348 real,
intent(inout),
dimension(:) :: dsv_ds
350 integer,
intent(in) :: start
351 integer,
intent(in) :: npts
354 real :: al0, p0, lambda, i_denom
357 do j=start,start+npts-1
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))
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))
371 end subroutine calculate_specvol_derivs_wright
379 subroutine calculate_compress_wright(T, S, pressure, rho, drho_dp, start, npts)
380 real,
intent(in),
dimension(:) :: t
381 real,
intent(in),
dimension(:) :: s
382 real,
intent(in),
dimension(:) :: pressure
383 real,
intent(inout),
dimension(:) :: rho
384 real,
intent(inout),
dimension(:) :: drho_dp
387 integer,
intent(in) :: start
388 integer,
intent(in) :: npts
392 real :: al0, p0, lambda, i_denom
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))
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
404 end subroutine calculate_compress_wright
409 subroutine int_density_dz_wright(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, &
410 dpa, intz_dpa, intx_dpa, inty_dpa, &
411 bathyT, dz_neglect, useMassWghtInterp, rho_scale, pres_scale)
413 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
416 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
418 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
420 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
422 real,
intent(in) :: rho_ref
425 real,
intent(in) :: rho_0
428 real,
intent(in) :: g_e
430 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
433 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
434 optional,
intent(inout) :: intz_dpa
437 real,
dimension(HI%IsdB:HI%IedB,HI%jsd:HI%jed), &
438 optional,
intent(inout) :: intx_dpa
441 real,
dimension(HI%isd:HI%ied,HI%JsdB:HI%JedB), &
442 optional,
intent(inout) :: inty_dpa
445 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
446 optional,
intent(in) :: bathyt
447 real,
optional,
intent(in) :: dz_neglect
448 logical,
optional,
intent(in) :: usemasswghtinterp
450 real,
optional,
intent(in) :: rho_scale
452 real,
optional,
intent(in) :: pres_scale
456 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed) :: al0_2d, p0_2d, lambda_2d
457 real :: al0, p0, lambda
459 real :: eps, eps2, rem
464 real :: p_ave, i_al0, i_lzz
469 real :: hwt_ll, hwt_lr
470 real :: hwt_rl, hwt_rr
477 logical :: do_massweight
478 real,
parameter :: c1_3 = 1.0/3.0, c1_7 = 1.0/7.0
479 real,
parameter :: c1_9 = 1.0/9.0, c1_90 = 1.0/90.0
480 integer :: is, ie, js, je, isq, ieq, jsq, jeq, i, j, m
484 isq = hi%IscB ; ieq = hi%IecB
485 jsq = hi%JscB ; jeq = hi%JecB
486 is = hi%isc ; ie = hi%iec
487 js = hi%jsc ; je = hi%jec
489 if (
present(pres_scale))
then
490 gxrho = pres_scale * g_e * rho_0 ; g_earth = pres_scale * g_e
491 pa_to_rl2_t2 = 1.0 / pres_scale
493 gxrho = g_e * rho_0 ; g_earth = g_e
496 if (
present(rho_scale))
then
497 g_earth = g_earth * rho_scale
498 rho_ref_mks = rho_ref / rho_scale ; i_rho = rho_scale / rho_0
500 rho_ref_mks = rho_ref ; i_rho = 1.0 / rho_0
503 do_massweight = .false.
504 if (
present(usemasswghtinterp))
then ;
if (usemasswghtinterp)
then
505 do_massweight = .true.
512 do j=jsq,jeq+1 ;
do i=isq,ieq+1
513 al0_2d(i,j) = (a0 + a1*t(i,j)) + a2*s(i,j)
514 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))
515 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))
517 al0 = al0_2d(i,j) ; p0 = p0_2d(i,j) ; lambda = lambda_2d(i,j)
519 dz = z_t(i,j) - z_b(i,j)
520 p_ave = -0.5*gxrho*(z_t(i,j)+z_b(i,j))
523 i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
524 eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
528 rho_anom = (p0 + p_ave)*(i_lzz*i_al0) - rho_ref_mks
529 rem = i_rho * (lambda * i_al0**2) * eps2 * &
530 (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
531 dpa(i,j) = pa_to_rl2_t2 * (g_earth*rho_anom*dz - 2.0*eps*rem)
532 if (
present(intz_dpa)) &
533 intz_dpa(i,j) = pa_to_rl2_t2 * (0.5*g_earth*rho_anom*dz**2 - dz*(1.0+eps)*rem)
536 if (
present(intx_dpa))
then ;
do j=js,je ;
do i=isq,ieq
542 hwght = max(0., -bathyt(i,j)-z_t(i+1,j), -bathyt(i+1,j)-z_t(i,j))
544 hl = (z_t(i,j) - z_b(i,j)) + dz_neglect
545 hr = (z_t(i+1,j) - z_b(i+1,j)) + dz_neglect
546 hwght = hwght * ( (hl-hr)/(hl+hr) )**2
547 idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
548 hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
549 hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
551 hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
554 intz(1) = dpa(i,j) ; intz(5) = dpa(i+1,j)
556 wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
557 wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
559 al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i+1,j)
560 p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i+1,j)
561 lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i+1,j)
563 dz = wt_l*(z_t(i,j) - z_b(i,j)) + wt_r*(z_t(i+1,j) - z_b(i+1,j))
564 p_ave = -0.5*gxrho*(wt_l*(z_t(i,j)+z_b(i,j)) + &
565 wt_r*(z_t(i+1,j)+z_b(i+1,j)))
568 i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
569 eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
571 intz(m) = pa_to_rl2_t2 * ( g_earth*dz*((p0 + p_ave)*(i_lzz*i_al0) - rho_ref_mks) - 2.0*eps * &
572 i_rho * (lambda * i_al0**2) * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2))) )
575 intx_dpa(i,j) = c1_90*(7.0*(intz(1)+intz(5)) + 32.0*(intz(2)+intz(4)) + 12.0*intz(3))
576 enddo ;
enddo ;
endif
578 if (
present(inty_dpa))
then ;
do j=jsq,jeq ;
do i=is,ie
584 hwght = max(0., -bathyt(i,j)-z_t(i,j+1), -bathyt(i,j+1)-z_t(i,j))
586 hl = (z_t(i,j) - z_b(i,j)) + dz_neglect
587 hr = (z_t(i,j+1) - z_b(i,j+1)) + dz_neglect
588 hwght = hwght * ( (hl-hr)/(hl+hr) )**2
589 idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
590 hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
591 hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
593 hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
596 intz(1) = dpa(i,j) ; intz(5) = dpa(i,j+1)
598 wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
599 wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
601 al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i,j+1)
602 p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i,j+1)
603 lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i,j+1)
605 dz = wt_l*(z_t(i,j) - z_b(i,j)) + wt_r*(z_t(i,j+1) - z_b(i,j+1))
606 p_ave = -0.5*gxrho*(wt_l*(z_t(i,j)+z_b(i,j)) + &
607 wt_r*(z_t(i,j+1)+z_b(i,j+1)))
610 i_lzz = 1.0 / (p0 + (lambda * i_al0) + p_ave)
611 eps = 0.5*gxrho*dz*i_lzz ; eps2 = eps*eps
613 intz(m) = pa_to_rl2_t2 * ( g_earth*dz*((p0 + p_ave)*(i_lzz*i_al0) - rho_ref_mks) - 2.0*eps * &
614 i_rho * (lambda * i_al0**2) * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2))) )
617 inty_dpa(i,j) = c1_90*(7.0*(intz(1)+intz(5)) + 32.0*(intz(2)+intz(4)) + 12.0*intz(3))
618 enddo ;
enddo ;
endif
620 end subroutine int_density_dz_wright
628 subroutine int_spec_vol_dp_wright(T, S, p_t, p_b, spv_ref, HI, dza, &
629 intp_dza, intx_dza, inty_dza, halo_size, &
630 bathyP, dP_neglect, useMassWghtInterp, SV_scale, pres_scale)
632 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
635 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
637 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
639 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
641 real,
intent(in) :: spv_ref
645 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
648 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
649 optional,
intent(inout) :: intp_dza
653 real,
dimension(HI%IsdB:HI%IedB,HI%jsd:HI%jed), &
654 optional,
intent(inout) :: intx_dza
658 real,
dimension(HI%isd:HI%ied,HI%JsdB:HI%JedB), &
659 optional,
intent(inout) :: inty_dza
663 integer,
optional,
intent(in) :: halo_size
665 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
666 optional,
intent(in) :: bathyp
667 real,
optional,
intent(in) :: dp_neglect
669 logical,
optional,
intent(in) :: usemasswghtinterp
671 real,
optional,
intent(in) :: sv_scale
673 real,
optional,
intent(in) :: pres_scale
677 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed) :: al0_2d, p0_2d, lambda_2d
692 real :: hwt_ll, hwt_lr
693 real :: hwt_rl, hwt_rr
698 logical :: do_massweight
699 real,
parameter :: c1_3 = 1.0/3.0, c1_7 = 1.0/7.0
700 real,
parameter :: c1_9 = 1.0/9.0, c1_90 = 1.0/90.0
701 integer :: isq, ieq, jsq, jeq, ish, ieh, jsh, jeh, i, j, m, halo
703 isq = hi%IscB ; ieq = hi%IecB ; jsq = hi%JscB ; jeq = hi%JecB
704 halo = 0 ;
if (
present(halo_size)) halo = max(halo_size,0)
705 ish = hi%isc-halo ; ieh = hi%iec+halo ; jsh = hi%jsc-halo ; jeh = hi%jec+halo
706 if (
present(intx_dza))
then ; ish = min(isq,ish) ; ieh = max(ieq+1,ieh);
endif
707 if (
present(inty_dza))
then ; jsh = min(jsq,jsh) ; jeh = max(jeq+1,jeh);
endif
710 al0_scale = 1.0 ;
if (
present(sv_scale)) al0_scale = sv_scale
712 if (
present(pres_scale))
then ;
if (pres_scale /= 1.0)
then
713 p0_scale = 1.0 / pres_scale
715 lam_scale = al0_scale * p0_scale
717 do_massweight = .false.
718 if (
present(usemasswghtinterp))
then ;
if (usemasswghtinterp)
then
719 do_massweight = .true.
727 do j=jsh,jeh ;
do i=ish,ieh
728 al0_2d(i,j) = al0_scale * ( (a0 + a1*t(i,j)) + a2*s(i,j) )
729 p0_2d(i,j) = p0_scale * ( (b0 + b4*s(i,j)) + t(i,j) * (b1 + t(i,j)*((b2 + b3*t(i,j))) + b5*s(i,j)) )
730 lambda_2d(i,j) = lam_scale * ( (c0 + c4*s(i,j)) + t(i,j) * (c1 + t(i,j)*((c2 + c3*t(i,j))) + c5*s(i,j)) )
732 al0 = al0_2d(i,j) ; p0 = p0_2d(i,j) ; lambda = lambda_2d(i,j)
733 dp = p_b(i,j) - p_t(i,j)
734 p_ave = 0.5*(p_t(i,j)+p_b(i,j))
736 eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
737 alpha_anom = al0 + lambda / (p0 + p_ave) - spv_ref
738 rem = lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
739 dza(i,j) = alpha_anom*dp + 2.0*eps*rem
740 if (
present(intp_dza)) &
741 intp_dza(i,j) = 0.5*alpha_anom*dp**2 - dp*(1.0-eps)*rem
744 if (
present(intx_dza))
then ;
do j=hi%jsc,hi%jec ;
do i=isq,ieq
750 hwght = max(0., bathyp(i,j)-p_t(i+1,j), bathyp(i+1,j)-p_t(i,j))
752 hl = (p_b(i,j) - p_t(i,j)) + dp_neglect
753 hr = (p_b(i+1,j) - p_t(i+1,j)) + dp_neglect
754 hwght = hwght * ( (hl-hr)/(hl+hr) )**2
755 idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
756 hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
757 hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
759 hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
762 intp(1) = dza(i,j) ; intp(5) = dza(i+1,j)
764 wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
765 wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
769 al0 = wtt_l*al0_2d(i,j) + wtt_r*al0_2d(i+1,j)
770 p0 = wtt_l*p0_2d(i,j) + wtt_r*p0_2d(i+1,j)
771 lambda = wtt_l*lambda_2d(i,j) + wtt_r*lambda_2d(i+1,j)
773 dp = wt_l*(p_b(i,j) - p_t(i,j)) + wt_r*(p_b(i+1,j) - p_t(i+1,j))
774 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)))
776 eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
777 intp(m) = (al0 + lambda / (p0 + p_ave) - spv_ref)*dp + 2.0*eps* &
778 lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
781 intx_dza(i,j) = c1_90*(7.0*(intp(1)+intp(5)) + 32.0*(intp(2)+intp(4)) + &
783 enddo ;
enddo ;
endif
785 if (
present(inty_dza))
then ;
do j=jsq,jeq ;
do i=hi%isc,hi%iec
791 hwght = max(0., bathyp(i,j)-p_t(i,j+1), bathyp(i,j+1)-p_t(i,j))
793 hl = (p_b(i,j) - p_t(i,j)) + dp_neglect
794 hr = (p_b(i,j+1) - p_t(i,j+1)) + dp_neglect
795 hwght = hwght * ( (hl-hr)/(hl+hr) )**2
796 idenom = 1.0 / ( hwght*(hr + hl) + hl*hr )
797 hwt_ll = (hwght*hl + hr*hl) * idenom ; hwt_lr = (hwght*hr) * idenom
798 hwt_rr = (hwght*hr + hr*hl) * idenom ; hwt_rl = (hwght*hl) * idenom
800 hwt_ll = 1.0 ; hwt_lr = 0.0 ; hwt_rr = 1.0 ; hwt_rl = 0.0
803 intp(1) = dza(i,j) ; intp(5) = dza(i,j+1)
805 wt_l = 0.25*real(5-m) ; wt_r = 1.0-wt_l
806 wtt_l = wt_l*hwt_ll + wt_r*hwt_rl ; wtt_r = wt_l*hwt_lr + wt_r*hwt_rr
810 al0 = wt_l*al0_2d(i,j) + wt_r*al0_2d(i,j+1)
811 p0 = wt_l*p0_2d(i,j) + wt_r*p0_2d(i,j+1)
812 lambda = wt_l*lambda_2d(i,j) + wt_r*lambda_2d(i,j+1)
814 dp = wt_l*(p_b(i,j) - p_t(i,j)) + wt_r*(p_b(i,j+1) - p_t(i,j+1))
815 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)))
817 eps = 0.5 * dp / (p0 + p_ave) ; eps2 = eps*eps
818 intp(m) = (al0 + lambda / (p0 + p_ave) - spv_ref)*dp + 2.0*eps* &
819 lambda * eps2 * (c1_3 + eps2*(0.2 + eps2*(c1_7 + c1_9*eps2)))
822 inty_dza(i,j) = c1_90*(7.0*(intp(1)+intp(5)) + 32.0*(intp(2)+intp(4)) + &
824 enddo ;
enddo ;
endif
825 end subroutine int_spec_vol_dp_wright