12 use mom_open_boundary,
only : obc_direction_e, obc_direction_w, obc_direction_n, obc_direction_s
14 implicit none ;
private
16 #include <MOM_memory.h>
18 public calc_isoneutral_slopes, vert_fill_ts
28 subroutine calc_isoneutral_slopes(G, GV, US, h, e, tv, dt_kappa_smooth, &
29 slope_x, slope_y, N2_u, N2_v, halo, OBC)
33 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
34 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: e
38 real,
intent(in) :: dt_kappa_smooth
40 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: slope_x
41 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)+1),
intent(inout) :: slope_y
42 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), &
43 optional,
intent(inout) :: n2_u
45 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)+1), &
46 optional,
intent(inout) :: n2_v
48 integer,
optional,
intent(in) :: halo
54 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)) :: &
60 real,
dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
62 real,
dimension(SZIB_(G)) :: &
65 real,
dimension(SZI_(G)) :: &
68 real,
dimension(SZIB_(G)) :: &
72 real,
dimension(SZI_(G)) :: &
82 real :: haa, hab, hal, har
84 real :: wta, wtb, wtl, wtr
104 logical :: present_n2_u, present_n2_v
105 integer,
dimension(2) :: eosdom_u, eosdom_v
106 integer :: is, ie, js, je, nz, isdb
109 logical :: local_open_u_bc, local_open_v_bc
111 if (
present(halo))
then
112 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
114 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
116 nz = g%ke ; isdb = g%IsdB
118 h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect**2
119 z_to_l = us%Z_to_L ; h_to_z = gv%H_to_Z
123 l_to_z = 1.0 / z_to_l
124 dz_neglect = gv%H_subroundoff * h_to_z
126 local_open_u_bc = .false.
127 local_open_v_bc = .false.
128 if (
present(obc))
then ;
if (
associated(obc))
then
129 local_open_u_bc = obc%open_u_BCs_exist_globally
130 local_open_v_bc = obc%open_v_BCs_exist_globally
133 use_eos =
associated(tv%eqn_of_state)
135 present_n2_u =
PRESENT(n2_u)
136 present_n2_v =
PRESENT(n2_v)
137 g_rho0 = (us%L_to_Z*l_to_z*gv%g_Earth) / gv%Rho0
138 if (present_n2_u)
then
139 do j=js,je ;
do i=is-1,ie
144 if (present_n2_v)
then
145 do j=js-1,je ;
do i=is,ie
152 if (
present(halo))
then
153 call vert_fill_ts(h, tv%T, tv%S, dt_kappa_smooth, t, s, g, gv, halo+1)
155 call vert_fill_ts(h, tv%T, tv%S, dt_kappa_smooth, t, s, g, gv, 1)
160 if (
associated(tv%p_surf))
then
162 do j=js-1,je+1 ;
do i=is-1,ie+1
163 pres(i,j,1) = tv%p_surf(i,j)
167 do j=js-1,je+1 ;
do i=is-1,ie+1
173 do k=1,nz ;
do i=is-1,ie+1
174 pres(i,j,k+1) = pres(i,j,k) + gv%g_Earth * gv%H_to_RZ * h(i,j,k)
178 eosdom_u(1) = is-1 - (g%IsdB-1) ; eosdom_u(2) = ie - (g%IsdB-1)
188 do j=js,je ;
do k=nz,2,-1
189 if (.not.(use_eos))
then
190 drdia = 0.0 ; drdib = 0.0
191 drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = gv%Rlay(k)-gv%Rlay(k-1)
197 pres_u(i) = 0.5*(pres(i,j,k) + pres(i+1,j,k))
198 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)))
199 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)))
202 tv%eqn_of_state, eosdom_u)
208 drdia = drho_dt_u(i) * (t(i+1,j,k-1)-t(i,j,k-1)) + &
209 drho_ds_u(i) * (s(i+1,j,k-1)-s(i,j,k-1))
210 drdib = drho_dt_u(i) * (t(i+1,j,k)-t(i,j,k)) + &
211 drho_ds_u(i) * (s(i+1,j,k)-s(i,j,k))
214 drdkl = (drho_dt_u(i) * (t(i,j,k)-t(i,j,k-1)) + &
215 drho_ds_u(i) * (s(i,j,k)-s(i,j,k-1)))
216 drdkr = (drho_dt_u(i) * (t(i+1,j,k)-t(i+1,j,k-1)) + &
217 drho_ds_u(i) * (s(i+1,j,k)-s(i+1,j,k-1)))
222 hg2a = h(i,j,k-1)*h(i+1,j,k-1) + h_neglect2
223 hg2b = h(i,j,k)*h(i+1,j,k) + h_neglect2
224 hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
225 hg2r = h(i+1,j,k-1)*h(i+1,j,k) + h_neglect2
226 haa = 0.5*(h(i,j,k-1) + h(i+1,j,k-1))
227 hab = 0.5*(h(i,j,k) + h(i+1,j,k)) + h_neglect
228 hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
229 har = 0.5*(h(i+1,j,k-1) + h(i+1,j,k)) + h_neglect
230 if (gv%Boussinesq)
then
231 dzal = hal * h_to_z ; dzar = har * h_to_z
233 dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
234 dzar = 0.5*(e(i+1,j,k-1) - e(i+1,j,k+1)) + dz_neglect
238 wta = hg2a*hab ; wtb = hg2b*haa
239 wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
241 drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
246 drdx = ((wta * drdia + wtb * drdib) / (wta + wtb) - &
247 drdz * (e(i,j,k)-e(i+1,j,k))) * g%IdxCu(i,j)
251 mag_grad2 = drdx**2 + (l_to_z*drdz)**2
252 if (mag_grad2 > 0.0)
then
253 slope_x(i,j,k) = drdx / sqrt(mag_grad2)
258 if (present_n2_u) n2_u(i,j,k) = g_rho0 * drdz * g%mask2dCu(i,j)
261 slope_x(i,j,k) = (z_to_l*(e(i,j,k)-e(i+1,j,k))) * g%IdxCu(i,j)
263 if (local_open_u_bc)
then
264 l_seg = obc%segnum_u(i,j)
265 if (l_seg /= obc_none)
then
266 if (obc%segment(l_seg)%open)
then
277 slope_x(i,j,k) = slope_x(i,j,k) * max(g%mask2dT(i,j),g%mask2dT(i+1,j))
283 eosdom_v(1) = is - (g%isd-1) ; eosdom_v(2) = ie - (g%isd-1)
294 do j=js-1,je ;
do k=nz,2,-1
295 if (.not.(use_eos))
then
296 drdja = 0.0 ; drdjb = 0.0
297 drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = gv%Rlay(k)-gv%Rlay(k-1)
302 pres_v(i) = 0.5*(pres(i,j,k) + pres(i,j+1,k))
303 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)))
304 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)))
312 drdja = drho_dt_v(i) * (t(i,j+1,k-1)-t(i,j,k-1)) + &
313 drho_ds_v(i) * (s(i,j+1,k-1)-s(i,j,k-1))
314 drdjb = drho_dt_v(i) * (t(i,j+1,k)-t(i,j,k)) + &
315 drho_ds_v(i) * (s(i,j+1,k)-s(i,j,k))
318 drdkl = (drho_dt_v(i) * (t(i,j,k)-t(i,j,k-1)) + &
319 drho_ds_v(i) * (s(i,j,k)-s(i,j,k-1)))
320 drdkr = (drho_dt_v(i) * (t(i,j+1,k)-t(i,j+1,k-1)) + &
321 drho_ds_v(i) * (s(i,j+1,k)-s(i,j+1,k-1)))
325 hg2a = h(i,j,k-1)*h(i,j+1,k-1) + h_neglect2
326 hg2b = h(i,j,k)*h(i,j+1,k) + h_neglect2
327 hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
328 hg2r = h(i,j+1,k-1)*h(i,j+1,k) + h_neglect2
329 haa = 0.5*(h(i,j,k-1) + h(i,j+1,k-1)) + h_neglect
330 hab = 0.5*(h(i,j,k) + h(i,j+1,k)) + h_neglect
331 hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
332 har = 0.5*(h(i,j+1,k-1) + h(i,j+1,k)) + h_neglect
333 if (gv%Boussinesq)
then
334 dzal = hal * h_to_z ; dzar = har * h_to_z
336 dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
337 dzar = 0.5*(e(i,j+1,k-1) - e(i,j+1,k+1)) + dz_neglect
341 wta = hg2a*hab ; wtb = hg2b*haa
342 wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
344 drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
349 drdy = ((wta * drdja + wtb * drdjb) / (wta + wtb) - &
350 drdz * (e(i,j,k)-e(i,j+1,k))) * g%IdyCv(i,j)
354 mag_grad2 = drdy**2 + (l_to_z*drdz)**2
355 if (mag_grad2 > 0.0)
then
356 slope_y(i,j,k) = drdy / sqrt(mag_grad2)
361 if (present_n2_v) n2_v(i,j,k) = g_rho0 * drdz * g%mask2dCv(i,j)
364 slope_y(i,j,k) = (z_to_l*(e(i,j,k)-e(i,j+1,k))) * g%IdyCv(i,j)
366 if (local_open_v_bc)
then
367 l_seg = obc%segnum_v(i,j)
368 if (l_seg /= obc_none)
then
369 if (obc%segment(l_seg)%open)
then
380 slope_y(i,j,k) = slope_y(i,j,k) * max(g%mask2dT(i,j),g%mask2dT(i,j+1))
386 end subroutine calc_isoneutral_slopes
390 subroutine vert_fill_ts(h, T_in, S_in, kappa_dt, T_f, S_f, G, GV, halo_here, larger_h_denom)
393 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
394 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: t_in
395 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: s_in
396 real,
intent(in) :: kappa_dt
398 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: t_f
399 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: s_f
400 integer,
optional,
intent(in) :: halo_here
402 logical,
optional,
intent(in) :: larger_h_denom
408 real :: ent(szi_(g),szk_(g)+1)
410 real :: b1(szi_(g)), d1(szi_(g))
411 real :: c1(szi_(g),szk_(g))
419 integer :: i, j, k, is, ie, js, je, nz, halo
421 halo=0 ;
if (
present(halo_here)) halo = max(halo_here,0)
423 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo ; nz = gv%ke
425 h_neglect = gv%H_subroundoff
426 kap_dt_x2 = (2.0*kappa_dt)*gv%Z_to_H**2
428 if (
present(larger_h_denom))
then
429 if (larger_h_denom) h0 = 1.0e-16*sqrt(kappa_dt)*gv%Z_to_H
432 if (kap_dt_x2 <= 0.0)
then
434 do k=1,nz ;
do j=js,je ;
do i=is,ie
435 t_f(i,j,k) = t_in(i,j,k) ; s_f(i,j,k) = s_in(i,j,k)
436 enddo ;
enddo ;
enddo
441 ent(i,2) = kap_dt_x2 / ((h(i,j,1)+h(i,j,2)) + h0)
442 h_tr = h(i,j,1) + h_neglect
443 b1(i) = 1.0 / (h_tr + ent(i,2))
445 t_f(i,j,1) = (b1(i)*h_tr)*t_in(i,j,1)
446 s_f(i,j,1) = (b1(i)*h_tr)*s_in(i,j,1)
448 do k=2,nz-1 ;
do i=is,ie
449 ent(i,k+1) = kap_dt_x2 / ((h(i,j,k)+h(i,j,k+1)) + h0)
450 h_tr = h(i,j,k) + h_neglect
451 c1(i,k) = ent(i,k) * b1(i)
452 b1(i) = 1.0 / ((h_tr + d1(i)*ent(i,k)) + ent(i,k+1))
453 d1(i) = b1(i) * (h_tr + d1(i)*ent(i,k))
454 t_f(i,j,k) = b1(i) * (h_tr*t_in(i,j,k) + ent(i,k)*t_f(i,j,k-1))
455 s_f(i,j,k) = b1(i) * (h_tr*s_in(i,j,k) + ent(i,k)*s_f(i,j,k-1))
458 c1(i,nz) = ent(i,nz) * b1(i)
459 h_tr = h(i,j,nz) + h_neglect
460 b1(i) = 1.0 / (h_tr + d1(i)*ent(i,nz))
461 t_f(i,j,nz) = b1(i) * (h_tr*t_in(i,j,nz) + ent(i,nz)*t_f(i,j,nz-1))
462 s_f(i,j,nz) = b1(i) * (h_tr*s_in(i,j,nz) + ent(i,nz)*s_f(i,j,nz-1))
464 do k=nz-1,1,-1 ;
do i=is,ie
465 t_f(i,j,k) = t_f(i,j,k) + c1(i,k+1)*t_f(i,j,k+1)
466 s_f(i,j,k) = s_f(i,j,k) + c1(i,k+1)*s_f(i,j,k+1)
471 end subroutine vert_fill_ts