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.
Changeset 5758 for branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_triad.F90 – NEMO

Ignore:
Timestamp:
2015-09-24T08:31:40+02:00 (9 years ago)
Author:
gm
Message:

#1593: LDF-ADV, step II.1: phasing the improvements/simplifications of diffusive trend (see wiki)

File:
1 moved

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_triad.F90

    r5722 r5758  
    1 MODULE traldf_iso_grif 
     1MODULE traldf_triad 
    22   !!====================================================================== 
    3    !!                   ***  MODULE  traldf_iso_grif  *** 
     3   !!                   ***  MODULE  traldf_triad  *** 
    44   !! Ocean  tracers:  horizontal component of the lateral tracer mixing trend 
    55   !!====================================================================== 
    6    !! History : 3.3  ! 2010-10  (G. Nurser, C. Harris, G. Madec) 
    7    !!                !          Griffies operator version adapted from traldf_iso.F90 
     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 
    88   !!---------------------------------------------------------------------- 
    9 #if   defined key_ldfslp   ||   defined key_esopa 
     9 
    1010   !!---------------------------------------------------------------------- 
    11    !!   'key_ldfslp'               slope of the lateral diffusive direction 
    12    !!---------------------------------------------------------------------- 
    13    !!   tra_ldf_iso_grif  : update the tracer trend with the horizontal component 
    14    !!                       of the Griffies iso-neutral laplacian operator 
     11   !!   tra_ldf_triad : update the tracer trend with the iso-neutral laplacian triad-operator 
    1512   !!---------------------------------------------------------------------- 
    1613   USE oce             ! ocean dynamics and active tracers 
     
    1916   USE trc_oce         ! share passive tracers/Ocean variables 
    2017   USE zdf_oce         ! ocean vertical physics 
    21    USE ldftra_oce      ! ocean active tracers: lateral physics 
    22    USE ldfslp          ! iso-neutral slopes 
     18   USE ldftra          ! lateral physics: eddy diffusivity 
     19   USE ldfslp          ! lateral physics: iso-neutral slopes 
     20   USE traldf_iso      ! lateral diffusion (Madec operator)         (tra_ldf_iso routine) 
    2321   USE diaptr          ! poleward transport diagnostics 
     22   USE zpshde          ! partial step: hor. derivative     (zps_hde routine) 
     23   ! 
    2424   USE in_out_manager  ! I/O manager 
    2525   USE iom             ! I/O library 
     
    2929   USE timing          ! Timing 
    3030 
    31  
    3231   IMPLICIT NONE 
    3332   PRIVATE 
    3433 
    35    PUBLIC   tra_ldf_iso_grif   ! routine called by traldf.F90 
    36  
    37    REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   psix_eiv, psiy_eiv   !: eiv stream function (diag only) 
    38    REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   ah_wslp2             !: aeiv*w-slope^2 
    39    REAL(wp),         DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   zdkt3d               !: vertical tracer gradient at 2 levels 
     34   PUBLIC   tra_ldf_triad   ! routine called by traldf.F90 
     35 
     36   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   zdkt3d   !: vertical tracer gradient at 2 levels 
    4037 
    4138   !! * Substitutions 
    4239#  include "domzgr_substitute.h90" 
    43 #  include "ldftra_substitute.h90" 
    4440#  include "vectopt_loop_substitute.h90" 
    45 #  include "ldfeiv_substitute.h90" 
    4641   !!---------------------------------------------------------------------- 
    47    !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     42   !! NEMO/OPA 3.7 , NEMO Consortium (2015) 
    4843   !! $Id$ 
    4944   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    5146CONTAINS 
    5247 
    53   SUBROUTINE tra_ldf_iso_grif( kt, kit000, cdtype, pgu, pgv,              & 
    54        &                                   ptb, pta, kjpt, pahtb0 ) 
     48  SUBROUTINE tra_ldf_triad( kt, kit000, cdtype, pahu, pahv, pgu , pgv ,   & 
     49      &                                                     pgui, pgvi,   & 
     50      &                                         ptb , ptbb, pta , kjpt, kpass ) 
    5551      !!---------------------------------------------------------------------- 
    56       !!                  ***  ROUTINE tra_ldf_iso_grif  *** 
     52      !!                  ***  ROUTINE tra_ldf_triad  *** 
    5753      !! 
    5854      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive 
     
    6662      !!      nal or geopotential slopes computed in routine ldfslp. 
    6763      !! 
    68       !!      1st part :  masked horizontal derivative of T  ( di[ t ] ) 
    69       !!      ========    with partial cell update if ln_zps=T. 
     64      !!      see documentation for the desciption 
    7065      !! 
    71       !!      2nd part :  horizontal fluxes of the lateral mixing operator 
    72       !!      ======== 
    73       !!         zftu = (aht+ahtb0) e2u*e3u/e1u di[ tb ] 
    74       !!               - aht       e2u*uslp    dk[ mi(mk(tb)) ] 
    75       !!         zftv = (aht+ahtb0) e1v*e3v/e2v dj[ tb ] 
    76       !!               - aht       e2u*vslp    dk[ mj(mk(tb)) ] 
    77       !!      take the horizontal divergence of the fluxes: 
    78       !!         difft = 1/(e1t*e2t*e3t) {  di-1[ zftu ] +  dj-1[ zftv ]  } 
    79       !!      Add this trend to the general trend (ta,sa): 
    80       !!         ta = ta + difft 
    81       !! 
    82       !!      3rd part: vertical trends of the lateral mixing operator 
    83       !!      ========  (excluding the vertical flux proportional to dk[t] ) 
    84       !!      vertical fluxes associated with the rotated lateral mixing: 
    85       !!         zftw =-aht {  e2t*wslpi di[ mi(mk(tb)) ] 
    86       !!                     + e1t*wslpj dj[ mj(mk(tb)) ]  } 
    87       !!      take the horizontal divergence of the fluxes: 
    88       !!         difft = 1/(e1t*e2t*e3t) dk[ zftw ] 
    89       !!      Add this trend to the general trend (ta,sa): 
    90       !!         pta = pta + difft 
    91       !! 
    92       !! ** Action :   Update pta arrays with the before rotated diffusion 
     66      !! ** Action :   pta   updated with the before rotated diffusion 
     67      !!               ah_wslp2 .... 
     68      !!               akz   stabilizing vertical diffusivity coefficient (used in trazdf_imp) 
    9369      !!---------------------------------------------------------------------- 
    94       USE oce     , ONLY:   zftu => ua       , zftv => va            ! (ua,va) used as 3D workspace 
    95       ! 
    9670      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index 
    9771      INTEGER                              , INTENT(in   ) ::   kit000     ! first time step index 
    9872      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
    9973      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers 
     74      INTEGER                              , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
     75      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
    10076      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
    101       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb        ! before and now tracer fields 
     77      REAL(wp), DIMENSION(jpi,jpj,    kjpt), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
     78      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb        ! tracer (kpass=1) or laplacian of tracer (kpass=2) 
     79      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptbb       ! tracer (only used in kpass=2) 
    10280      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend 
    103       REAL(wp)                             , INTENT(in   ) ::   pahtb0     ! background diffusion coef 
    104       ! 
    105       INTEGER  ::  ji, jj, jk,jn   ! dummy loop indices 
    106       INTEGER  ::  ip,jp,kp        ! dummy loop indices 
    107       INTEGER  ::  ierr            ! temporary integer 
    108       REAL(wp) ::  zmsku, zabe1, zcof1, zcoef3   ! local scalars 
    109       REAL(wp) ::  zmskv, zabe2, zcof2, zcoef4   !   -      - 
    110       REAL(wp) ::  zcoef0, zbtr                  !   -      - 
     81      ! 
     82      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices 
     83      INTEGER  ::  ip,jp,kp         ! dummy loop indices 
     84      INTEGER  ::  ierr            ! local integer 
     85      REAL(wp) ::  zmsku, zabe1, zcof1, zcoef3          ! local scalars 
     86      REAL(wp) ::  zmskv, zabe2, zcof2, zcoef4          !   -      - 
     87      REAL(wp) ::  zcoef0, ze3w_2, zsign, z2dt, z1_2dt  !   -      - 
    11188      ! 
    11289      REAL(wp) ::   zslope_skew, zslope_iso, zslope2, zbu, zbv 
    113       REAL(wp) ::   ze1ur, zdxt, ze2vr, ze3wr, zdyt, zdzt 
     90      REAL(wp) ::   ze1ur, ze2vr, ze3wr, zdxt, zdyt, zdzt 
    11491      REAL(wp) ::   zah, zah_slp, zaei_slp 
    115       REAL(wp), POINTER, DIMENSION(:,:  ) :: z2d 
    116       REAL(wp), POINTER, DIMENSION(:,:,:) :: zdit, zdjt, ztfw  
    117       REAL(wp), POINTER, DIMENSION(:,:,:) :: zw3d       ! 3D workspace 
     92#if defined key_diaar5 
     93      REAL(wp) ::   zztmp              ! local scalar 
     94#endif 
     95      REAL(wp), POINTER, DIMENSION(:,:  ) :: z2d                                            ! 2D workspace 
     96      REAL(wp), POINTER, DIMENSION(:,:,:) :: zdit, zdjt, zftu, zftv, ztfw, zpsi_uw, zpsi_vw   ! 3D     - 
    11897      !!---------------------------------------------------------------------- 
    11998      ! 
    120       IF( nn_timing == 1 )  CALL timing_start('tra_ldf_iso_grif') 
    121       ! 
    122       CALL wrk_alloc( jpi, jpj,      z2d )  
    123       CALL wrk_alloc( jpi, jpj, jpk, zdit, zdjt, ztfw  )  
    124       ! 
    125  
    126       IF( kt == kit000 .AND. .NOT.ALLOCATED(ah_wslp2) )  THEN 
     99      IF( nn_timing == 1 )  CALL timing_start('tra_ldf_triad') 
     100      ! 
     101      CALL wrk_alloc( jpi,jpj,       z2d )  
     102      CALL wrk_alloc( jpi,jpj,jpk,   zdit, zdjt, zftu, zftv, ztfw, zpsi_uw, zpsi_vw  )  
     103      ! 
     104      IF( .NOT.ALLOCATED(zdkt3d) )  THEN 
     105         ALLOCATE( zdkt3d(jpi,jpj,0:1) , STAT=ierr ) 
     106         IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
     107         IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_triad: unable to allocate arrays') 
     108      ENDIF 
     109     ! 
     110      IF( kpass == 1 .AND. kt == kit000 )  THEN 
    127111         IF(lwp) WRITE(numout,*) 
    128          IF(lwp) WRITE(numout,*) 'tra_ldf_iso_grif : rotated laplacian diffusion operator on ', cdtype 
    129          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
    130          ALLOCATE( ah_wslp2(jpi,jpj,jpk) , zdkt3d(jpi,jpj,0:1), STAT=ierr ) 
    131          IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
    132          IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_iso_grif: unable to allocate arrays') 
    133          IF( ln_traldf_gdia ) THEN 
    134             IF (.NOT. ALLOCATED(psix_eiv))THEN 
    135                 ALLOCATE( psix_eiv(jpi,jpj,jpk) , psiy_eiv(jpi,jpj,jpk) , STAT=ierr ) 
    136                 IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
    137                 IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_iso_grif: unable to allocate diagnostics') 
    138             ENDIF 
     112         IF(lwp) WRITE(numout,*) 'tra_ldf_triad : rotated laplacian diffusion operator on ', cdtype 
     113         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~' 
     114      ENDIF 
     115      ! 
     116      !                                               ! set time step size (Euler/Leapfrog) 
     117      IF( neuler == 0 .AND. kt == kit000 ) THEN   ;   z2dt =     rdttra(1)      ! at nit000   (Euler) 
     118      ELSE                                        ;   z2dt = 2.* rdttra(1)      !             (Leapfrog) 
     119      ENDIF 
     120      z1_2dt = 1._wp / z2dt 
     121      ! 
     122      IF( kpass == 1 ) THEN   ;   zsign =  1._wp      ! bilaplacian operator require a minus sign (eddy diffusivity >0) 
     123      ELSE                    ;   zsign = -1._wp 
     124      ENDIF 
     125                   
     126      !!---------------------------------------------------------------------- 
     127      !!   0 - calculate  ah_wslp2, akz, and optionally zpsi_uw, zpsi_vw 
     128      !!---------------------------------------------------------------------- 
     129      ! 
     130      IF( kpass == 1 ) THEN         !==  first pass only  and whatever the tracer is  ==! 
     131         ! 
     132         akz     (:,:,:) = 0._wp       
     133         ah_wslp2(:,:,:) = 0._wp 
     134         IF( ln_ldfeiv_dia ) THEN 
     135            zpsi_uw(:,:,:) = 0._wp 
     136            zpsi_vw(:,:,:) = 0._wp 
    139137         ENDIF 
    140       ENDIF 
    141  
    142       !!---------------------------------------------------------------------- 
    143       !!   0 - calculate  ah_wslp2, psix_eiv, psiy_eiv 
    144       !!---------------------------------------------------------------------- 
    145  
    146       !!gm Future development: consider using Ah defined at T-points and attached to the 4 t-point triads 
    147  
    148       ah_wslp2(:,:,:) = 0._wp 
    149       IF( ln_traldf_gdia ) THEN 
    150          psix_eiv(:,:,:) = 0._wp 
    151          psiy_eiv(:,:,:) = 0._wp 
    152       ENDIF 
    153  
    154       DO ip = 0, 1 
    155          DO kp = 0, 1 
    156             DO jk = 1, jpkm1 
    157                DO jj = 1, jpjm1 
    158                   DO ji = 1, fs_jpim1 
    159                      ze1ur = 1._wp / e1u(ji,jj) 
    160                      ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
    161                      zbu   = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
    162                      zah   = fsahtu(ji,jj,jk)                                  ! fsaht(ji+ip,jj,jk) 
    163                      zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
    164                      ! Subtract s-coordinate slope at t-points to give slope rel to s surfaces 
    165                      ! (do this by *adding* gradient of depth) 
    166                      zslope2 = zslope_skew + ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * ze1ur * umask(ji,jj,jk+kp) 
    167                      zslope2 = zslope2 *zslope2 
    168                      ah_wslp2(ji+ip,jj,jk+kp) = ah_wslp2(ji+ip,jj,jk+kp)    & 
    169                         &                     + zah * ( zbu * ze3wr / ( e1t(ji+ip,jj) * e2t(ji+ip,jj) ) ) * zslope2 
    170                      IF( ln_traldf_gdia ) THEN 
    171                         zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew           ! fsaeit(ji+ip,jj,jk)*zslope_skew 
    172                         psix_eiv(ji,jj,jk+kp) = psix_eiv(ji,jj,jk+kp) + 0.25_wp * zaei_slp 
    173                      ENDIF 
     138         ! 
     139         DO ip = 0, 1                            ! i-k triads 
     140            DO kp = 0, 1 
     141               DO jk = 1, jpkm1 
     142                  DO jj = 1, jpjm1 
     143                     DO ji = 1, fs_jpim1 
     144                        ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
     145                        zbu   = e1e2u(ji,jj) * fse3u(ji,jj,jk) 
     146                        zah   = 0.25_wp * pahu(ji,jj,jk) 
     147                        zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
     148                        ! Subtract s-coordinate slope at t-points to give slope rel to s-surfaces (do this by *adding* gradient of depth) 
     149                        zslope2 = zslope_skew + ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) * r1_e1u(ji,jj) * umask(ji,jj,jk+kp) 
     150                        zslope2 = zslope2 *zslope2 
     151                        ah_wslp2(ji+ip,jj,jk+kp) = ah_wslp2(ji+ip,jj,jk+kp) + zah * zbu * ze3wr * r1_e1e2t(ji+ip,jj) * zslope2 
     152                        akz     (ji+ip,jj,jk+kp) = akz     (ji+ip,jj,jk+kp) + zah * r1_e1u(ji,jj)       & 
     153                           &                                                      * r1_e1u(ji,jj) * umask(ji,jj,jk+kp) 
     154                        ! 
     155                       IF( ln_ldfeiv_dia )   zpsi_uw(ji,jj,jk+kp) = zpsi_uw(ji,jj,jk+kp)   & 
     156                           &                                       + 0.25_wp * aeiu(ji,jj,jk) * e2u(ji,jj) * zslope_skew 
     157                     END DO 
    174158                  END DO 
    175159               END DO 
    176160            END DO 
    177161         END DO 
    178       END DO 
    179       ! 
    180       DO jp = 0, 1 
    181          DO kp = 0, 1 
    182             DO jk = 1, jpkm1 
    183                DO jj = 1, jpjm1 
    184                   DO ji=1,fs_jpim1 
    185                      ze2vr = 1._wp / e2v(ji,jj) 
    186                      ze3wr = 1.0_wp / fse3w(ji,jj+jp,jk+kp) 
    187                      zbv   = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
    188                      zah   = fsahtv(ji,jj,jk)                                  ! fsaht(ji,jj+jp,jk) 
    189                      zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
    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 = zslope_skew + ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) * ze2vr * vmask(ji,jj,jk+kp) 
    193                      zslope2 = zslope2 * zslope2 
    194                      ah_wslp2(ji,jj+jp,jk+kp) = ah_wslp2(ji,jj+jp,jk+kp)   & 
    195                         &                     + zah * ( zbv * ze3wr / ( e1t(ji,jj+jp) * e2t(ji,jj+jp) ) ) * zslope2 
    196                      IF( ln_traldf_gdia ) THEN 
    197                         zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew           ! fsaeit(ji,jj+jp,jk)*zslope_skew 
    198                         psiy_eiv(ji,jj,jk+kp) = psiy_eiv(ji,jj,jk+kp) + 0.25_wp * zaei_slp 
    199                      ENDIF 
     162         ! 
     163         DO jp = 0, 1                            ! j-k triads  
     164            DO kp = 0, 1 
     165               DO jk = 1, jpkm1 
     166                  DO jj = 1, jpjm1 
     167                     DO ji = 1, fs_jpim1 
     168                        ze3wr = 1.0_wp / fse3w(ji,jj+jp,jk+kp) 
     169                        zbv   = e1e2v(ji,jj) * fse3v(ji,jj,jk) 
     170                        zah   = 0.25_wp * pahv(ji,jj,jk) 
     171                        zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
     172                        ! Subtract s-coordinate slope at t-points to give slope rel to s surfaces 
     173                        !    (do this by *adding* gradient of depth) 
     174                        zslope2 = zslope_skew + ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk+kp) 
     175                        zslope2 = zslope2 * zslope2 
     176                        ah_wslp2(ji,jj+jp,jk+kp) = ah_wslp2(ji,jj+jp,jk+kp) + zah * zbv * ze3wr * r1_e1e2t(ji,jj+jp) * zslope2 
     177                        akz     (ji,jj+jp,jk+kp) = akz     (ji,jj+jp,jk+kp) + zah * r1_e2v(ji,jj)     & 
     178                           &                                                      * r1_e2v(ji,jj) * vmask(ji,jj,jk+kp) 
     179                        ! 
     180                        IF( ln_ldfeiv_dia )   zpsi_vw(ji,jj,jk+kp) = zpsi_vw(ji,jj,jk+kp)   & 
     181                           &                                       + 0.25 * aeiv(ji,jj,jk) * e1v(ji,jj) * zslope_skew 
     182                     END DO 
    200183                  END DO 
    201184               END DO 
    202185            END DO 
    203186         END DO 
    204       END DO 
    205       ! 
    206       IF( iom_use("uoce_eiv") .OR. iom_use("voce_eiv") .OR. iom_use("woce_eiv") )  THEN 
    207          ! 
    208          IF( ln_traldf_gdia .AND. cdtype == 'TRA' ) THEN 
    209             CALL wrk_alloc( jpi , jpj , jpk  , zw3d ) 
    210             DO jk=1,jpkm1 
    211                zw3d(:,:,jk) = (psix_eiv(:,:,jk+1) - psix_eiv(:,:,jk))/fse3u(:,:,jk)  ! u_eiv = -dpsix/dz 
    212             END DO 
    213             zw3d(:,:,jpk) = 0._wp 
    214             CALL iom_put( "uoce_eiv", zw3d )    ! i-eiv current 
    215  
    216             DO jk=1,jpk-1 
    217                zw3d(:,:,jk) = (psiy_eiv(:,:,jk+1) - psiy_eiv(:,:,jk))/fse3v(:,:,jk)  ! v_eiv = -dpsiy/dz 
    218             END DO 
    219             zw3d(:,:,jpk) = 0._wp 
    220             CALL iom_put( "voce_eiv", zw3d )    ! j-eiv current 
    221  
    222             DO jk=1,jpk-1 
    223                DO jj = 2, jpjm1 
    224                   DO ji = fs_2, fs_jpim1  ! vector opt. 
    225                      zw3d(ji,jj,jk) = (psiy_eiv(ji,jj,jk) - psiy_eiv(ji,jj-1,jk))/e2t(ji,jj) + & 
    226                           &    (psix_eiv(ji,jj,jk) - psix_eiv(ji-1,jj,jk))/e1t(ji,jj) ! w_eiv = dpsiy/dy + dpsiy/dx 
    227                   END DO 
    228                END DO 
    229             END DO 
    230             zw3d(:,:,jpk) = 0._wp 
    231             CALL iom_put( "woce_eiv", zw3d )    ! vert. eiv current 
    232             CALL wrk_dealloc( jpi , jpj , jpk  , zw3d ) 
     187         ! 
     188         IF( ln_traldf_msc ) THEN                ! stabilizing vertical diffusivity coefficient 
     189            ! 
     190            IF( ln_traldf_blp ) THEN                ! bilaplacian operator 
     191               DO jk = 2, jpkm1 
     192                  DO jj = 1, jpjm1 
     193                     DO ji = 1, fs_jpim1 
     194                        akz(ji,jj,jk) = 16._wp * ah_wslp2(ji,jj,jk)   & 
     195                           &          * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ( fse3w(ji,jj,jk) * fse3w(ji,jj,jk) )  ) 
     196                     END DO 
     197                  END DO 
     198               END DO 
     199            ELSEIF( ln_traldf_lap ) THEN              ! laplacian operator 
     200               DO jk = 2, jpkm1 
     201                  DO jj = 1, jpjm1 
     202                     DO ji = 1, fs_jpim1 
     203                        ze3w_2 = fse3w(ji,jj,jk) * fse3w(ji,jj,jk) 
     204                        zcoef0 = z2dt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
     205                        akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * z1_2dt 
     206                     END DO 
     207                  END DO 
     208               END DO 
     209           ENDIF 
     210           ! 
     211         ELSE                                    ! 33 flux set to zero with akz=ah_wslp2 ==>> computed in full implicit 
     212            akz(:,:,:) = ah_wslp2(:,:,:)       
    233213         ENDIF 
    234214         ! 
    235       ENDIF 
    236       !                                                          ! =========== 
    237       DO jn = 1, kjpt                                            ! tracer loop 
    238          !                                                       ! =========== 
     215         IF( ln_ldfeiv_dia .AND. cdtype == 'TRA' )   CALL ldf_eiv_dia( zpsi_uw, zpsi_vw ) 
     216         ! 
     217      ENDIF                                  !==  end 1st pass only  ==! 
     218      ! 
     219      !                                                           ! =========== 
     220      DO jn = 1, kjpt                                             ! tracer loop 
     221         !                                                        ! =========== 
    239222         ! Zero fluxes for each tracer 
     223!!gm  this should probably be done outside the jn loop 
    240224         ztfw(:,:,:) = 0._wp 
    241225         zftu(:,:,:) = 0._wp 
    242226         zftv(:,:,:) = 0._wp 
    243227         ! 
    244          DO jk = 1, jpkm1                          !==  before lateral T & S gradients at T-level jk  ==! 
     228         DO jk = 1, jpkm1        !==  before lateral T & S gradients at T-level jk  ==! 
    245229            DO jj = 1, jpjm1 
    246230               DO ji = 1, fs_jpim1   ! vector opt. 
     
    250234            END DO 
    251235         END DO 
    252          IF( ln_zps.and.l_grad_zps ) THEN              ! partial steps: correction at the last level 
    253             DO jj = 1, jpjm1 
    254                DO ji = 1, jpim1 
     236         IF( ln_zps .AND. l_grad_zps ) THEN    ! partial steps: correction at top/bottom ocean level 
     237            DO jj = 1, jpjm1                       ! bottom level 
     238               DO ji = 1, fs_jpim1   ! vector opt. 
    255239                  zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn) 
    256240                  zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn) 
    257241               END DO 
    258242            END DO 
     243            IF( ln_isfcav ) THEN                   ! top level (ocean cavities only) 
     244               DO jj = 1, jpjm1 
     245                  DO ji = 1, fs_jpim1   ! vector opt. 
     246                     IF( miku(ji,jj)  > 1 )   zdit(ji,jj,miku(ji,jj) ) = pgui(ji,jj,jn)  
     247                     IF( mikv(ji,jj)  > 1 )   zdjt(ji,jj,mikv(ji,jj) ) = pgvi(ji,jj,jn)  
     248                  END DO 
     249               END DO 
     250            ENDIF 
    259251         ENDIF 
    260252 
     
    272264            ELSE                 ;   zdkt3d(:,:,0) = ( ptb(:,:,jk-1,jn) - ptb(:,:,jk,jn) ) * tmask(:,:,jk) 
    273265            ENDIF 
    274  
    275  
    276             IF (ln_botmix_grif) THEN 
     266            ! 
     267            zaei_slp = 0._wp 
     268            ! 
     269            IF( ln_botmix_triad ) THEN 
    277270               DO ip = 0, 1              !==  Horizontal & vertical fluxes 
    278271                  DO kp = 0, 1 
    279272                     DO jj = 1, jpjm1 
    280273                        DO ji = 1, fs_jpim1 
    281                            ze1ur = 1._wp / e1u(ji,jj) 
     274                           ze1ur = r1_e1u(ji,jj) 
     275                           zdxt  = zdit(ji,jj,jk) * ze1ur 
     276                           ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
     277                           zdzt  = zdkt3d(ji+ip,jj,kp) * ze3wr 
     278                           zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
     279                           zslope_iso  = triadi  (ji+ip,jj,jk,1-ip,kp) 
     280 
     281                           zbu = 0.25_wp * e1e2u(ji,jj) * fse3u(ji,jj,jk) 
     282                           ! ln_botmix_triad is .T. don't mask zah for bottom half cells    !!gm ?????   ahu is masked.... 
     283                           zah = pahu(ji,jj,jk) 
     284                           zah_slp  = zah * zslope_iso 
     285                           IF( ln_ldfeiv )   zaei_slp = aeiu(ji,jj,jk) * zslope_skew 
     286                           zftu(ji   ,jj,jk   ) = zftu(ji   ,jj,jk   ) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
     287                           ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - ( zah_slp + zaei_slp) * zdxt                 * zbu * ze3wr 
     288                        END DO 
     289                     END DO 
     290                  END DO 
     291               END DO 
     292 
     293               DO jp = 0, 1 
     294                  DO kp = 0, 1 
     295                     DO jj = 1, jpjm1 
     296                        DO ji = 1, fs_jpim1 
     297                           ze2vr = r1_e2v(ji,jj) 
     298                           zdyt  = zdjt(ji,jj,jk) * ze2vr 
     299                           ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
     300                           zdzt  = zdkt3d(ji,jj+jp,kp) * ze3wr 
     301                           zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
     302                           zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
     303                           zbv = 0.25_wp * e1e2v(ji,jj) * fse3v(ji,jj,jk) 
     304                           ! ln_botmix_triad is .T. don't mask zah for bottom half cells    !!gm ?????  ahv is masked... 
     305                           zah = pahv(ji,jj,jk) 
     306                           zah_slp = zah * zslope_iso 
     307                           IF( ln_ldfeiv )   zaei_slp = aeiv(ji,jj,jk) * zslope_skew 
     308                           zftv(ji,jj   ,jk   ) = zftv(ji,jj   ,jk   ) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
     309                           ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - ( zah_slp + zaei_slp ) * zdyt                * zbv * ze3wr 
     310                        END DO 
     311                     END DO 
     312                  END DO 
     313               END DO 
     314 
     315            ELSE 
     316 
     317               DO ip = 0, 1               !==  Horizontal & vertical fluxes 
     318                  DO kp = 0, 1 
     319                     DO jj = 1, jpjm1 
     320                        DO ji = 1, fs_jpim1 
     321                           ze1ur = r1_e1u(ji,jj) 
    282322                           zdxt  = zdit(ji,jj,jk) * ze1ur 
    283323                           ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
     
    286326                           zslope_iso  = triadi(ji+ip,jj,jk,1-ip,kp) 
    287327 
    288                            zbu = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
    289                            ! ln_botmix_grif is .T. don't mask zah for bottom half cells 
    290                            zah = fsahtu(ji,jj,jk)   !*umask(ji,jj,jk+kp)         !fsaht(ji+ip,jj,jk)           ===>>  ???? 
     328                           zbu = 0.25_wp * e1e2u(ji,jj) * fse3u(ji,jj,jk) 
     329                           ! ln_botmix_triad is .F. mask zah for bottom half cells 
     330                           zah = pahu(ji,jj,jk) * umask(ji,jj,jk+kp)         ! pahu(ji+ip,jj,jk)   ===>>  ???? 
    291331                           zah_slp  = zah * zslope_iso 
    292                            zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew    !fsaeit(ji+ip,jj,jk)*zslope_skew 
    293                            zftu(ji,jj,jk) = zftu(ji,jj,jk) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
     332                           IF( ln_ldfeiv )   zaei_slp = aeiu(ji,jj,jk) * zslope_skew        ! fsaeit(ji+ip,jj,jk)*zslope_skew 
     333                           zftu(ji   ,jj,jk   ) = zftu(ji   ,jj,jk   ) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
    294334                           ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - (zah_slp + zaei_slp) * zdxt * zbu * ze3wr 
    295335                        END DO 
     
    302342                     DO jj = 1, jpjm1 
    303343                        DO ji = 1, fs_jpim1 
    304                            ze2vr = 1._wp / e2v(ji,jj) 
     344                           ze2vr = r1_e2v(ji,jj) 
    305345                           zdyt  = zdjt(ji,jj,jk) * ze2vr 
    306346                           ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
     
    308348                           zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
    309349                           zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
    310                            zbv = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
    311                            ! ln_botmix_grif is .T. don't mask zah for bottom half cells 
    312                            zah = fsahtv(ji,jj,jk)        !*vmask(ji,jj,jk+kp)  ! fsaht(ji,jj+jp,jk) 
     350                           zbv = 0.25_wp * e1e2v(ji,jj) * fse3v(ji,jj,jk) 
     351                           ! ln_botmix_triad is .F. mask zah for bottom half cells 
     352                           zah = pahv(ji,jj,jk) * vmask(ji,jj,jk+kp)         ! pahv(ji,jj+jp,jk)  ???? 
    313353                           zah_slp = zah * zslope_iso 
    314                            zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew        ! fsaeit(ji,jj+jp,jk)*zslope_skew 
     354                           IF( ln_ldfeiv )   zaei_slp = aeiv(ji,jj,jk) * zslope_skew        ! fsaeit(ji,jj+jp,jk)*zslope_skew 
    315355                           zftv(ji,jj,jk) = zftv(ji,jj,jk) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
    316356                           ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - (zah_slp + zaei_slp) * zdyt * zbv * ze3wr 
     
    319359                  END DO 
    320360               END DO 
    321             ELSE 
    322                DO ip = 0, 1              !==  Horizontal & vertical fluxes 
    323                   DO kp = 0, 1 
    324                      DO jj = 1, jpjm1 
    325                         DO ji = 1, fs_jpim1 
    326                            ze1ur = 1._wp / e1u(ji,jj) 
    327                            zdxt  = zdit(ji,jj,jk) * ze1ur 
    328                            ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
    329                            zdzt  = zdkt3d(ji+ip,jj,kp) * ze3wr 
    330                            zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
    331                            zslope_iso  = triadi(ji+ip,jj,jk,1-ip,kp) 
    332  
    333                            zbu = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
    334                            ! ln_botmix_grif is .F. mask zah for bottom half cells 
    335                            zah = fsahtu(ji,jj,jk) * umask(ji,jj,jk+kp)         ! fsaht(ji+ip,jj,jk)   ===>>  ???? 
    336                            zah_slp  = zah * zslope_iso 
    337                            zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew        ! fsaeit(ji+ip,jj,jk)*zslope_skew 
    338                            zftu(ji,jj,jk) = zftu(ji,jj,jk) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
    339                            ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - (zah_slp + zaei_slp) * zdxt * zbu * ze3wr 
    340                         END DO 
    341                      END DO 
    342                   END DO 
    343                END DO 
    344  
    345                DO jp = 0, 1 
    346                   DO kp = 0, 1 
    347                      DO jj = 1, jpjm1 
    348                         DO ji = 1, fs_jpim1 
    349                            ze2vr = 1._wp / e2v(ji,jj) 
    350                            zdyt  = zdjt(ji,jj,jk) * ze2vr 
    351                            ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
    352                            zdzt  = zdkt3d(ji,jj+jp,kp) * ze3wr 
    353                            zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
    354                            zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
    355                            zbv = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
    356                            ! ln_botmix_grif is .F. mask zah for bottom half cells 
    357                            zah = fsahtv(ji,jj,jk) * vmask(ji,jj,jk+kp)         ! fsaht(ji,jj+jp,jk) 
    358                            zah_slp = zah * zslope_iso 
    359                            zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew        ! fsaeit(ji,jj+jp,jk)*zslope_skew 
    360                            zftv(ji,jj,jk) = zftv(ji,jj,jk) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
    361                            ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - (zah_slp + zaei_slp) * zdyt * zbv * ze3wr 
    362                         END DO 
    363                      END DO 
    364                   END DO 
    365                END DO 
    366             END IF 
    367             !                          !==  divergence and add to the general trend  ==! 
     361            ENDIF 
     362            !                             !==  horizontal divergence and add to the general trend  ==! 
    368363            DO jj = 2 , jpjm1 
    369364               DO ji = fs_2, fs_jpim1  ! vector opt. 
    370                   zbtr = 1._wp / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    371                   pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zbtr * (   zftu(ji-1,jj,jk) - zftu(ji,jj,jk)   & 
    372                      &                                           + zftv(ji,jj-1,jk) - zftv(ji,jj,jk)   ) 
     365                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zsign * (  zftu(ji-1,jj,jk) - zftu(ji,jj,jk)       & 
     366                     &                                           + zftv(ji,jj-1,jk) - zftv(ji,jj,jk)   )   & 
     367                     &                                        / (  e1e2t(ji,jj) * fse3t(ji,jj,jk)  ) 
    373368               END DO 
    374369            END DO 
     
    376371         END DO 
    377372         ! 
    378          DO jk = 1, jpkm1              !== Divergence of vertical fluxes added to the general tracer trend 
     373         !                                !==  add the vertical 33 flux  ==! 
     374         IF( ln_traldf_lap ) THEN               ! laplacian case: eddy coef = ah_wslp2 - akz 
     375            DO jk = 2, jpkm1        
     376               DO jj = 1, jpjm1 
     377                  DO ji = fs_2, fs_jpim1 
     378                     ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / fse3w(ji,jj,jk) * tmask(ji,jj,jk)   & 
     379                        &                            * ( ah_wslp2(ji,jj,jk) - akz(ji,jj,jk) )             & 
     380                        &                            * ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn) ) 
     381                  END DO 
     382               END DO 
     383            END DO 
     384         ELSE                                   ! bilaplacian  
     385            SELECT CASE( kpass ) 
     386            CASE(  1  )                            ! 1st pass : eddy coef = ah_wslp2 
     387               DO jk = 2, jpkm1  
     388                  DO jj = 1, jpjm1 
     389                     DO ji = fs_2, fs_jpim1 
     390                        ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / fse3w(ji,jj,jk) * tmask(ji,jj,jk)             & 
     391                           &                            * ah_wslp2(ji,jj,jk) * ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn) ) 
     392                     END DO 
     393                  END DO 
     394               END DO  
     395            CASE(  2  )                            ! 2nd pass : eddy flux = ah_wslp2 and akz applied on ptb  and ptbb gradients, resp. 
     396               DO jk = 2, jpkm1  
     397                  DO jj = 1, jpjm1 
     398                     DO ji = fs_2, fs_jpim1 
     399                        ztfw(ji,jj,jk) = ztfw(ji,jj,jk) - e1e2t(ji,jj) / fse3w(ji,jj,jk) * tmask(ji,jj,jk)                      & 
     400                           &                            * (  ah_wslp2(ji,jj,jk) * ( ptb (ji,jj,jk-1,jn) - ptb (ji,jj,jk,jn) )   & 
     401                           &                               + akz     (ji,jj,jk) * ( ptbb(ji,jj,jk-1,jn) - ptbb(ji,jj,jk,jn) )   ) 
     402                     END DO 
     403                  END DO 
     404               END DO 
     405            END SELECT  
     406         ENDIF 
     407         ! 
     408         DO jk = 1, jpkm1                 !==  Divergence of vertical fluxes added to pta  ==! 
    379409            DO jj = 2, jpjm1 
    380410               DO ji = fs_2, fs_jpim1  ! vector opt. 
    381                   pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + (  ztfw(ji,jj,jk+1) - ztfw(ji,jj,jk)  )   & 
    382                      &                                / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
     411                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zsign * (  ztfw(ji,jj,jk+1) - ztfw(ji,jj,jk)  )   & 
     412                     &                                        / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    383413               END DO 
    384414            END DO 
    385415         END DO 
    386416         ! 
    387          !                             ! "Poleward" diffusive heat or salt transports (T-S case only) 
    388          IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
    389             IF( jn == jp_tem)   htr_ldf(:) = ptr_sj( zftv(:,:,:) )        ! 3.3  names 
    390             IF( jn == jp_sal)   str_ldf(:) = ptr_sj( zftv(:,:,:) ) 
    391          ENDIF 
    392  
    393          IF( iom_use("udiff_heattr") .OR. iom_use("vdiff_heattr") ) THEN 
    394            ! 
    395            IF( cdtype == 'TRA' .AND. jn == jp_tem  ) THEN 
    396                z2d(:,:) = 0._wp  
    397                DO jk = 1, jpkm1 
    398                   DO jj = 2, jpjm1 
    399                      DO ji = fs_2, fs_jpim1   ! vector opt. 
    400                         z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk)  
    401                      END DO 
    402                   END DO 
    403                END DO 
    404                z2d(:,:) = rau0_rcp * z2d(:,:)  
    405                CALL lbc_lnk( z2d, 'U', -1. ) 
    406                CALL iom_put( "udiff_heattr", z2d )                  ! heat transport in i-direction 
     417         IF( ( kpass == 1 .AND. ln_traldf_lap ) .OR.  &     !==  first pass only (  laplacian)  ==! 
     418             ( kpass == 2 .AND. ln_traldf_blp ) ) THEN      !==  2nd   pass      (bilaplacian)  ==! 
     419            ! 
     420            !                          ! "Poleward" diffusive heat or salt transports (T-S case only) 
     421            IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
     422               IF( jn == jp_tem)   htr_ldf(:) = ptr_sj( zftv(:,:,:) )        ! 3.3  names 
     423               IF( jn == jp_sal)   str_ldf(:) = ptr_sj( zftv(:,:,:) ) 
     424            ENDIF 
     425            ! 
     426            IF( iom_use("udiff_heattr") .OR. iom_use("vdiff_heattr") ) THEN 
     427              ! 
     428              IF( cdtype == 'TRA' .AND. jn == jp_tem  ) THEN 
     429                  z2d(:,:) = zftu(ji,jj,1)  
     430                  DO jk = 2, jpkm1 
     431                     DO jj = 2, jpjm1 
     432                        DO ji = fs_2, fs_jpim1   ! vector opt. 
     433                           z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk)  
     434                        END DO 
     435                     END DO 
     436                  END DO 
     437                  z2d(:,:) = rau0_rcp * z2d(:,:)  
     438                  CALL lbc_lnk( z2d, 'U', -1. ) 
     439                  CALL iom_put( "udiff_heattr", z2d )                  ! heat i-transport 
     440                  ! 
     441                  z2d(:,:) = zftv(ji,jj,1)  
     442                  DO jk = 2, jpkm1 
     443                     DO jj = 2, jpjm1 
     444                        DO ji = fs_2, fs_jpim1   ! vector opt. 
     445                           z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk)  
     446                        END DO 
     447                     END DO 
     448                  END DO 
     449                  z2d(:,:) = rau0_rcp * z2d(:,:)      
     450                  CALL lbc_lnk( z2d, 'V', -1. ) 
     451                  CALL iom_put( "vdiff_heattr", z2d )                  !  heat j-transport 
     452               ENDIF 
    407453               ! 
    408                z2d(:,:) = 0._wp  
    409                DO jk = 1, jpkm1 
    410                   DO jj = 2, jpjm1 
    411                      DO ji = fs_2, fs_jpim1   ! vector opt. 
    412                         z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk)  
    413                      END DO 
    414                   END DO 
    415                END DO 
    416                z2d(:,:) = rau0_rcp * z2d(:,:)      
    417                CALL lbc_lnk( z2d, 'V', -1. ) 
    418                CALL iom_put( "vdiff_heattr", z2d )                  !  heat transport in i-direction 
    419             END IF 
    420             ! 
    421          ENDIF 
    422          ! 
    423       END DO 
    424       ! 
    425       CALL wrk_dealloc( jpi, jpj,      z2d )  
    426       CALL wrk_dealloc( jpi, jpj, jpk, zdit, zdjt, ztfw  )  
    427       ! 
    428       IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_iso_grif') 
    429       ! 
    430   END SUBROUTINE tra_ldf_iso_grif 
    431  
    432 #else 
    433    !!---------------------------------------------------------------------- 
    434    !!   default option :   Dummy code   NO rotation of the diffusive tensor 
    435    !!---------------------------------------------------------------------- 
    436    REAL, PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   psix_eiv, psiy_eiv   !: eiv stream function (diag only) 
    437 CONTAINS 
    438    SUBROUTINE tra_ldf_iso_grif( kt, kit000, cdtype, pgu, pgv,              & 
    439        &                                   ptb, pta, kjpt, pahtb0 ) 
    440       CHARACTER(len=3) ::   cdtype 
    441       INTEGER          ::   kit000     ! first time step index 
    442       REAL, DIMENSION(:,:,:) ::   pgu, pgv   ! tracer gradient at pstep levels 
    443       REAL, DIMENSION(:,:,:,:) ::   ptb, pta 
    444       WRITE(*,*) 'tra_ldf_iso_grif: You should not have seen this print! error?', kt, cdtype,    & 
    445          &                  pgu(1,1,1), pgv(1,1,1), ptb(1,1,1,1), pta(1,1,1,1), kjpt, pahtb0 
    446    END SUBROUTINE tra_ldf_iso_grif 
    447 #endif 
     454            ENDIF 
     455            ! 
     456         ENDIF                                                    !== end pass selection  ==! 
     457         ! 
     458         !                                                        ! =============== 
     459      END DO                                                      ! end tracer loop 
     460      !                                                           ! =============== 
     461      ! 
     462      CALL wrk_dealloc( jpi,jpj,       z2d )  
     463      CALL wrk_dealloc( jpi,jpj,jpk,   zdit, zdjt, zftu, zftv, ztfw, zpsi_uw, zpsi_vw  )  
     464      ! 
     465      IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_triad') 
     466      ! 
     467   END SUBROUTINE tra_ldf_triad 
    448468 
    449469   !!============================================================================== 
    450 END MODULE traldf_iso_grif 
     470END MODULE traldf_triad 
Note: See TracChangeset for help on using the changeset viewer.