12 implicit none ;
private
14 #include <MOM_memory.h>
27 slope_x, slope_y, N2_u, N2_v, halo)
31 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
32 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: e
36 real,
intent(in) :: dt_kappa_smooth
38 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: slope_x
39 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)+1),
intent(inout) :: slope_y
40 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), &
41 optional,
intent(inout) :: n2_u
43 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)+1), &
44 optional,
intent(inout) :: n2_v
46 integer,
optional,
intent(in) :: halo
51 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)) :: &
57 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
59 real,
dimension(SZIB_(G)) :: &
62 real,
dimension(SZI_(G)) :: &
65 real,
dimension(SZIB_(G)) :: &
69 real,
dimension(SZI_(G)) :: &
79 real :: haa, hab, hal, har
81 real :: wta, wtb, wtl, wtr
101 logical :: present_n2_u, present_n2_v
102 integer :: is, ie, js, je, nz, isdb
105 if (
present(halo))
then
106 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
108 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
110 nz = g%ke ; isdb = g%IsdB
112 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect**2
113 z_to_l = us%Z_to_L ; h_to_z = gv%H_to_Z
117 l_to_z = 1.0 / z_to_l
118 dz_neglect = gv%H_subroundoff * h_to_z
120 use_eos =
associated(tv%eqn_of_state)
122 present_n2_u =
PRESENT(n2_u)
123 present_n2_v =
PRESENT(n2_v)
124 g_rho0 = (us%L_to_Z*l_to_z*gv%g_Earth) / gv%Rho0
125 if (present_n2_u)
then
126 do j=js,je ;
do i=is-1,ie
131 if (present_n2_v)
then
132 do j=js-1,je ;
do i=is,ie
139 if (
present(halo))
then
140 call vert_fill_ts(h, tv%T, tv%S, dt_kappa_smooth, t, s, g, gv, halo+1)
142 call vert_fill_ts(h, tv%T, tv%S, dt_kappa_smooth, t, s, g, gv, 1)
148 do j=js-1,je+1 ;
do i=is-1,ie+1
150 pres(i,j,2) = pres(i,j,1) + gv%H_to_Pa*h(i,j,1)
154 do k=2,nz ;
do i=is-1,ie+1
155 pres(i,j,k+1) = pres(i,j,k) + gv%H_to_Pa*h(i,j,k)
166 do j=js,je ;
do k=nz,2,-1
167 if (.not.(use_eos))
then
168 drdia = 0.0 ; drdib = 0.0
169 drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = gv%Rlay(k)-gv%Rlay(k-1)
175 pres_u(i) = 0.5*(pres(i,j,k) + pres(i+1,j,k))
176 t_u(i) = 0.25*((t(i,j,k) + t(i+1,j,k)) + (t(i,j,k-1) + t(i+1,j,k-1)))
177 s_u(i) = 0.25*((s(i,j,k) + s(i+1,j,k)) + (s(i,j,k-1) + s(i+1,j,k-1)))
180 drho_ds_u, (is-isdb+1)-1, ie-is+2, tv%eqn_of_state, scale=us%kg_m3_to_R)
186 drdia = drho_dt_u(i) * (t(i+1,j,k-1)-t(i,j,k-1)) + &
187 drho_ds_u(i) * (s(i+1,j,k-1)-s(i,j,k-1))
188 drdib = drho_dt_u(i) * (t(i+1,j,k)-t(i,j,k)) + &
189 drho_ds_u(i) * (s(i+1,j,k)-s(i,j,k))
192 drdkl = (drho_dt_u(i) * (t(i,j,k)-t(i,j,k-1)) + &
193 drho_ds_u(i) * (s(i,j,k)-s(i,j,k-1)))
194 drdkr = (drho_dt_u(i) * (t(i+1,j,k)-t(i+1,j,k-1)) + &
195 drho_ds_u(i) * (s(i+1,j,k)-s(i+1,j,k-1)))
200 hg2a = h(i,j,k-1)*h(i+1,j,k-1) + h_neglect2
201 hg2b = h(i,j,k)*h(i+1,j,k) + h_neglect2
202 hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
203 hg2r = h(i+1,j,k-1)*h(i+1,j,k) + h_neglect2
204 haa = 0.5*(h(i,j,k-1) + h(i+1,j,k-1))
205 hab = 0.5*(h(i,j,k) + h(i+1,j,k)) + h_neglect
206 hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
207 har = 0.5*(h(i+1,j,k-1) + h(i+1,j,k)) + h_neglect
208 if (gv%Boussinesq)
then
209 dzal = hal * h_to_z ; dzar = har * h_to_z
211 dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
212 dzar = 0.5*(e(i+1,j,k-1) - e(i+1,j,k+1)) + dz_neglect
216 wta = hg2a*hab ; wtb = hg2b*haa
217 wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
219 drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
224 drdx = ((wta * drdia + wtb * drdib) / (wta + wtb) - &
225 drdz * (e(i,j,k)-e(i+1,j,k))) * g%IdxCu(i,j)
229 mag_grad2 = drdx**2 + (l_to_z*drdz)**2
230 if (mag_grad2 > 0.0)
then
231 slope_x(i,j,k) = drdx / sqrt(mag_grad2)
236 if (present_n2_u) n2_u(i,j,k) = g_rho0 * drdz * g%mask2dCu(i,j)
239 slope_x(i,j,k) = (z_to_l*(e(i,j,k)-e(i+1,j,k))) * g%IdxCu(i,j)
253 do j=js-1,je ;
do k=nz,2,-1
254 if (.not.(use_eos))
then
255 drdja = 0.0 ; drdjb = 0.0
256 drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = gv%Rlay(k)-gv%Rlay(k-1)
261 pres_v(i) = 0.5*(pres(i,j,k) + pres(i,j+1,k))
262 t_v(i) = 0.25*((t(i,j,k) + t(i,j+1,k)) + (t(i,j,k-1) + t(i,j+1,k-1)))
263 s_v(i) = 0.25*((s(i,j,k) + s(i,j+1,k)) + (s(i,j,k-1) + s(i,j+1,k-1)))
266 drho_ds_v, is, ie-is+1, tv%eqn_of_state, scale=us%kg_m3_to_R)
271 drdja = drho_dt_v(i) * (t(i,j+1,k-1)-t(i,j,k-1)) + &
272 drho_ds_v(i) * (s(i,j+1,k-1)-s(i,j,k-1))
273 drdjb = drho_dt_v(i) * (t(i,j+1,k)-t(i,j,k)) + &
274 drho_ds_v(i) * (s(i,j+1,k)-s(i,j,k))
277 drdkl = (drho_dt_v(i) * (t(i,j,k)-t(i,j,k-1)) + &
278 drho_ds_v(i) * (s(i,j,k)-s(i,j,k-1)))
279 drdkr = (drho_dt_v(i) * (t(i,j+1,k)-t(i,j+1,k-1)) + &
280 drho_ds_v(i) * (s(i,j+1,k)-s(i,j+1,k-1)))
284 hg2a = h(i,j,k-1)*h(i,j+1,k-1) + h_neglect2
285 hg2b = h(i,j,k)*h(i,j+1,k) + h_neglect2
286 hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
287 hg2r = h(i,j+1,k-1)*h(i,j+1,k) + h_neglect2
288 haa = 0.5*(h(i,j,k-1) + h(i,j+1,k-1)) + h_neglect
289 hab = 0.5*(h(i,j,k) + h(i,j+1,k)) + h_neglect
290 hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
291 har = 0.5*(h(i,j+1,k-1) + h(i,j+1,k)) + h_neglect
292 if (gv%Boussinesq)
then
293 dzal = hal * h_to_z ; dzar = har * h_to_z
295 dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
296 dzar = 0.5*(e(i,j+1,k-1) - e(i,j+1,k+1)) + dz_neglect
300 wta = hg2a*hab ; wtb = hg2b*haa
301 wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
303 drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
308 drdy = ((wta * drdja + wtb * drdjb) / (wta + wtb) - &
309 drdz * (e(i,j,k)-e(i,j+1,k))) * g%IdyCv(i,j)
313 mag_grad2 = drdy**2 + (l_to_z*drdz)**2
314 if (mag_grad2 > 0.0)
then
315 slope_y(i,j,k) = drdy / sqrt(mag_grad2)
320 if (present_n2_v) n2_v(i,j,k) = g_rho0 * drdz * g%mask2dCv(i,j)
323 slope_y(i,j,k) = (z_to_l*(e(i,j,k)-e(i,j+1,k))) * g%IdyCv(i,j)
333 subroutine vert_fill_ts(h, T_in, S_in, kappa_dt, T_f, S_f, G, GV, halo_here, larger_h_denom)
336 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
337 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: t_in
338 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: s_in
339 real,
intent(in) :: kappa_dt
341 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: t_f
342 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: s_f
343 integer,
optional,
intent(in) :: halo_here
345 logical,
optional,
intent(in) :: larger_h_denom
351 real :: ent(szi_(g),szk_(g)+1)
353 real :: b1(szi_(g)), d1(szi_(g))
354 real :: c1(szi_(g),szk_(g))
362 integer :: i, j, k, is, ie, js, je, nz, halo
364 halo=0 ;
if (
present(halo_here)) halo = max(halo_here,0)
366 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo ; nz = gv%ke
368 h_neglect = gv%H_subroundoff
369 kap_dt_x2 = (2.0*kappa_dt)*gv%Z_to_H**2
371 if (
present(larger_h_denom))
then
372 if (larger_h_denom) h0 = 1.0e-16*sqrt(kappa_dt)*gv%Z_to_H
375 if (kap_dt_x2 <= 0.0)
then
377 do k=1,nz ;
do j=js,je ;
do i=is,ie
378 t_f(i,j,k) = t_in(i,j,k) ; s_f(i,j,k) = s_in(i,j,k)
379 enddo ;
enddo ;
enddo
384 ent(i,2) = kap_dt_x2 / ((h(i,j,1)+h(i,j,2)) + h0)
385 h_tr = h(i,j,1) + h_neglect
386 b1(i) = 1.0 / (h_tr + ent(i,2))
388 t_f(i,j,1) = (b1(i)*h_tr)*t_in(i,j,1)
389 s_f(i,j,1) = (b1(i)*h_tr)*s_in(i,j,1)
391 do k=2,nz-1 ;
do i=is,ie
392 ent(i,k+1) = kap_dt_x2 / ((h(i,j,k)+h(i,j,k+1)) + h0)
393 h_tr = h(i,j,k) + h_neglect
394 c1(i,k) = ent(i,k) * b1(i)
395 b1(i) = 1.0 / ((h_tr + d1(i)*ent(i,k)) + ent(i,k+1))
396 d1(i) = b1(i) * (h_tr + d1(i)*ent(i,k))
397 t_f(i,j,k) = b1(i) * (h_tr*t_in(i,j,k) + ent(i,k)*t_f(i,j,k-1))
398 s_f(i,j,k) = b1(i) * (h_tr*s_in(i,j,k) + ent(i,k)*s_f(i,j,k-1))
401 c1(i,nz) = ent(i,nz) * b1(i)
402 h_tr = h(i,j,nz) + h_neglect
403 b1(i) = 1.0 / (h_tr + d1(i)*ent(i,nz))
404 t_f(i,j,nz) = b1(i) * (h_tr*t_in(i,j,nz) + ent(i,nz)*t_f(i,j,nz-1))
405 s_f(i,j,nz) = b1(i) * (h_tr*s_in(i,j,nz) + ent(i,nz)*s_f(i,j,nz-1))
407 do k=nz-1,1,-1 ;
do i=is,ie
408 t_f(i,j,k) = t_f(i,j,k) + c1(i,k+1)*t_f(i,j,k+1)
409 s_f(i,j,k) = s_f(i,j,k) + c1(i,k+1)*s_f(i,j,k+1)