New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
traldf_iso.F90 in NEMO/branches/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE/TRA – NEMO

source: NEMO/branches/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE/TRA/traldf_iso.F90 @ 14632

Last change on this file since 14632 was 14632, checked in by hadcv, 3 years ago

#2600: Take out some unnecessary tags

  • Property svn:keywords set to Id
File size: 23.9 KB
Line 
1MODULE traldf_iso
2   !!======================================================================
3   !!                   ***  MODULE  traldf_iso  ***
4   !! Ocean  tracers:  horizontal component of the lateral tracer mixing trend
5   !!======================================================================
6   !! History :  OPA  ! 1994-08  (G. Madec, M. Imbard)
7   !!            8.0  ! 1997-05  (G. Madec)  split into traldf and trazdf
8   !!            NEMO ! 2002-08  (G. Madec)  Free form, F90
9   !!            1.0  ! 2005-11  (G. Madec)  merge traldf and trazdf :-)
10   !!            3.3  ! 2010-09  (C. Ethe, G. Madec) Merge TRA-TRC
11   !!            3.7  ! 2014-01  (G. Madec, S. Masson)  restructuration/simplification of aht/aeiv specification
12   !!             -   ! 2014-02  (F. Lemarie, G. Madec)  triad operator (Griffies) + Method of Stabilizing Correction
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   tra_ldf_iso   : update the tracer trend with the horizontal component of a iso-neutral laplacian operator
17   !!                   and with the vertical part of the isopycnal or geopotential s-coord. operator
18   !!----------------------------------------------------------------------
19   USE oce            ! ocean dynamics and active tracers
20   USE dom_oce        ! ocean space and time domain
21   USE domutl, ONLY : is_tile
22   USE trc_oce        ! share passive tracers/Ocean variables
23   USE zdf_oce        ! ocean vertical physics
24   USE ldftra         ! lateral diffusion: tracer eddy coefficients
25   USE ldfslp         ! iso-neutral slopes
26   USE diaptr         ! poleward transport diagnostics
27   USE diaar5         ! AR5 diagnostics
28   !
29   USE in_out_manager ! I/O manager
30   USE iom            ! I/O library
31   USE phycst         ! physical constants
32   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
33
34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC   tra_ldf_iso   ! routine called by step.F90
38
39   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
40   LOGICAL  ::   l_hst   ! flag to compute heat transport
41
42   !! * Substitutions
43#  include "do_loop_substitute.h90"
44#  include "domzgr_substitute.h90"
45   !!----------------------------------------------------------------------
46   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
47   !! $Id$
48   !! Software governed by the CeCILL license (see ./LICENSE)
49   !!----------------------------------------------------------------------
50CONTAINS
51
52   SUBROUTINE tra_ldf_iso( kt, Kmm, kit000, cdtype, pahu, pahv,             &
53      &                                             pgu , pgv , pgui, pgvi, &
54      &                                             pt, pt2, pt_rhs, kjpt, kpass )
55      !!
56      INTEGER                     , INTENT(in   ) ::   kt         ! ocean time-step index
57      INTEGER                     , INTENT(in   ) ::   kit000     ! first time step index
58      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
59      INTEGER                     , INTENT(in   ) ::   kjpt       ! number of tracers
60      INTEGER                     , INTENT(in   ) ::   kpass      ! =1/2 first or second passage
61      INTEGER                     , INTENT(in   ) ::   Kmm        ! ocean time level index
62      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s]
63      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels
64      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels
65      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2)
66      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2)
67      REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) ::   pt_rhs     ! tracer trend
68      !!
69      CALL tra_ldf_iso_t( kt, Kmm, kit000, cdtype, pahu, pahv, is_tile(pahu),                             &
70         &                                         pgu , pgv , is_tile(pgu) , pgui, pgvi, is_tile(pgui),  &
71         &                                         pt, is_tile(pt), pt2, is_tile(pt2), pt_rhs, is_tile(pt_rhs), kjpt, kpass )
72   END SUBROUTINE tra_ldf_iso
73
74
75  SUBROUTINE tra_ldf_iso_t( kt, Kmm, kit000, cdtype, pahu, pahv, ktah,                    &
76      &                                              pgu , pgv , ktg , pgui, pgvi, ktgi,  &
77      &                                              pt, ktt, pt2, ktt2, pt_rhs, ktt_rhs, kjpt, kpass )
78      !!----------------------------------------------------------------------
79      !!                  ***  ROUTINE tra_ldf_iso  ***
80      !!
81      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive
82      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and
83      !!      add it to the general trend of tracer equation.
84      !!
85      !! ** Method  :   The horizontal component of the lateral diffusive trends
86      !!      is provided by a 2nd order operator rotated along neural or geopo-
87      !!      tential surfaces to which an eddy induced advection can be added
88      !!      It is computed using before fields (forward in time) and isopyc-
89      !!      nal or geopotential slopes computed in routine ldfslp.
90      !!
91      !!      1st part :  masked horizontal derivative of T  ( di[ t ] )
92      !!      ========    with partial cell update if ln_zps=T
93      !!                  with top     cell update if ln_isfcav
94      !!
95      !!      2nd part :  horizontal fluxes of the lateral mixing operator
96      !!      ========
97      !!         zftu =  pahu e2u*e3u/e1u di[ tb ]
98      !!               - pahu e2u*uslp    dk[ mi(mk(tb)) ]
99      !!         zftv =  pahv e1v*e3v/e2v dj[ tb ]
100      !!               - pahv e2u*vslp    dk[ mj(mk(tb)) ]
101      !!      take the horizontal divergence of the fluxes:
102      !!         difft = 1/(e1e2t*e3t) {  di-1[ zftu ] +  dj-1[ zftv ]  }
103      !!      Add this trend to the general trend (ta,sa):
104      !!         ta = ta + difft
105      !!
106      !!      3rd part: vertical trends of the lateral mixing operator
107      !!      ========  (excluding the vertical flux proportional to dk[t] )
108      !!      vertical fluxes associated with the rotated lateral mixing:
109      !!         zftw = - {  mi(mk(pahu)) * e2t*wslpi di[ mi(mk(tb)) ]
110      !!                   + mj(mk(pahv)) * e1t*wslpj dj[ mj(mk(tb)) ]  }
111      !!      take the horizontal divergence of the fluxes:
112      !!         difft = 1/(e1e2t*e3t) dk[ zftw ]
113      !!      Add this trend to the general trend (ta,sa):
114      !!         pt_rhs = pt_rhs + difft
115      !!
116      !! ** Action :   Update pt_rhs arrays with the before rotated diffusion
117      !!----------------------------------------------------------------------
118      INTEGER                                   , INTENT(in   ) ::   kt         ! ocean time-step index
119      INTEGER                                   , INTENT(in   ) ::   kit000     ! first time step index
120      CHARACTER(len=3)                          , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
121      INTEGER                                   , INTENT(in   ) ::   kjpt       ! number of tracers
122      INTEGER                                   , INTENT(in   ) ::   kpass      ! =1/2 first or second passage
123      INTEGER                                   , INTENT(in   ) ::   Kmm        ! ocean time level index
124      INTEGER                                   , INTENT(in   ) ::   ktah, ktg, ktgi, ktt, ktt2, ktt_rhs
125      REAL(wp), DIMENSION(A2D_T(ktah)   ,JPK)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s]
126      REAL(wp), DIMENSION(A2D_T(ktg)        ,KJPT), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels
127      REAL(wp), DIMENSION(A2D_T(ktgi)       ,KJPT), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels
128      REAL(wp), DIMENSION(A2D_T(ktt)    ,JPK,KJPT), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2)
129      REAL(wp), DIMENSION(A2D_T(ktt2)   ,JPK,KJPT), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2)
130      REAL(wp), DIMENSION(A2D_T(ktt_rhs),JPK,KJPT), INTENT(inout) ::   pt_rhs     ! tracer trend
131      !
132      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
133      INTEGER  ::  ikt
134      INTEGER  ::  ierr, iij        ! local integer
135      REAL(wp) ::  zmsku, zahu_w, zabe1, zcof1, zcoef3   ! local scalars
136      REAL(wp) ::  zmskv, zahv_w, zabe2, zcof2, zcoef4   !   -      -
137      REAL(wp) ::  zcoef0, ze3w_2, zsign                 !   -      -
138      REAL(wp), DIMENSION(A2D(nn_hls))     ::   zdkt, zdk1t, z2d
139      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zdit, zdjt, zftu, zftv, ztfw
140      !!----------------------------------------------------------------------
141      !
142      IF( kpass == 1 .AND. kt == kit000 )  THEN
143         IF( .NOT. l_istiled .OR. ntile == 1 )  THEN                       ! Do only on the first tile
144            IF(lwp) WRITE(numout,*)
145            IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype
146            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
147         ENDIF
148         !
149         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpk )
150         DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpk )
151            akz     (ji,jj,jk) = 0._wp
152            ah_wslp2(ji,jj,jk) = 0._wp
153         END_3D
154      ENDIF
155      !
156      IF( .NOT. l_istiled .OR. ntile == 1 )  THEN                           ! Do only on the first tile
157         l_hst = .FALSE.
158         l_ptr = .FALSE.
159         IF( cdtype == 'TRA' .AND. ( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf' ) ) )     l_ptr = .TRUE.
160         IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. &
161            &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )   l_hst = .TRUE.
162      ENDIF
163      !
164      ! Define pt_rhs halo points for multi-point haloes in bilaplacian case
165      IF( nldf_tra == np_blp_i .AND. kpass == 1 ) THEN ; iij = nn_hls
166      ELSE                                             ; iij = 1
167      ENDIF
168
169      !
170      IF( kpass == 1 ) THEN   ;   zsign =  1._wp      ! bilaplacian operator require a minus sign (eddy diffusivity >0)
171      ELSE                    ;   zsign = -1._wp
172      ENDIF
173
174      !!----------------------------------------------------------------------
175      !!   0 - calculate  ah_wslp2 and akz
176      !!----------------------------------------------------------------------
177      !
178      IF( kpass == 1 ) THEN                  !==  first pass only  ==!
179         !
180         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
181         DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
182            !
183            zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          &
184               &                           + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk) , 1._wp  )
185            zmskv = wmask(ji,jj,jk) / MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)          &
186               &                           + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk) , 1._wp  )
187               !
188            ! NOTE: [halo1-halo2]
189            ! Extra brackets required to ensure consistent floating point arithmetic for different nn_hls for bilaplacian
190            zahu_w = (   (pahu(ji  ,jj,jk-1) + pahu(ji-1,jj,jk))    &
191               &       + (pahu(ji-1,jj,jk-1) + pahu(ji  ,jj,jk))  ) * zmsku
192            zahv_w = (   (pahv(ji,jj  ,jk-1) + pahv(ji,jj-1,jk))    &
193               &       + (pahv(ji,jj-1,jk-1) + pahv(ji,jj  ,jk))  ) * zmskv
194               !
195            ah_wslp2(ji,jj,jk) = zahu_w * wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
196               &               + zahv_w * wslpj(ji,jj,jk) * wslpj(ji,jj,jk)
197         END_3D
198         !
199         IF( ln_traldf_msc ) THEN                ! stabilizing vertical diffusivity coefficient
200            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
201            DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
202               ! NOTE: [halo1-halo2]
203               ! Extra brackets required to ensure consistent floating point arithmetic for different nn_hls for bilaplacian
204               akz(ji,jj,jk) = 0.25_wp * (                                                                       &
205                  &              (( pahu(ji  ,jj,jk) + pahu(ji  ,jj,jk-1) ) / ( e1u(ji  ,jj) * e1u(ji  ,jj) )    &
206                  &            + (  pahu(ji-1,jj,jk) + pahu(ji-1,jj,jk-1) ) / ( e1u(ji-1,jj) * e1u(ji-1,jj) ))   &
207                  &            + (( pahv(ji,jj  ,jk) + pahv(ji,jj  ,jk-1) ) / ( e2v(ji,jj  ) * e2v(ji,jj  ) )    &
208                  &            + (  pahv(ji,jj-1,jk) + pahv(ji,jj-1,jk-1) ) / ( e2v(ji,jj-1) * e2v(ji,jj-1) ))   )
209            END_3D
210            !
211            IF( ln_traldf_blp ) THEN                ! bilaplacian operator
212               ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
213               DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
214                  akz(ji,jj,jk) = 16._wp   &
215                     &   * ah_wslp2   (ji,jj,jk)   &
216                     &   * (  akz     (ji,jj,jk)   &
217                     &      + ah_wslp2(ji,jj,jk)   &
218                     &        / ( e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm) )  )
219               END_3D
220            ELSEIF( ln_traldf_lap ) THEN              ! laplacian operator
221               ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
222               DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
223                  ze3w_2 = e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm)
224                  zcoef0 = rDt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  )
225                  akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_Dt
226               END_3D
227           ENDIF
228           !
229         ELSE                                    ! 33 flux set to zero with akz=ah_wslp2 ==>> computed in full implicit
230            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpk )
231            DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpk )
232               akz(ji,jj,jk) = ah_wslp2(ji,jj,jk)
233            END_3D
234         ENDIF
235      ENDIF
236      !
237      !                                                          ! ===========
238      DO jn = 1, kjpt                                            ! tracer loop
239         !                                                       ! ===========
240         !
241         !!----------------------------------------------------------------------
242         !!   I - masked horizontal derivative
243         !!----------------------------------------------------------------------
244         zdit(:,:,:) = 0._wp
245         zdjt(:,:,:) = 0._wp
246
247         ! Horizontal tracer gradient
248         ! [comm_cleanup] ! DO_3D( 1, 0, 1, 0, 1, jpkm1 )
249         DO_3D( iij, iij-1, iij, iij-1, 1, jpkm1 )
250            zdit(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn) - pt(ji,jj,jk,jn) ) * umask(ji,jj,jk)
251            zdjt(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn) - pt(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
252         END_3D
253         IF( ln_zps ) THEN      ! botton and surface ocean correction of the horizontal gradient
254            ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )           ! bottom correction (partial bottom cell)
255            DO_2D( iij, iij-1, iij, iij-1 )            ! bottom correction (partial bottom cell)
256               zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)
257               zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)
258            END_2D
259            IF( ln_isfcav ) THEN      ! first wet level beneath a cavity
260               ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
261               DO_2D( iij, iij-1, iij, iij-1 )
262                  IF( miku(ji,jj) > 1 )   zdit(ji,jj,miku(ji,jj)) = pgui(ji,jj,jn)
263                  IF( mikv(ji,jj) > 1 )   zdjt(ji,jj,mikv(ji,jj)) = pgvi(ji,jj,jn)
264               END_2D
265            ENDIF
266         ENDIF
267         !
268         !!----------------------------------------------------------------------
269         !!   II - horizontal trend  (full)
270         !!----------------------------------------------------------------------
271         !
272         DO jk = 1, jpkm1                                 ! Horizontal slab
273            !
274            ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
275            DO_2D( iij, iij, iij, iij )
276               !                             !== Vertical tracer gradient
277               zdk1t(ji,jj) = ( pt(ji,jj,jk,jn) - pt(ji,jj,jk+1,jn) ) * wmask(ji,jj,jk+1)     ! level jk+1
278               !
279               IF( jk == 1 ) THEN   ;   zdkt(ji,jj) = zdk1t(ji,jj)                            ! surface: zdkt(jk=1)=zdkt(jk=2)
280               ELSE                 ;   zdkt(ji,jj) = ( pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) ) * wmask(ji,jj,jk)
281               ENDIF
282            END_2D
283            !
284            ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )           !==  Horizontal fluxes
285            DO_2D( iij, iij-1, iij, iij-1 )           !==  Horizontal fluxes
286               zabe1 = pahu(ji,jj,jk) * e2_e1u(ji,jj) * e3u(ji,jj,jk,Kmm)
287               zabe2 = pahv(ji,jj,jk) * e1_e2v(ji,jj) * e3v(ji,jj,jk,Kmm)
288               !
289               zmsku = 1. / MAX(  wmask(ji+1,jj,jk  ) + wmask(ji,jj,jk+1)   &
290                  &             + wmask(ji+1,jj,jk+1) + wmask(ji,jj,jk  ), 1. )
291               !
292               zmskv = 1. / MAX(  wmask(ji,jj+1,jk  ) + wmask(ji,jj,jk+1)   &
293                  &             + wmask(ji,jj+1,jk+1) + wmask(ji,jj,jk  ), 1. )
294               !
295               zcof1 = - pahu(ji,jj,jk) * e2u(ji,jj) * uslp(ji,jj,jk) * zmsku
296               zcof2 = - pahv(ji,jj,jk) * e1v(ji,jj) * vslp(ji,jj,jk) * zmskv
297               !
298               ! NOTE: [halo1-halo2]
299               ! Extra brackets required to ensure consistent floating point arithmetic for different nn_hls for bilaplacian
300               zftu(ji,jj,jk) = (  zabe1 * zdit(ji,jj,jk)                          &
301                  &              + zcof1 * (  (zdkt (ji+1,jj) + zdk1t(ji,jj))      &
302                  &                         + (zdk1t(ji+1,jj) + zdkt (ji,jj))  )  ) * umask(ji,jj,jk)
303               zftv(ji,jj,jk) = (  zabe2 * zdjt(ji,jj,jk)                          &
304                  &              + zcof2 * (  (zdkt (ji,jj+1) + zdk1t(ji,jj))      &
305                  &                         + (zdk1t(ji,jj+1) + zdkt (ji,jj))  )  ) * vmask(ji,jj,jk)
306            END_2D
307            !
308            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )           !== horizontal divergence and add to pta
309            DO_2D( iij-1, iij-1, iij-1, iij-1 )           !== horizontal divergence and add to pta
310               ! NOTE: [halo1-halo2]
311               ! Extra brackets required to ensure consistent floating point arithmetic for different nn_hls for bilaplacian
312               pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)    &
313                  &       + zsign * (  (zftu(ji,jj,jk) - zftu(ji-1,jj,jk)) + &
314                  &                    (zftv(ji,jj,jk) - zftv(ji,jj-1,jk))  ) * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
315            END_2D
316         END DO                                        !   End of slab
317
318         !!----------------------------------------------------------------------
319         !!   III - vertical trend (full)
320         !!----------------------------------------------------------------------
321         !
322         ! Vertical fluxes
323         ! ---------------
324         !                          ! Surface and bottom vertical fluxes set to zero
325         ztfw(:,:, 1 ) = 0._wp      ;      ztfw(:,:,jpk) = 0._wp
326
327         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )    ! interior (2=<jk=<jpk-1)
328         DO_3D( iij-1, iij-1, iij-1, iij-1, 2, jpkm1 )    ! interior (2=<jk=<jpk-1)
329            !
330            zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          &
331               &                           + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk) , 1._wp  )
332            zmskv = wmask(ji,jj,jk) / MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)          &
333               &                           + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk) , 1._wp  )
334               !
335            zahu_w = (   pahu(ji  ,jj,jk-1) + pahu(ji-1,jj,jk)    &
336               &       + pahu(ji-1,jj,jk-1) + pahu(ji  ,jj,jk)  ) * zmsku
337            zahv_w = (   pahv(ji,jj  ,jk-1) + pahv(ji,jj-1,jk)    &
338               &       + pahv(ji,jj-1,jk-1) + pahv(ji,jj  ,jk)  ) * zmskv
339               !
340            zcoef3 = - zahu_w * e2t(ji,jj) * zmsku * wslpi (ji,jj,jk)   !wslpi & j are already w-masked
341            zcoef4 = - zahv_w * e1t(ji,jj) * zmskv * wslpj (ji,jj,jk)
342            !
343            ! NOTE: [halo1-halo2]
344            ! Extra brackets required to ensure consistent floating point arithmetic for different nn_hls for bilaplacian
345            ztfw(ji,jj,jk) = zcoef3 * (   (zdit(ji  ,jj  ,jk-1) + zdit(ji-1,jj  ,jk))      &
346               &                        + (zdit(ji-1,jj  ,jk-1) + zdit(ji  ,jj  ,jk))  )   &
347               &           + zcoef4 * (   (zdjt(ji  ,jj  ,jk-1) + zdjt(ji  ,jj-1,jk))      &
348               &                        + (zdjt(ji  ,jj-1,jk-1) + zdjt(ji  ,jj  ,jk))  )
349         END_3D
350         !                                !==  add the vertical 33 flux  ==!
351         IF( ln_traldf_lap ) THEN               ! laplacian case: eddy coef = ah_wslp2 - akz
352            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
353            DO_3D( iij-1, iij-1, iij-1, iij-1, 2, jpkm1 )
354               ztfw(ji,jj,jk) = ztfw(ji,jj,jk) + e1e2t(ji,jj) / e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)   &
355                  &                            * ( ah_wslp2(ji,jj,jk) - akz(ji,jj,jk) )               &
356                  &                            * (  pt(ji,jj,jk-1,jn) -  pt(ji,jj,jk,jn) )
357            END_3D
358            !
359         ELSE                                   ! bilaplacian
360            SELECT CASE( kpass )
361            CASE(  1  )                            ! 1st pass : eddy coef = ah_wslp2
362               ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
363               DO_3D( iij-1, iij-1, iij-1, iij-1, 2, jpkm1 )
364                  ztfw(ji,jj,jk) =   &
365                     &  ztfw(ji,jj,jk) + ah_wslp2(ji,jj,jk) * e1e2t(ji,jj)   &
366                     &           * ( pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) ) / e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
367               END_3D
368            CASE(  2  )                         ! 2nd pass : eddy flux = ah_wslp2 and akz applied on pt  and pt2 gradients, resp.
369               DO_3D( 0, 0, 0, 0, 2, jpkm1 )
370                  ztfw(ji,jj,jk) = ztfw(ji,jj,jk) + e1e2t(ji,jj) / e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)                  &
371                     &                            * (  ah_wslp2(ji,jj,jk) * ( pt (ji,jj,jk-1,jn) - pt (ji,jj,jk,jn) )   &
372                     &                            +         akz(ji,jj,jk) * ( pt2(ji,jj,jk-1,jn) - pt2(ji,jj,jk,jn) )   )
373               END_3D
374            END SELECT
375         ENDIF
376         !
377         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !==  Divergence of vertical fluxes added to pta  ==!
378         DO_3D( iij-1, iij-1, iij-1, iij-1, 1, jpkm1 )    !==  Divergence of vertical fluxes added to pta  ==!
379            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn) + zsign * (  ztfw (ji,jj,jk) - ztfw(ji,jj,jk+1)  ) * r1_e1e2t(ji,jj)   &
380               &                                             / e3t(ji,jj,jk,Kmm)
381         END_3D
382         !
383         IF( ( kpass == 1 .AND. ln_traldf_lap ) .OR.  &     !==  first pass only (  laplacian)  ==!
384             ( kpass == 2 .AND. ln_traldf_blp ) ) THEN      !==  2nd   pass      (bilaplacian)  ==!
385            !
386            !                             ! "Poleward" diffusive heat or salt transports (T-S case only)
387               ! note sign is reversed to give down-gradient diffusive transports )
388            IF( l_ptr )  CALL dia_ptr_hst( jn, 'ldf', -zftv(:,:,:)  )
389            !                          ! Diffusive heat transports
390            IF( l_hst )  CALL dia_ar5_hst( jn, 'ldf', -zftu(:,:,:), -zftv(:,:,:) )
391            !
392         ENDIF                                                    !== end pass selection  ==!
393         !
394         !                                                        ! ===============
395      END DO                                                      ! end tracer loop
396      !
397   END SUBROUTINE tra_ldf_iso_t
398
399   !!==============================================================================
400END MODULE traldf_iso
Note: See TracBrowser for help on using the repository browser.