5#if defined TANGENT && defined FOUR_DVAR
61 & LBi, UBi, LBj, UBj, &
62 & IminS, ImaxS, JminS, JmaxS, &
63 & Nghost, NHsteps, DTsizeH, &
65 & pm, pn, pmon_u, pnom_v, &
67 & rmask, umask, vmask, &
82 integer,
intent(in) :: ng, tile, model
83 integer,
intent(in) :: LBi, UBi, LBj, UBj
84 integer,
intent(in) :: IminS, ImaxS, JminS, JmaxS
85 integer,
intent(in) :: Nghost, NHsteps
87 real(r8),
intent(in) :: DTsizeH
90 real(r8),
intent(in) :: pm(LBi:,LBj:)
91 real(r8),
intent(in) :: pn(LBi:,LBj:)
92 real(r8),
intent(in) :: pmon_u(LBi:,LBj:)
93 real(r8),
intent(in) :: pnom_v(LBi:,LBj:)
95 real(r8),
intent(in) :: rmask(LBi:,LBj:)
96 real(r8),
intent(in) :: umask(LBi:,LBj:)
97 real(r8),
intent(in) :: vmask(LBi:,LBj:)
99 real(r8),
intent(in) :: Kh(LBi:,LBj:)
100 real(r8),
intent(inout) :: tl_A(LBi:,LBj:)
102 real(r8),
intent(in) :: pm(LBi:UBi,LBj:UBj)
103 real(r8),
intent(in) :: pn(LBi:UBi,LBj:UBj)
104 real(r8),
intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
105 real(r8),
intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
107 real(r8),
intent(in) :: rmask(LBi:UBi,LBj:UBj)
108 real(r8),
intent(in) :: umask(LBi:UBi,LBj:UBj)
109 real(r8),
intent(in) :: vmask(LBi:UBi,LBj:UBj)
111 real(r8),
intent(in) :: Kh(LBi:UBi,LBj:UBj)
112 real(r8),
intent(inout) :: tl_A(LBi:UBi,LBj:UBj)
117 integer :: Nnew, Nold, Nsav, i, j, step
119 real(r8),
dimension(LBi:UBi,LBj:UBj,2) :: tl_Awrk
121 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
122 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
123 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac
125# include "set_bounds.h"
136 hfac(i,j)=dtsizeh*pm(i,j)*pn(i,j)
149 & lbi, ubi, lbj, ubj, &
159 & lbi, ubi, lbj, ubj, &
168 tl_awrk(i,j,nold)=tl_a(i,j)
185 tl_fx(i,j)=pmon_u(i,j)*0.5_r8*(kh(i-1,j)+kh(i,j))* &
186 & (tl_awrk(i,j,nold)-tl_awrk(i-1,j,nold))
190 tl_fx(i,j)=tl_fx(i,j)*umask(i,j)
199 tl_fe(i,j)=pnom_v(i,j)*0.5_r8*(kh(i,j-1)+kh(i,j))* &
200 & (tl_awrk(i,j,nold)-tl_awrk(i,j-1,nold))
204 tl_fe(i,j)=tl_fe(i,j)*vmask(i,j)
218 tl_awrk(i,j,nnew)=tl_awrk(i,j,nold)+ &
220 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
221 & tl_fe(i,j+1)-tl_fe(i,j))
232 & lbi, ubi, lbj, ubj, &
242 & lbi, ubi, lbj, ubj, &
263 tl_a(i,j)=tl_awrk(i,j,nold)
271 & lbi, ubi, lbj, ubj, &
281 & lbi, ubi, lbj, ubj, &
292 & LBi, UBi, LBj, UBj, &
293 & IminS, ImaxS, JminS, JmaxS, &
294 & Nghost, NHsteps, DTsizeH, &
296 & pm, pn, pmon_r, pnom_p, &
313 integer,
intent(in) :: ng, tile, model
314 integer,
intent(in) :: LBi, UBi, LBj, UBj
315 integer,
intent(in) :: IminS, ImaxS, JminS, JmaxS
316 integer,
intent(in) :: Nghost, NHsteps
318 real(r8),
intent(in) :: DTsizeH
321 real(r8),
intent(in) :: pm(LBi:,LBj:)
322 real(r8),
intent(in) :: pn(LBi:,LBj:)
323 real(r8),
intent(in) :: pmon_r(LBi:,LBj:)
324 real(r8),
intent(in) :: pnom_p(LBi:,LBj:)
326 real(r8),
intent(in) :: umask(LBi:,LBj:)
327 real(r8),
intent(in) :: pmask(LBi:,LBj:)
329 real(r8),
intent(in) :: Kh(LBi:,LBj:)
330 real(r8),
intent(inout) :: tl_A(LBi:,LBj:)
332 real(r8),
intent(in) :: pm(LBi:UBi,LBj:UBj)
333 real(r8),
intent(in) :: pn(LBi:UBi,LBj:UBj)
334 real(r8),
intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
335 real(r8),
intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
337 real(r8),
intent(in) :: umask(LBi:UBi,LBj:UBj)
338 real(r8),
intent(in) :: pmask(LBi:UBi,LBj:UBj)
340 real(r8),
intent(in) :: Kh(LBi:UBi,LBj:UBj)
341 real(r8),
intent(inout) :: tl_A(LBi:UBi,LBj:UBj)
346 integer :: Nnew, Nold, Nsav, i, j, step
350 real(r8),
dimension(LBi:UBi,LBj:UBj,2) :: tl_Awrk
352 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
353 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
354 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac
356# include "set_bounds.h"
368 hfac(i,j)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
381 & lbi, ubi, lbj, ubj, &
391 & lbi, ubi, lbj, ubj, &
400 tl_awrk(i,j,nold)=tl_a(i,j)
417 tl_fx(i,j)=pmon_r(i,j)*kh(i,j)* &
418 & (tl_awrk(i+1,j,nold)-tl_awrk(i,j,nold))
427 tl_fe(i,j)=pnom_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
428 & kh(i-1,j-1)+kh(i,j-1))* &
429 & (tl_awrk(i,j,nold)-tl_awrk(i,j-1,nold))
433 tl_fe(i,j)=tl_fe(i,j)*pmask(i,j)
447 tl_awrk(i,j,nnew)=tl_awrk(i,j,nold)+ &
449 & (tl_fx(i,j)-tl_fx(i-1,j)+ &
450 & tl_fe(i,j+1)-tl_fe(i,j))
461 & lbi, ubi, lbj, ubj, &
471 & lbi, ubi, lbj, ubj, &
492 tl_a(i,j)=tl_awrk(i,j,nold)
500 & lbi, ubi, lbj, ubj, &
510 & lbi, ubi, lbj, ubj, &
521 & LBi, UBi, LBj, UBj, &
522 & IminS, ImaxS, JminS, JmaxS, &
523 & Nghost, NHsteps, DTsizeH, &
525 & pm, pn, pmon_p, pnom_r, &
542 integer,
intent(in) :: ng, tile, model
543 integer,
intent(in) :: LBi, UBi, LBj, UBj
544 integer,
intent(in) :: IminS, ImaxS, JminS, JmaxS
545 integer,
intent(in) :: Nghost, NHsteps
547 real(r8),
intent(in) :: DTsizeH
550 real(r8),
intent(in) :: pm(LBi:,LBj:)
551 real(r8),
intent(in) :: pn(LBi:,LBj:)
552 real(r8),
intent(in) :: pmon_p(LBi:,LBj:)
553 real(r8),
intent(in) :: pnom_r(LBi:,LBj:)
555 real(r8),
intent(in) :: vmask(LBi:,LBj:)
556 real(r8),
intent(in) :: pmask(LBi:,LBj:)
558 real(r8),
intent(in) :: Kh(LBi:,LBj:)
559 real(r8),
intent(inout) :: tl_A(LBi:,LBj:)
561 real(r8),
intent(in) :: pm(LBi:UBi,LBj:UBj)
562 real(r8),
intent(in) :: pn(LBi:UBi,LBj:UBj)
563 real(r8),
intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
564 real(r8),
intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
566 real(r8),
intent(in) :: vmask(LBi:UBi,LBj:UBj)
567 real(r8),
intent(in) :: pmask(LBi:UBi,LBj:UBj)
569 real(r8),
intent(in) :: Kh(LBi:UBi,LBj:UBj)
570 real(r8),
intent(inout) :: tl_A(LBi:UBi,LBj:UBj)
575 integer :: Nnew, Nold, Nsav, i, j, step
579 real(r8),
dimension(LBi:UBi,LBj:UBj,2) :: tl_Awrk
581 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
582 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
583 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac
585# include "set_bounds.h"
597 hfac(i,j)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
610 & lbi, ubi, lbj, ubj, &
620 & lbi, ubi, lbj, ubj, &
629 tl_awrk(i,j,nold)=tl_a(i,j)
647 tl_fx(i,j)=pmon_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
648 & kh(i-1,j-1)+kh(i,j-1))* &
649 & (tl_awrk(i,j,nold)-tl_awrk(i-1,j,nold))
654 tl_fx(i,j)=tl_fx(i,j)*pmask(i,j)
663 tl_fe(i,j)=pnom_r(i,j)*kh(i,j)* &
664 & (tl_awrk(i,j+1,nold)-tl_awrk(i,j,nold))
677 tl_awrk(i,j,nnew)=tl_awrk(i,j,nold)+ &
679 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
680 & tl_fe(i,j)-tl_fe(i,j-1))
691 & lbi, ubi, lbj, ubj, &
701 & lbi, ubi, lbj, ubj, &
722 tl_a(i,j)=tl_awrk(i,j,nold)
730 & lbi, ubi, lbj, ubj, &
740 & lbi, ubi, lbj, ubj, &
subroutine dabc_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine dabc_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine dabc_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine tl_conv_u2d_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_r, pnom_p, umask, pmask, tl_a)
subroutine tl_conv_r2d_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_u, pnom_v, rmask, umask, vmask, tl_a)
subroutine tl_conv_v2d_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_p, pnom_r, vmask, pmask, tl_a)