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_triad.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_triad.F90 @ 14713

Last change on this file since 14713 was 14713, checked in by hadcv, 4 years ago

#2600: Resolve uninitialised arrays in traldf_triad.F90 from [14667]

  • Property svn:keywords set to Id
File size: 35.2 KB
Line 
1MODULE traldf_triad
2   !!======================================================================
3   !!                   ***  MODULE  traldf_triad  ***
4   !! Ocean  tracers:  horizontal component of the lateral tracer mixing trend
5   !!======================================================================
6   !! History :  3.3  ! 2010-10  (G. Nurser, C. Harris, G. Madec)  Griffies operator (original code)
7   !!            3.7  ! 2013-12  (F. Lemarie, G. Madec)  triad operator (Griffies) + Method of Stabilizing Correction
8   !!----------------------------------------------------------------------
9
10   !!----------------------------------------------------------------------
11   !!   tra_ldf_triad : update the tracer trend with the iso-neutral laplacian triad-operator
12   !!----------------------------------------------------------------------
13   USE oce            ! ocean dynamics and active tracers
14   USE dom_oce        ! ocean space and time domain
15   USE domutl, ONLY : is_tile
16   USE phycst         ! physical constants
17   USE trc_oce        ! share passive tracers/Ocean variables
18   USE zdf_oce        ! ocean vertical physics
19   USE ldftra         ! lateral physics: eddy diffusivity
20   USE ldfslp         ! lateral physics: iso-neutral slopes
21   USE traldf_iso     ! lateral diffusion (Madec operator)         (tra_ldf_iso routine)
22   USE diaptr         ! poleward transport diagnostics
23   USE diaar5         ! AR5 diagnostics
24   USE zpshde         ! partial step: hor. derivative     (zps_hde routine)
25   !
26   USE in_out_manager ! I/O manager
27   USE iom            ! I/O library
28   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
29   USE lib_mpp        ! MPP library
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC   tra_ldf_triad   ! routine called by traldf.F90
35
36   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
37   LOGICAL  ::   l_hst   ! flag to compute heat transport
38
39
40   !! * Substitutions
41#  include "do_loop_substitute.h90"
42#  include "domzgr_substitute.h90"
43   !!----------------------------------------------------------------------
44   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
45   !! $Id$
46   !! Software governed by the CeCILL license (see ./LICENSE)
47   !!----------------------------------------------------------------------
48CONTAINS
49
50   SUBROUTINE tra_ldf_triad( kt, Kmm, kit000, cdtype, pahu, pahv,             &
51      &                                               pgu , pgv , pgui, pgvi, &
52      &                                               pt, pt2, pt_rhs, kjpt, kpass )
53      !!
54      INTEGER                     , INTENT(in   ) ::   kt         ! ocean time-step index
55      INTEGER                     , INTENT(in   ) ::   kit000     ! first time step index
56      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
57      INTEGER                     , INTENT(in   ) ::   kjpt       ! number of tracers
58      INTEGER                     , INTENT(in   ) ::   kpass      ! =1/2 first or second passage
59      INTEGER                     , INTENT(in   ) ::   Kmm        ! ocean time level indices
60      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s]
61      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgu , pgv  ! tracer gradient at pstep levels
62      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels
63      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2)
64      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2)
65      REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) ::   pt_rhs     ! tracer trend
66      !!
67      CALL tra_ldf_triad_t( kt, Kmm, kit000, cdtype, pahu, pahv, is_tile(pahu),                            &
68      &                                              pgu , pgv , is_tile(pgu) , pgui, pgvi, is_tile(pgui), &
69      &                                              pt, is_tile(pt), pt2, is_tile(pt2), pt_rhs, is_tile(pt_rhs), kjpt, kpass )
70   END SUBROUTINE tra_ldf_triad
71
72
73  SUBROUTINE tra_ldf_triad_t( kt, Kmm, kit000, cdtype, pahu, pahv, ktah,                   &
74      &                                                pgu , pgv , ktg , pgui, pgvi, ktgi, &
75      &                                                pt, ktt, pt2, ktt2, pt_rhs, ktt_rhs, kjpt, kpass )
76      !!----------------------------------------------------------------------
77      !!                  ***  ROUTINE tra_ldf_triad  ***
78      !!
79      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive
80      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and
81      !!      add it to the general trend of tracer equation.
82      !!
83      !! ** Method  :   The horizontal component of the lateral diffusive trends
84      !!      is provided by a 2nd order operator rotated along neural or geopo-
85      !!      tential surfaces to which an eddy induced advection can be added
86      !!      It is computed using before fields (forward in time) and isopyc-
87      !!      nal or geopotential slopes computed in routine ldfslp.
88      !!
89      !!      see documentation for the desciption
90      !!
91      !! ** Action :   pt_rhs   updated with the before rotated diffusion
92      !!               ah_wslp2 ....
93      !!               akz   stabilizing vertical diffusivity coefficient (used in trazdf_imp)
94      !!----------------------------------------------------------------------
95      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index
96      INTEGER                              , INTENT(in   ) ::   kit000     ! first time step index
97      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
98      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers
99      INTEGER                              , INTENT(in   ) ::   kpass      ! =1/2 first or second passage
100      INTEGER                              , INTENT(in)    ::   Kmm        ! ocean time level indices
101      INTEGER                              , INTENT(in   ) ::   ktah, ktg, ktgi, ktt, ktt2, ktt_rhs
102      REAL(wp), DIMENSION(A2D_T(ktah),   JPK)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s]
103      REAL(wp), DIMENSION(A2D_T(ktg),        KJPT), INTENT(in   ) ::   pgu , pgv  ! tracer gradient at pstep levels
104      REAL(wp), DIMENSION(A2D_T(ktgi),       KJPT), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels
105      REAL(wp), DIMENSION(A2D_T(ktt),    JPK,KJPT), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2)
106      REAL(wp), DIMENSION(A2D_T(ktt2),   JPK,KJPT), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2)
107      REAL(wp), DIMENSION(A2D_T(ktt_rhs),JPK,KJPT), INTENT(inout) ::   pt_rhs     ! tracer trend
108      !
109      INTEGER  ::  ji, jj, jk, jn, kp, iij   ! dummy loop indices
110      REAL(wp) ::  zcoef0, ze3w_2, zsign          !   -      -
111      !
112      REAL(wp) ::   zslope2, zbu, zbv, zbu1, zbv1, zslope21, zah, zah1, zah_ip1, zah_jp1, zbu_ip1, zbv_jp1
113      REAL(wp) ::   ze1ur, ze2vr, ze3wr, zdxt, zdyt, zdzt, zdyt_jp1, ze3wr_jp1, zdzt_jp1, zah_slp1, zah_slp_jp1, zaei_slp_jp1
114      REAL(wp) ::   zah_slp, zaei_slp, zdxt_ip1, ze3wr_ip1, zdzt_ip1, zah_slp_ip1, zaei_slp_ip1, zaei_slp1
115      REAL(wp), DIMENSION(A2D(nn_hls),0:1) ::   zdkt3d                                           ! vertical tracer gradient at 2 levels
116      REAL(wp), DIMENSION(A2D(nn_hls)    ) ::   z2d                                              ! 2D workspace
117      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zdit, zdjt, zftu, zftv, ztfw, zpsi_uw, zpsi_vw   ! 3D     -
118      !!----------------------------------------------------------------------
119      !
120      IF( .NOT. l_istiled .OR. ntile == 1 )  THEN                       ! Do only on the first tile
121         IF( kpass == 1 .AND. kt == kit000 )  THEN
122            IF(lwp) WRITE(numout,*)
123            IF(lwp) WRITE(numout,*) 'tra_ldf_triad : rotated laplacian diffusion operator on ', cdtype
124            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~'
125         ENDIF
126         !
127         l_hst = .FALSE.
128         l_ptr = .FALSE.
129         IF( cdtype == 'TRA' ) THEN
130            IF( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf') )      l_ptr = .TRUE.
131            IF( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.                   &
132            &   iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  )   l_hst = .TRUE.
133         ENDIF
134      ENDIF
135      !
136      ! Define pt_rhs halo points for multi-point haloes in bilaplacian case
137      IF( nldf_tra == np_blp_it .AND. kpass == 1 ) THEN ; iij = nn_hls
138      ELSE                                              ; iij = 1
139      ENDIF
140
141      !
142      IF( kpass == 1 ) THEN   ;   zsign =  1._wp      ! bilaplacian operator require a minus sign (eddy diffusivity >0)
143      ELSE                    ;   zsign = -1._wp
144      ENDIF
145      !
146      !!----------------------------------------------------------------------
147      !!   0 - calculate  ah_wslp2, akz, and optionally zpsi_uw, zpsi_vw
148      !!----------------------------------------------------------------------
149      !
150      IF( kpass == 1 ) THEN         !==  first pass only  and whatever the tracer is  ==!
151         !
152         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpk )
153         DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpk )
154            akz     (ji,jj,jk) = 0._wp
155            ah_wslp2(ji,jj,jk) = 0._wp
156         END_3D
157         !
158         DO kp = 0, 1                            ! i-k triads
159            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpkm1 )
160            DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )
161               ze3wr = 1._wp / e3w(ji,jj,jk+kp,Kmm)
162               zbu   = e1e2u(ji,jj) * e3u(ji,jj,jk,Kmm)
163               zbu1  = e1e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm)
164               zah   = 0.25_wp * pahu(ji,jj,jk)
165               zah1  = 0.25_wp * pahu(ji-1,jj,jk)
166               ! Subtract s-coordinate slope at t-points to give slope rel to s-surfaces (do this by *adding* gradient of depth)
167               zslope2 = triadi_g(ji,jj,jk,1,kp) + ( gdept(ji+1,jj,jk,Kmm) - gdept(ji,jj,jk,Kmm) ) * r1_e1u(ji,jj) * umask(ji,jj,jk+kp)
168               zslope2 = zslope2 *zslope2
169               zslope21 = triadi_g(ji,jj,jk,0,kp) + ( gdept(ji,jj,jk,Kmm) - gdept(ji-1,jj,jk,Kmm) ) * r1_e1u(ji-1,jj) * umask(ji-1,jj,jk+kp)
170               zslope21 = zslope21 *zslope21
171               ! round brackets added to fix the order of floating point operations
172               ! needed to ensure halo 1 - halo 2 compatibility
173               ah_wslp2(ji,jj,jk+kp) =  ah_wslp2(ji,jj,jk+kp) + ( zah * zbu * ze3wr * r1_e1e2t(ji,jj) * zslope2                    &
174                        &                                       + zah1 * zbu1 * ze3wr * r1_e1e2t(ji,jj) * zslope21                 &
175                        &                                       )                                                                  ! bracket for halo 1 - halo 2 compatibility
176               akz     (ji,jj,jk+kp) =  akz     (ji,jj,jk+kp) + ( zah * r1_e1u(ji,jj) * r1_e1u(ji,jj) * umask(ji,jj,jk+kp)         &
177                                                                + zah1 * r1_e1u(ji-1,jj) * r1_e1u(ji-1,jj) * umask(ji-1,jj,jk+kp)  &
178                        &                                       )                                                                  ! bracket for halo 1 - halo 2 compatibility
179            END_3D
180         END DO
181         !
182         DO kp = 0, 1                            ! j-k triads
183            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpkm1 )
184            DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )
185               ze3wr = 1.0_wp / e3w(ji,jj,jk+kp,Kmm)
186               zbv   = e1e2v(ji,jj) * e3v(ji,jj,jk,Kmm)
187               zbv1   = e1e2v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm)
188               zah   = 0.25_wp * pahv(ji,jj,jk)
189               zah1   = 0.25_wp * pahv(ji,jj-1,jk)
190               ! Subtract s-coordinate slope at t-points to give slope rel to s surfaces
191               !    (do this by *adding* gradient of depth)
192               zslope2 = triadj_g(ji,jj,jk,1,kp) + ( gdept(ji,jj+1,jk,Kmm) - gdept(ji,jj,jk,Kmm) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk+kp)
193               zslope2 = zslope2 * zslope2
194               zslope21 = triadj_g(ji,jj,jk,0,kp) + ( gdept(ji,jj,jk,Kmm) - gdept(ji,jj-1,jk,Kmm) ) * r1_e2v(ji,jj-1) * vmask(ji,jj-1,jk+kp)
195               zslope21 = zslope21 * zslope21
196               ! round brackets added to fix the order of floating point operations
197               ! needed to ensure halo 1 - halo 2 compatibility
198               ah_wslp2(ji,jj,jk+kp) = ah_wslp2(ji,jj,jk+kp) + ( zah * zbv * ze3wr * r1_e1e2t(ji,jj) * zslope2                     &
199                        &                                      + zah1 * zbv1 * ze3wr * r1_e1e2t(ji,jj) * zslope21                  &
200                        &                                      )                                                                   ! bracket for halo 1 - halo 2 compatibility
201               akz     (ji,jj,jk+kp) = akz     (ji,jj,jk+kp) + ( zah * r1_e2v(ji,jj) * r1_e2v(ji,jj) * vmask(ji,jj,jk+kp)          &
202                        &                                      + zah1 * r1_e2v(ji,jj-1) * r1_e2v(ji,jj-1) * vmask(ji,jj-1,jk+kp)   &
203                        &                                      )                                                                   ! bracket for halo 1 - halo 2 compatibility
204               !
205            END_3D
206         END DO
207         !
208         IF( ln_traldf_msc ) THEN                ! stabilizing vertical diffusivity coefficient
209            !
210            IF( ln_traldf_blp ) THEN                ! bilaplacian operator
211               ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
212               DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
213                  akz(ji,jj,jk) = 16._wp           &
214                     &   * ah_wslp2   (ji,jj,jk)   &
215                     &   * (  akz     (ji,jj,jk)   &
216                     &      + ah_wslp2(ji,jj,jk)   &
217                     &        / ( e3w (ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm) )  )
218               END_3D
219            ELSEIF( ln_traldf_lap ) THEN              ! laplacian operator
220               ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
221               DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
222                  ze3w_2 = e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm)
223                  zcoef0 = rDt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  )
224                  akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_Dt
225               END_3D
226           ENDIF
227           !
228         ELSE                                    ! 33 flux set to zero with akz=ah_wslp2 ==>> computed in full implicit
229            ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpk )
230            DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpk )
231               akz(ji,jj,jk) = ah_wslp2(ji,jj,jk)
232            END_3D
233         ENDIF
234         !
235         IF( ln_ldfeiv_dia .AND. cdtype == 'TRA' ) THEN
236            zpsi_uw(:,:,:) = 0._wp
237            zpsi_vw(:,:,:) = 0._wp
238
239            DO kp = 0, 1
240               DO_3D( 1, 0, 1, 0, 1, jpkm1 )
241                  ! round brackets added to fix the order of floating point operations
242                  ! needed to ensure halo 1 - halo 2 compatibility
243                  zpsi_uw(ji,jj,jk+kp) = zpsi_uw(ji,jj,jk+kp)                                     &
244                     & + ( 0.25_wp * aeiu(ji,jj,jk) * e2u(ji,jj) * triadi_g(ji,jj,jk,1,kp)        &
245                     &   + 0.25_wp * aeiu(ji,jj,jk) * e2u(ji,jj) * triadi_g(ji+1,jj,jk,0,kp)      &
246                     &   )                                                                        ! bracket for halo 1 - halo 2 compatibility
247                  zpsi_vw(ji,jj,jk+kp) = zpsi_vw(ji,jj,jk+kp)                                     &
248                     & + ( 0.25_wp * aeiv(ji,jj,jk) * e1v(ji,jj) * triadj_g(ji,jj,jk,1,kp)        &
249                     &   + 0.25_wp * aeiv(ji,jj,jk) * e1v(ji,jj) * triadj_g(ji,jj+1,jk,0,kp)      &
250                     &   )                                                                        ! bracket for halo 1 - halo 2 compatibility
251               END_3D
252            END DO
253            CALL ldf_eiv_dia( zpsi_uw, zpsi_vw, Kmm )
254         ENDIF
255         !
256      ENDIF                                  !==  end 1st pass only  ==!
257      !
258      !                                                           ! ===========
259      DO jn = 1, kjpt                                             ! tracer loop
260         !                                                        ! ===========
261         ! Zero fluxes for each tracer
262!!gm  this should probably be done outside the jn loop
263         ztfw(:,:,:) = 0._wp
264         zftu(:,:,:) = 0._wp
265         zftv(:,:,:) = 0._wp
266         ! NOTE: [tiling] these are zeroed to avoid floating point exceptions due to undefined values when calculating zdxt_ip1 & zdyt_jp1
267         zdit(:,:,:) = 0._wp
268         zdjt(:,:,:) = 0._wp
269         !
270         ! [comm_cleanup] ! DO_3D( 1, 0, 1, 0, 1, jpkm1 )    !==  before lateral T & S gradients at T-level jk  ==!
271         DO_3D( iij, iij-1, iij, iij-1, 1, jpkm1 )    !==  before lateral T & S gradients at T-level jk  ==!
272            zdit(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn) - pt(ji,jj,jk,jn) ) * umask(ji,jj,jk)
273            zdjt(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn) - pt(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
274         END_3D
275         IF( ln_zps .AND. l_grad_zps ) THEN    ! partial steps: correction at top/bottom ocean level
276            ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )                    ! bottom level
277            DO_2D( iij, iij-1, iij, iij-1 )                    ! bottom level
278               zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)
279               zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)
280            END_2D
281            IF( ln_isfcav ) THEN                   ! top level (ocean cavities only)
282               ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
283               DO_2D( iij, iij-1, iij, iij-1 )
284                  IF( miku(ji,jj)  > 1 )   zdit(ji,jj,miku(ji,jj) ) = pgui(ji,jj,jn)
285                  IF( mikv(ji,jj)  > 1 )   zdjt(ji,jj,mikv(ji,jj) ) = pgvi(ji,jj,jn)
286               END_2D
287            ENDIF
288         ENDIF
289         !
290         !!----------------------------------------------------------------------
291         !!   II - horizontal trend  (full)
292         !!----------------------------------------------------------------------
293         !
294         DO jk = 1, jpkm1
295            !                    !==  Vertical tracer gradient at level jk and jk+1
296            ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
297            DO_2D( iij, iij, iij, iij )
298               zdkt3d(ji,jj,1) = ( pt(ji,jj,jk,jn) - pt(ji,jj,jk+1,jn) ) * tmask(ji,jj,jk+1)
299            END_2D
300            !
301            !                    ! surface boundary condition: zdkt3d(jk=0)=zdkt3d(jk=1)
302            IF( jk == 1 ) THEN   ;   zdkt3d(:,:,0) = zdkt3d(:,:,1)
303            ELSE
304               ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
305               DO_2D( iij, iij, iij, iij )
306                  zdkt3d(ji,jj,0) = ( pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) ) * tmask(ji,jj,jk)
307               END_2D
308            ENDIF
309            !
310            zaei_slp = 0._wp
311            zaei_slp_ip1 = 0._wp
312            zaei_slp_jp1 = 0._wp
313            zaei_slp1 = 0._wp
314            !
315            IF( ln_botmix_triad ) THEN
316               DO kp = 0, 1              !==  Horizontal & vertical fluxes
317                  ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
318                  DO_2D( iij, iij-1, iij, iij-1 )
319                     ze1ur = r1_e1u(ji,jj)
320                     zdxt  = zdit(ji,jj,jk) * ze1ur
321                     zdxt_ip1  = zdit(ji+1,jj,jk) * r1_e1u(ji+1,jj)
322                     ze3wr = 1._wp / e3w(ji,jj,jk+kp,Kmm)
323                     ze3wr_ip1 = 1._wp / e3w(ji+1,jj,jk+kp,Kmm)
324                     zdzt  = zdkt3d(ji,jj,kp) * ze3wr
325                     zdzt_ip1  = zdkt3d(ji+1,jj,kp) * ze3wr_ip1
326                     !
327                     zbu = 0.25_wp * e1e2u(ji,jj) * e3u(ji,jj,jk,Kmm)
328                     zbu_ip1 = 0.25_wp * e1e2u(ji+1,jj) * e3u(ji+1,jj,jk,Kmm)
329                     ! ln_botmix_triad is .T. don't mask zah for bottom half cells    !!gm ?????   ahu is masked....
330                     zah = pahu(ji,jj,jk)
331                     zah_ip1 = pahu(ji+1,jj,jk)
332                     zah_slp  = zah * triadi(ji,jj,jk,1,kp)
333                     zah_slp_ip1  = zah_ip1 * triadi(ji+1,jj,jk,1,kp)
334                     zah_slp1  = zah * triadi(ji+1,jj,jk,0,kp)
335                     IF( ln_ldfeiv )   THEN
336                        zaei_slp = aeiu(ji,jj,jk) * triadi_g(ji,jj,jk,1,kp)
337                        zaei_slp_ip1 = aeiu(ji+1,jj,jk) * triadi_g(ji+1,jj,jk,1,kp)
338                        zaei_slp1 = aeiu(ji,jj,jk) * triadi_g(ji+1,jj,jk,0,kp)
339                     ENDIF
340                     ! round brackets added to fix the order of floating point operations
341                     ! needed to ensure halo 1 - halo 2 compatibility
342                     zftu(ji   ,jj,jk  ) =  zftu(ji   ,jj,jk )                                                               &
343                                         &    - ( ( zah * zdxt + ( zah_slp - zaei_slp ) * zdzt ) * zbu * ze1ur               &
344                                         &      + ( zah * zdxt + zah_slp1 * zdzt_ip1 - zaei_slp1 * zdzt_ip1 ) * zbu * ze1ur  &
345                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
346                     ztfw(ji+1,jj,jk+kp) =  ztfw(ji+1,jj,jk+kp)                                                              &
347                                         &    - ( (zah_slp_ip1 + zaei_slp_ip1) * zdxt_ip1 * zbu_ip1 * ze3wr_ip1              &
348                                         &      + ( zah_slp1 + zaei_slp1) * zdxt * zbu * ze3wr_ip1                           &
349                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
350                  END_2D
351               END DO
352               !
353               DO kp = 0, 1
354                  ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
355                  DO_2D( iij, iij-1, iij, iij-1 )
356                     ze2vr = r1_e2v(ji,jj)
357                     zdyt  = zdjt(ji,jj,jk) * ze2vr
358                     zdyt_jp1  = zdjt(ji,jj+1,jk) * r1_e2v(ji,jj+1)
359                     ze3wr = 1._wp / e3w(ji,jj,jk+kp,Kmm)
360                     ze3wr_jp1 = 1._wp / e3w(ji,jj+1,jk+kp,Kmm)
361                     zdzt  = zdkt3d(ji,jj,kp) * ze3wr
362                     zdzt_jp1  = zdkt3d(ji,jj+1,kp) * ze3wr_jp1
363                     zbv = 0.25_wp * e1e2v(ji,jj) * e3v(ji,jj,jk,Kmm)
364                     zbv_jp1 = 0.25_wp * e1e2v(ji,jj+1) * e3v(ji,jj+1,jk,Kmm)
365                     ! ln_botmix_triad is .T. don't mask zah for bottom half cells    !!gm ?????   ahu is masked....
366                     zah = pahv(ji,jj,jk)          ! pahv(ji,jj+jp,jk)  ????
367                     zah_jp1 = pahv(ji,jj+1,jk)
368                     zah_slp = zah * triadj(ji,jj,jk,1,kp)
369                     zah_slp1 = zah * triadj(ji,jj+1,jk,0,kp)
370                     zah_slp_jp1 = zah_jp1 * triadj(ji,jj+1,jk,1,kp)
371                     IF( ln_ldfeiv )   THEN
372                        zaei_slp = aeiv(ji,jj,jk) * triadj_g(ji,jj,jk,1,kp)
373                        zaei_slp_jp1 = aeiv(ji,jj+1,jk) * triadj_g(ji,jj+1,jk,1,kp)
374                        zaei_slp1 = aeiv(ji,jj,jk) * triadj_g(ji,jj+1,jk,0,kp)
375                     ENDIF
376                     ! round brackets added to fix the order of floating point operations
377                     ! needed to ensure halo 1 - halo 2 compatibility
378                     zftv(ji,jj  ,jk   ) =  zftv(ji,jj  ,jk   )                                                              &
379                                         &    - ( ( zah * zdyt + ( zah_slp - zaei_slp ) * zdzt ) * zbv * ze2vr               &
380                                         &      + ( zah * zdyt + zah_slp1 * zdzt_jp1 - zaei_slp1 * zdzt_jp1 ) * zbv * ze2vr  &
381                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
382                     ztfw(ji,jj+1,jk+kp) =  ztfw(ji,jj+1,jk+kp)                                                              &
383                                         &    - ( ( zah_slp_jp1 + zaei_slp_jp1) * zdyt_jp1 * zbv_jp1 * ze3wr_jp1             &
384                                         &      + ( zah_slp1 + zaei_slp1) * zdyt * zbv * ze3wr_jp1                           &
385                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
386                  END_2D
387               END DO
388               !
389            ELSE
390               !
391               DO kp = 0, 1               !==  Horizontal & vertical fluxes
392                  ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
393                  DO_2D( iij, iij-1, iij, iij-1 )
394                     ze1ur = r1_e1u(ji,jj)
395                     zdxt  = zdit(ji,jj,jk) * ze1ur
396                     zdxt_ip1  = zdit(ji+1,jj,jk) * r1_e1u(ji+1,jj)
397                     ze3wr = 1._wp / e3w(ji,jj,jk+kp,Kmm)
398                     ze3wr_ip1 = 1._wp / e3w(ji+1,jj,jk+kp,Kmm)
399                     zdzt  = zdkt3d(ji,jj,kp) * ze3wr
400                     zdzt_ip1  = zdkt3d(ji+1,jj,kp) * ze3wr_ip1
401                     !
402                     zbu = 0.25_wp * e1e2u(ji,jj) * e3u(ji,jj,jk,Kmm)
403                     zbu_ip1 = 0.25_wp * e1e2u(ji+1,jj) * e3u(ji+1,jj,jk,Kmm)
404                     ! ln_botmix_triad is .F. mask zah for bottom half cells
405                     zah = pahu(ji,jj,jk) * umask(ji,jj,jk+kp)         ! pahu(ji+ip,jj,jk)   ===>>  ????
406                     zah_ip1 = pahu(ji+1,jj,jk) * umask(ji+1,jj,jk+kp)
407                     zah_slp  = zah * triadi(ji,jj,jk,1,kp)
408                     zah_slp_ip1  = zah_ip1 * triadi(ji+1,jj,jk,1,kp)
409                     zah_slp1  = zah * triadi(ji+1,jj,jk,0,kp)
410                     IF( ln_ldfeiv )   THEN
411                        zaei_slp = aeiu(ji,jj,jk) * triadi_g(ji,jj,jk,1,kp)
412                        zaei_slp_ip1 = aeiu(ji+1,jj,jk) * triadi_g(ji+1,jj,jk,1,kp)
413                        zaei_slp1 = aeiu(ji,jj,jk) * triadi_g(ji+1,jj,jk,0,kp)
414                     ENDIF
415!                     zftu(ji   ,jj,jk  ) = zftu(ji   ,jj,jk ) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur - ( zah * zdxt + (zah_slp1 - zaei_slp1) * zdzt_ip1 ) * zbu * ze1ur
416!                     ztfw(ji+1,jj,jk+kp) = ztfw(ji+1,jj,jk+kp) - (zah_slp_ip1 + zaei_slp_ip1) * zdxt_ip1 * zbu_ip1 * ze3wr_ip1 - (zah_slp1 + zaei_slp1) * zdxt * zbu * ze3wr_ip1
417                     ! round brackets added to fix the order of floating point operations
418                     ! needed to ensure halo 1 - halo 2 compatibility
419                     zftu(ji   ,jj,jk  ) =  zftu(ji   ,jj,jk )                                                               &
420                                         &    - ( ( zah * zdxt + ( zah_slp - zaei_slp ) * zdzt ) * zbu * ze1ur               &
421                                         &      + ( zah * zdxt + zah_slp1 * zdzt_ip1 - zaei_slp1 * zdzt_ip1 ) * zbu * ze1ur  &
422                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
423                     ztfw(ji+1,jj,jk+kp) =  ztfw(ji+1,jj,jk+kp)                                                              &
424                                         &    - ( (zah_slp_ip1 + zaei_slp_ip1) * zdxt_ip1 * zbu_ip1 * ze3wr_ip1              &
425                                         &      + ( zah_slp1 + zaei_slp1) * zdxt * zbu * ze3wr_ip1                           &
426                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
427                  END_2D
428               END DO
429               !
430               DO kp = 0, 1
431                  ! [comm_cleanup] ! DO_2D( 1, 0, 1, 0 )
432                  DO_2D( iij, iij-1, iij, iij-1 )
433                     ze2vr = r1_e2v(ji,jj)
434                     zdyt  = zdjt(ji,jj,jk) * ze2vr
435                     zdyt_jp1  = zdjt(ji,jj+1,jk) * r1_e2v(ji,jj+1)
436                     ze3wr = 1._wp / e3w(ji,jj,jk+kp,Kmm)
437                     ze3wr_jp1 = 1._wp / e3w(ji,jj+1,jk+kp,Kmm)
438                     zdzt  = zdkt3d(ji,jj,kp) * ze3wr
439                     zdzt_jp1  = zdkt3d(ji,jj+1,kp) * ze3wr_jp1
440                     zbv = 0.25_wp * e1e2v(ji,jj) * e3v(ji,jj,jk,Kmm)
441                     zbv_jp1 = 0.25_wp * e1e2v(ji,jj+1) * e3v(ji,jj+1,jk,Kmm)
442                     ! ln_botmix_triad is .F. mask zah for bottom half cells
443                     zah = pahv(ji,jj,jk) * vmask(ji,jj,jk+kp)         ! pahv(ji,jj+jp,jk)  ????
444                     zah_jp1 = pahv(ji,jj+1,jk) * vmask(ji,jj+1,jk+kp)
445                     zah_slp = zah * triadj(ji,jj,jk,1,kp)
446                     zah_slp1 = zah * triadj(ji,jj+1,jk,0,kp)
447                     zah_slp_jp1 = zah_jp1 * triadj(ji,jj+1,jk,1,kp)
448                     IF( ln_ldfeiv )   THEN
449                        zaei_slp = aeiv(ji,jj,jk) * triadj_g(ji,jj,jk,1,kp)
450                        zaei_slp_jp1 = aeiv(ji,jj+1,jk) * triadj_g(ji,jj+1,jk,1,kp)
451                        zaei_slp1 = aeiv(ji,jj,jk) * triadj_g(ji,jj+1,jk,0,kp)
452                     ENDIF
453!                     zftv(ji,jj  ,jk   ) = zftv(ji,jj  ,jk   ) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr - ( zah * zdyt + (zah_slp1 - zaei_slp1) * zdzt_jp1 ) * zbv * ze2vr
454!                     ztfw(ji,jj+1,jk+kp) = ztfw(ji,jj+1,jk+kp) - ( zah_slp_jp1 + zaei_slp_jp1) * zdyt_jp1 * zbv_jp1 * ze3wr_jp1 - (zah_slp1 + zaei_slp1) * zdyt * zbv * ze3wr_jp1
455                     ! round brackets added to fix the order of floating point operations
456                     ! needed to ensure halo 1 - halo 2 compatibility
457                     zftv(ji,jj  ,jk   ) =  zftv(ji,jj  ,jk   )                                                              &
458                                         &    - ( ( zah * zdyt + ( zah_slp - zaei_slp ) * zdzt ) * zbv * ze2vr               &
459                                         &      + ( zah * zdyt + zah_slp1 * zdzt_jp1 - zaei_slp1 * zdzt_jp1 ) * zbv * ze2vr  &
460                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
461                     ztfw(ji,jj+1,jk+kp) =  ztfw(ji,jj+1,jk+kp)                                                              &
462                                         &    - ( ( zah_slp_jp1 + zaei_slp_jp1) * zdyt_jp1 * zbv_jp1 * ze3wr_jp1             &
463                                         &      + ( zah_slp1 + zaei_slp1) * zdyt * zbv * ze3wr_jp1                           &
464                                         &      )                                                                            ! bracket for halo 1 - halo 2 compatibility
465                  END_2D
466               END DO
467            ENDIF
468            !                             !==  horizontal divergence and add to the general trend  ==!
469            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
470            DO_2D( iij-1, iij-1, iij-1, iij-1 )
471               ! round brackets added to fix the order of floating point operations
472               ! needed to ensure halo 1 - halo 2 compatibility
473               pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)                                                &
474                  &                       + zsign * ( ( zftu(ji-1,jj  ,jk) - zftu(ji,jj,jk)             &
475                  &                                   )                                                 & ! bracket for halo 1 - halo 2 compatibility
476                  &                                 + ( zftv(ji,jj-1,jk) - zftv(ji,jj,jk)               &
477                  &                                   )                                                 & ! bracket for halo 1 - halo 2 compatibility
478                  &                                 ) / (  e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm)  )
479            END_2D
480            !
481         END DO
482         !
483         !                                !==  add the vertical 33 flux  ==!
484         IF( ln_traldf_lap ) THEN               ! laplacian case: eddy coef = ah_wslp2 - akz
485            ! [comm_cleanup] ! DO_3D( 0, 0, 1, 0, 2, jpkm1 )
486            DO_3D( iij-1, iij-1, iij-1, iij-1, 2, jpkm1 )
487               ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / e3w(ji,jj,jk,Kmm) * tmask(ji,jj,jk)   &
488                  &                            * ( ah_wslp2(ji,jj,jk) - akz(ji,jj,jk) )             &
489                  &                            * (  pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) )
490            END_3D
491         ELSE                                   ! bilaplacian
492            SELECT CASE( kpass )
493            CASE(  1  )                            ! 1st pass : eddy coef = ah_wslp2
494               ! [comm_cleanup] ! DO_3D( 0, 0, 1, 0, 2, jpkm1 )
495               DO_3D( iij-1, iij-1, iij-1, iij-1, 2, jpkm1 )
496                  ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / e3w(ji,jj,jk,Kmm) * tmask(ji,jj,jk)             &
497                     &                            * ah_wslp2(ji,jj,jk) * ( pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) )
498               END_3D
499            CASE(  2  )                            ! 2nd pass : eddy flux = ah_wslp2 and akz applied on pt  and pt2 gradients, resp.
500               ! [comm_cleanup] ! DO_3D( 0, 0, 1, 0, 2, jpkm1 )
501               DO_3D( 0, 0, 0, 0, 2, jpkm1 )
502                  ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / e3w(ji,jj,jk,Kmm) * tmask(ji,jj,jk)                      &
503                     &                            * (  ah_wslp2(ji,jj,jk) * ( pt (ji,jj,jk-1,jn) - pt (ji,jj,jk,jn) )   &
504                     &                               + akz     (ji,jj,jk) * ( pt2(ji,jj,jk-1,jn) - pt2(ji,jj,jk,jn) )   )
505               END_3D
506            END SELECT
507         ENDIF
508         !
509         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpkm1 )      !==  Divergence of vertical fluxes added to pta  ==!
510         DO_3D( iij-1, iij-1, iij-1, iij-1, 1, jpkm1 )      !==  Divergence of vertical fluxes added to pta  ==!
511            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)    &
512            &                                  + zsign * (  ztfw(ji,jj,jk+1) - ztfw(ji,jj,jk)  )   &
513               &                                              / ( e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm) )
514         END_3D
515         !
516         IF( ( kpass == 1 .AND. ln_traldf_lap ) .OR.  &     !==  first pass only (  laplacian)  ==!
517             ( kpass == 2 .AND. ln_traldf_blp ) ) THEN      !==  2nd   pass      (bilaplacian)  ==!
518            !
519            !                          ! "Poleward" diffusive heat or salt transports (T-S case only)
520            IF( l_ptr )  CALL dia_ptr_hst( jn, 'ldf', zftv(:,:,:)  )
521            !                          ! Diffusive heat transports
522            IF( l_hst )  CALL dia_ar5_hst( jn, 'ldf', zftu(:,:,:), zftv(:,:,:) )
523            !
524         ENDIF                                                    !== end pass selection  ==!
525         !
526         !                                                        ! ===============
527      END DO                                                      ! end tracer loop
528      !                                                           ! ===============
529   END SUBROUTINE tra_ldf_triad_t
530
531   !!==============================================================================
532END MODULE traldf_triad
Note: See TracBrowser for help on using the repository browser.