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 branches/UKMO/r8727_WAVE-2_Clementi_add_coupling/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/UKMO/r8727_WAVE-2_Clementi_add_coupling/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90 @ 8749

Last change on this file since 8749 was 8749, checked in by jcastill, 6 years ago

Remove svn keywords

File size: 21.8 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 trc_oce        ! share passive tracers/Ocean variables
22   USE zdf_oce        ! ocean vertical physics
23   USE ldftra         ! lateral diffusion: tracer eddy coefficients
24   USE ldfslp         ! iso-neutral slopes
25   USE diaptr         ! poleward transport diagnostics
26   USE diaar5         ! AR5 diagnostics
27   !
28   USE in_out_manager ! I/O manager
29   USE iom            ! I/O library
30   USE phycst         ! physical constants
31   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
32   USE wrk_nemo       ! Memory Allocation
33   USE timing         ! Timing
34
35   IMPLICIT NONE
36   PRIVATE
37
38   PUBLIC   tra_ldf_iso   ! routine called by step.F90
39
40   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
41   LOGICAL  ::   l_hst   ! flag to compute heat transport
42
43   !! * Substitutions
44#  include "vectopt_loop_substitute.h90"
45   !!----------------------------------------------------------------------
46   !! NEMO/OPA 3.7 , NEMO Consortium (2015)
47   !! $Id$
48   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
49   !!----------------------------------------------------------------------
50CONTAINS
51
52  SUBROUTINE tra_ldf_iso( kt, kit000, cdtype, pahu, pahv, pgu , pgv ,   &
53      &                                                   pgui, pgvi,   &
54      &                                       ptb , ptbb, pta , kjpt, kpass )
55      !!----------------------------------------------------------------------
56      !!                  ***  ROUTINE tra_ldf_iso  ***
57      !!
58      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive
59      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and
60      !!      add it to the general trend of tracer equation.
61      !!
62      !! ** Method  :   The horizontal component of the lateral diffusive trends
63      !!      is provided by a 2nd order operator rotated along neural or geopo-
64      !!      tential surfaces to which an eddy induced advection can be added
65      !!      It is computed using before fields (forward in time) and isopyc-
66      !!      nal or geopotential slopes computed in routine ldfslp.
67      !!
68      !!      1st part :  masked horizontal derivative of T  ( di[ t ] )
69      !!      ========    with partial cell update if ln_zps=T
70      !!                  with top     cell update if ln_isfcav
71      !!
72      !!      2nd part :  horizontal fluxes of the lateral mixing operator
73      !!      ========   
74      !!         zftu =  pahu e2u*e3u/e1u di[ tb ]
75      !!               - pahu e2u*uslp    dk[ mi(mk(tb)) ]
76      !!         zftv =  pahv e1v*e3v/e2v dj[ tb ]
77      !!               - pahv e2u*vslp    dk[ mj(mk(tb)) ]
78      !!      take the horizontal divergence of the fluxes:
79      !!         difft = 1/(e1e2t*e3t) {  di-1[ zftu ] +  dj-1[ zftv ]  }
80      !!      Add this trend to the general trend (ta,sa):
81      !!         ta = ta + difft
82      !!
83      !!      3rd part: vertical trends of the lateral mixing operator
84      !!      ========  (excluding the vertical flux proportional to dk[t] )
85      !!      vertical fluxes associated with the rotated lateral mixing:
86      !!         zftw = - {  mi(mk(pahu)) * e2t*wslpi di[ mi(mk(tb)) ]
87      !!                   + mj(mk(pahv)) * e1t*wslpj dj[ mj(mk(tb)) ]  }
88      !!      take the horizontal divergence of the fluxes:
89      !!         difft = 1/(e1e2t*e3t) dk[ zftw ]
90      !!      Add this trend to the general trend (ta,sa):
91      !!         pta = pta + difft
92      !!
93      !! ** Action :   Update pta arrays with the before rotated diffusion
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      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s]
101      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels
102      REAL(wp), DIMENSION(jpi,jpj,    kjpt), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels
103      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb        ! tracer (kpass=1) or laplacian of tracer (kpass=2)
104      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptbb       ! tracer (only used in kpass=2)
105      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend
106      !
107      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
108      INTEGER  ::  ikt
109      INTEGER  ::  ierr             ! local integer
110      REAL(wp) ::  zmsku, zahu_w, zabe1, zcof1, zcoef3   ! local scalars
111      REAL(wp) ::  zmskv, zahv_w, zabe2, zcof2, zcoef4   !   -      -
112      REAL(wp) ::  zcoef0, ze3w_2, zsign, z2dt, z1_2dt   !   -      -
113      REAL(wp), POINTER, DIMENSION(:,:)   ::   zdkt, zdk1t, z2d
114      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zdit, zdjt, zftu, zftv, ztfw 
115      !!----------------------------------------------------------------------
116      !
117      IF( nn_timing == 1 )  CALL timing_start('tra_ldf_iso')
118      !
119      CALL wrk_alloc( jpi,jpj,       zdkt, zdk1t, z2d ) 
120      CALL wrk_alloc( jpi,jpj,jpk,   zdit, zdjt , zftu, zftv, ztfw  ) 
121      !
122      IF( kt == kit000 )  THEN
123         IF(lwp) WRITE(numout,*)
124         IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype
125         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
126         !
127         akz     (:,:,:) = 0._wp     
128         ah_wslp2(:,:,:) = 0._wp
129      ENDIF
130      !   
131      l_hst = .FALSE.
132      l_ptr = .FALSE.
133      IF( cdtype == 'TRA' .AND. ln_diaptr )                                                 l_ptr = .TRUE. 
134      IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. &
135         &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )   l_hst = .TRUE.
136      !
137      !                                            ! set time step size (Euler/Leapfrog)
138      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   z2dt =     rdt      ! at nit000   (Euler)
139      ELSE                                        ;   z2dt = 2.* rdt      !             (Leapfrog)
140      ENDIF
141      z1_2dt = 1._wp / z2dt
142      !
143      IF( kpass == 1 ) THEN   ;   zsign =  1._wp      ! bilaplacian operator require a minus sign (eddy diffusivity >0)
144      ELSE                    ;   zsign = -1._wp
145      ENDIF
146         
147      !!----------------------------------------------------------------------
148      !!   0 - calculate  ah_wslp2 and akz
149      !!----------------------------------------------------------------------
150      !
151      IF( kpass == 1 ) THEN                  !==  first pass only  ==!
152         !
153         DO jk = 2, jpkm1
154            DO jj = 2, jpjm1
155               DO ji = fs_2, fs_jpim1   ! vector opt.
156                  !
157                  zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          &
158                     &                           + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk) , 1._wp  )
159                  zmskv = wmask(ji,jj,jk) / MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)          &
160                     &                           + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk) , 1._wp  )
161                     !
162                  zahu_w = (   pahu(ji  ,jj,jk-1) + pahu(ji-1,jj,jk)    &
163                     &       + pahu(ji-1,jj,jk-1) + pahu(ji  ,jj,jk)  ) * zmsku
164                  zahv_w = (   pahv(ji,jj  ,jk-1) + pahv(ji,jj-1,jk)    &
165                     &       + pahv(ji,jj-1,jk-1) + pahv(ji,jj  ,jk)  ) * zmskv
166                     !
167                  ah_wslp2(ji,jj,jk) = zahu_w * wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
168                     &               + zahv_w * wslpj(ji,jj,jk) * wslpj(ji,jj,jk)
169               END DO
170            END DO
171         END DO
172         !
173         IF( ln_traldf_msc ) THEN                ! stabilizing vertical diffusivity coefficient
174            DO jk = 2, jpkm1
175               DO jj = 2, jpjm1
176                  DO ji = fs_2, fs_jpim1
177                     akz(ji,jj,jk) = 0.25_wp * (                                                                     &
178                        &              ( pahu(ji  ,jj,jk) + pahu(ji  ,jj,jk-1) ) / ( e1u(ji  ,jj) * e1u(ji  ,jj) )   &
179                        &            + ( pahu(ji-1,jj,jk) + pahu(ji-1,jj,jk-1) ) / ( e1u(ji-1,jj) * e1u(ji-1,jj) )   &
180                        &            + ( pahv(ji,jj  ,jk) + pahv(ji,jj  ,jk-1) ) / ( e2v(ji,jj  ) * e2v(ji,jj  ) )   &
181                        &            + ( pahv(ji,jj-1,jk) + pahv(ji,jj-1,jk-1) ) / ( e2v(ji,jj-1) * e2v(ji,jj-1) )   )
182                  END DO
183               END DO
184            END DO
185            !
186            IF( ln_traldf_blp ) THEN                ! bilaplacian operator
187               DO jk = 2, jpkm1
188                  DO jj = 1, jpjm1
189                     DO ji = 1, fs_jpim1
190                        akz(ji,jj,jk) = 16._wp * ah_wslp2(ji,jj,jk)   &
191                           &          * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ( e3w_n(ji,jj,jk) * e3w_n(ji,jj,jk) )  )
192                     END DO
193                  END DO
194               END DO
195            ELSEIF( ln_traldf_lap ) THEN              ! laplacian operator
196               DO jk = 2, jpkm1
197                  DO jj = 1, jpjm1
198                     DO ji = 1, fs_jpim1
199                        ze3w_2 = e3w_n(ji,jj,jk) * e3w_n(ji,jj,jk)
200                        zcoef0 = z2dt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  )
201                        akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * z1_2dt
202                     END DO
203                  END DO
204               END DO
205           ENDIF
206           !
207         ELSE                                    ! 33 flux set to zero with akz=ah_wslp2 ==>> computed in full implicit
208            akz(:,:,:) = ah_wslp2(:,:,:)     
209         ENDIF
210      ENDIF
211      !
212      !                                                          ! ===========
213      DO jn = 1, kjpt                                            ! tracer loop
214         !                                                       ! ===========
215         !                                               
216         !!----------------------------------------------------------------------
217         !!   I - masked horizontal derivative
218         !!----------------------------------------------------------------------
219!!gm : bug.... why (x,:,:)?   (1,jpj,:) and (jpi,1,:) should be sufficient....
220         zdit (1,:,:) = 0._wp     ;     zdit (jpi,:,:) = 0._wp
221         zdjt (1,:,:) = 0._wp     ;     zdjt (jpi,:,:) = 0._wp
222         !!end
223
224         ! Horizontal tracer gradient
225         DO jk = 1, jpkm1
226            DO jj = 1, jpjm1
227               DO ji = 1, fs_jpim1   ! vector opt.
228                  zdit(ji,jj,jk) = ( ptb(ji+1,jj  ,jk,jn) - ptb(ji,jj,jk,jn) ) * umask(ji,jj,jk)
229                  zdjt(ji,jj,jk) = ( ptb(ji  ,jj+1,jk,jn) - ptb(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
230               END DO
231            END DO
232         END DO
233         IF( ln_zps ) THEN      ! botton and surface ocean correction of the horizontal gradient
234            DO jj = 1, jpjm1              ! bottom correction (partial bottom cell)
235               DO ji = 1, fs_jpim1   ! vector opt.
236                  zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)         
237                  zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)
238               END DO
239            END DO
240            IF( ln_isfcav ) THEN      ! first wet level beneath a cavity
241               DO jj = 1, jpjm1
242                  DO ji = 1, fs_jpim1   ! vector opt.
243                     IF( miku(ji,jj) > 1 )   zdit(ji,jj,miku(ji,jj)) = pgui(ji,jj,jn)         
244                     IF( mikv(ji,jj) > 1 )   zdjt(ji,jj,mikv(ji,jj)) = pgvi(ji,jj,jn)     
245                  END DO
246               END DO
247            ENDIF
248         ENDIF
249         !
250         !!----------------------------------------------------------------------
251         !!   II - horizontal trend  (full)
252         !!----------------------------------------------------------------------
253         !
254         DO jk = 1, jpkm1                                 ! Horizontal slab
255            !
256            !                             !== Vertical tracer gradient
257            zdk1t(:,:) = ( ptb(:,:,jk,jn) - ptb(:,:,jk+1,jn) ) * wmask(:,:,jk+1)     ! level jk+1
258            !
259            IF( jk == 1 ) THEN   ;   zdkt(:,:) = zdk1t(:,:)                          ! surface: zdkt(jk=1)=zdkt(jk=2)
260            ELSE                 ;   zdkt(:,:) = ( ptb(:,:,jk-1,jn) - ptb(:,:,jk,jn) ) * wmask(:,:,jk)
261            ENDIF
262            DO jj = 1 , jpjm1            !==  Horizontal fluxes
263               DO ji = 1, fs_jpim1   ! vector opt.
264                  zabe1 = pahu(ji,jj,jk) * e2_e1u(ji,jj) * e3u_n(ji,jj,jk)
265                  zabe2 = pahv(ji,jj,jk) * e1_e2v(ji,jj) * e3v_n(ji,jj,jk)
266                  !
267                  zmsku = 1. / MAX(  wmask(ji+1,jj,jk  ) + wmask(ji,jj,jk+1)   &
268                     &             + wmask(ji+1,jj,jk+1) + wmask(ji,jj,jk  ), 1. )
269                  !
270                  zmskv = 1. / MAX(  wmask(ji,jj+1,jk  ) + wmask(ji,jj,jk+1)   &
271                     &             + wmask(ji,jj+1,jk+1) + wmask(ji,jj,jk  ), 1. )
272                  !
273                  zcof1 = - pahu(ji,jj,jk) * e2u(ji,jj) * uslp(ji,jj,jk) * zmsku
274                  zcof2 = - pahv(ji,jj,jk) * e1v(ji,jj) * vslp(ji,jj,jk) * zmskv
275                  !
276                  zftu(ji,jj,jk ) = (  zabe1 * zdit(ji,jj,jk)   &
277                     &               + zcof1 * (  zdkt (ji+1,jj) + zdk1t(ji,jj)      &
278                     &                          + zdk1t(ji+1,jj) + zdkt (ji,jj)  )  ) * umask(ji,jj,jk)
279                  zftv(ji,jj,jk) = (  zabe2 * zdjt(ji,jj,jk)   &
280                     &               + zcof2 * (  zdkt (ji,jj+1) + zdk1t(ji,jj)      &
281                     &                          + zdk1t(ji,jj+1) + zdkt (ji,jj)  )  ) * vmask(ji,jj,jk)                 
282               END DO
283            END DO
284            !
285            DO jj = 2 , jpjm1          !== horizontal divergence and add to pta
286               DO ji = fs_2, fs_jpim1   ! vector opt.
287                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zsign * (  zftu(ji,jj,jk) - zftu(ji-1,jj,jk)      &
288                     &                                           + zftv(ji,jj,jk) - zftv(ji,jj-1,jk)  )   &
289                     &                                        * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
290               END DO
291            END DO
292         END DO                                        !   End of slab 
293
294         !!----------------------------------------------------------------------
295         !!   III - vertical trend (full)
296         !!----------------------------------------------------------------------
297         !
298         ztfw(1,:,:) = 0._wp     ;     ztfw(jpi,:,:) = 0._wp
299         !
300         ! Vertical fluxes
301         ! ---------------
302         !                          ! Surface and bottom vertical fluxes set to zero
303         ztfw(:,:, 1 ) = 0._wp      ;      ztfw(:,:,jpk) = 0._wp
304         
305         DO jk = 2, jpkm1           ! interior (2=<jk=<jpk-1)
306            DO jj = 2, jpjm1
307               DO ji = fs_2, fs_jpim1   ! vector opt.
308                  !
309                  zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          &
310                     &                           + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk) , 1._wp  )
311                  zmskv = wmask(ji,jj,jk) / MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)          &
312                     &                           + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk) , 1._wp  )
313                     !
314                  zahu_w = (   pahu(ji  ,jj,jk-1) + pahu(ji-1,jj,jk)    &
315                     &       + pahu(ji-1,jj,jk-1) + pahu(ji  ,jj,jk)  ) * zmsku
316                  zahv_w = (   pahv(ji,jj  ,jk-1) + pahv(ji,jj-1,jk)    &
317                     &       + pahv(ji,jj-1,jk-1) + pahv(ji,jj  ,jk)  ) * zmskv
318                     !
319                  zcoef3 = - zahu_w * e2t(ji,jj) * zmsku * wslpi (ji,jj,jk)   !wslpi & j are already w-masked
320                  zcoef4 = - zahv_w * e1t(ji,jj) * zmskv * wslpj (ji,jj,jk)
321                  !
322                  ztfw(ji,jj,jk) = zcoef3 * (   zdit(ji  ,jj  ,jk-1) + zdit(ji-1,jj  ,jk)      &
323                     &                        + zdit(ji-1,jj  ,jk-1) + zdit(ji  ,jj  ,jk)  )   &
324                     &           + zcoef4 * (   zdjt(ji  ,jj  ,jk-1) + zdjt(ji  ,jj-1,jk)      &
325                     &                        + zdjt(ji  ,jj-1,jk-1) + zdjt(ji  ,jj  ,jk)  )
326               END DO
327            END DO
328         END DO
329         !                                !==  add the vertical 33 flux  ==!
330         IF( ln_traldf_lap ) THEN               ! laplacian case: eddy coef = ah_wslp2 - akz
331            DO jk = 2, jpkm1       
332               DO jj = 1, jpjm1
333                  DO ji = fs_2, fs_jpim1
334                     ztfw(ji,jj,jk) = ztfw(ji,jj,jk) + e1e2t(ji,jj) / e3w_n(ji,jj,jk) * wmask(ji,jj,jk)   &
335                        &                            * ( ah_wslp2(ji,jj,jk) - akz(ji,jj,jk) )             &
336                        &                            * ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn) )
337                  END DO
338               END DO
339            END DO
340            !
341         ELSE                                   ! bilaplacian
342            SELECT CASE( kpass )
343            CASE(  1  )                            ! 1st pass : eddy coef = ah_wslp2
344               DO jk = 2, jpkm1 
345                  DO jj = 1, jpjm1
346                     DO ji = fs_2, fs_jpim1
347                        ztfw(ji,jj,jk) = ztfw(ji,jj,jk)    &
348                           &           + ah_wslp2(ji,jj,jk) * e1e2t(ji,jj)   &
349                           &           * ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn) ) / e3w_n(ji,jj,jk) * wmask(ji,jj,jk)
350                     END DO
351                  END DO
352               END DO
353            CASE(  2  )                         ! 2nd pass : eddy flux = ah_wslp2 and akz applied on ptb  and ptbb gradients, resp.
354               DO jk = 2, jpkm1 
355                  DO jj = 1, jpjm1
356                     DO ji = fs_2, fs_jpim1
357                        ztfw(ji,jj,jk) = ztfw(ji,jj,jk) + e1e2t(ji,jj) / e3w_n(ji,jj,jk) * wmask(ji,jj,jk)                      &
358                           &                            * (  ah_wslp2(ji,jj,jk) * ( ptb (ji,jj,jk-1,jn) - ptb (ji,jj,jk,jn) )   &
359                           &                               + akz     (ji,jj,jk) * ( ptbb(ji,jj,jk-1,jn) - ptbb(ji,jj,jk,jn) )   )
360                     END DO
361                  END DO
362               END DO
363            END SELECT
364         ENDIF
365         !         
366         DO jk = 1, jpkm1                 !==  Divergence of vertical fluxes added to pta  ==!
367            DO jj = 2, jpjm1
368               DO ji = fs_2, fs_jpim1   ! vector opt.
369                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zsign * (  ztfw (ji,jj,jk) - ztfw(ji,jj,jk+1)  )   &
370                     &                                        * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
371               END DO
372            END DO
373         END DO
374         !
375         IF( ( kpass == 1 .AND. ln_traldf_lap ) .OR.  &     !==  first pass only (  laplacian)  ==!
376             ( kpass == 2 .AND. ln_traldf_blp ) ) THEN      !==  2nd   pass      (bilaplacian)  ==!
377            !
378            !                             ! "Poleward" diffusive heat or salt transports (T-S case only)
379               ! note sign is reversed to give down-gradient diffusive transports )
380            IF( l_ptr )  CALL dia_ptr_hst( jn, 'ldf', -zftv(:,:,:)  )
381            !                          ! Diffusive heat transports
382            IF( l_hst )  CALL dia_ar5_hst( jn, 'ldf', -zftu(:,:,:), -zftv(:,:,:) )
383            !
384         ENDIF                                                    !== end pass selection  ==!
385         !
386         !                                                        ! ===============
387      END DO                                                      ! end tracer loop
388      !                                                           ! ===============
389      !
390      CALL wrk_dealloc( jpi, jpj,      zdkt, zdk1t, z2d ) 
391      CALL wrk_dealloc( jpi, jpj, jpk, zdit, zdjt , zftu, zftv, ztfw  ) 
392      !
393      IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_iso')
394      !
395   END SUBROUTINE tra_ldf_iso
396
397   !!==============================================================================
398END MODULE traldf_iso
Note: See TracBrowser for help on using the repository browser.