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 3294 for trunk/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso_grif.F90 – NEMO

Ignore:
Timestamp:
2012-01-28T17:44:18+01:00 (12 years ago)
Author:
rblod
Message:

Merge of 3.4beta into the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso_grif.F90

    r2715 r3294  
    44   !! Ocean  tracers:  horizontal component of the lateral tracer mixing trend 
    55   !!====================================================================== 
    6    !! History : 3.3  ! 2010-10  (G. Nurser, C. Harris, G. Madec)   
     6   !! History : 3.3  ! 2010-10  (G. Nurser, C. Harris, G. Madec) 
    77   !!                !          Griffies operator version adapted from traldf_iso.F90 
    88   !!---------------------------------------------------------------------- 
     
    1111   !!   'key_ldfslp'               slope of the lateral diffusive direction 
    1212   !!---------------------------------------------------------------------- 
    13    !!   tra_ldf_iso_grif  : update the tracer trend with the horizontal component   
    14    !!                       of the Griffies iso-neutral laplacian operator  
     13   !!   tra_ldf_iso_grif  : update the tracer trend with the horizontal component 
     14   !!                       of the Griffies iso-neutral laplacian operator 
    1515   !!---------------------------------------------------------------------- 
    1616   USE oce             ! ocean dynamics and active tracers 
     
    2626   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    2727   USE lib_mpp         ! MPP library 
     28   USE wrk_nemo        ! Memory Allocation 
     29   USE timing          ! Timing 
     30 
    2831 
    2932   IMPLICIT NONE 
     
    3437   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   psix_eiv, psiy_eiv   !: eiv stream function (diag only) 
    3538   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   ah_wslp2             !: aeiv*w-slope^2 
    36    REAL(wp),         DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   zdkt                 !  atypic workspace 
     39   REAL(wp),         DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   zdkt3d               !: vertical tracer gradient at 2 levels 
    3740 
    3841   !! * Substitutions 
     
    4851CONTAINS 
    4952 
    50   SUBROUTINE tra_ldf_iso_grif( kt, cdtype, pgu, pgv,              & 
     53  SUBROUTINE tra_ldf_iso_grif( kt, kit000, cdtype, pgu, pgv,              & 
    5154       &                                   ptb, pta, kjpt, pahtb0 ) 
    5255      !!---------------------------------------------------------------------- 
    5356      !!                  ***  ROUTINE tra_ldf_iso_grif  *** 
    5457      !! 
    55       !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive  
    56       !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and  
     58      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive 
     59      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and 
    5760      !!      add it to the general trend of tracer equation. 
    5861      !! 
    59       !! ** Method  :   The horizontal component of the lateral diffusive trends  
     62      !! ** Method  :   The horizontal component of the lateral diffusive trends 
    6063      !!      is provided by a 2nd order operator rotated along neural or geopo- 
    6164      !!      tential surfaces to which an eddy induced advection can be added 
     
    6770      !! 
    6871      !!      2nd part :  horizontal fluxes of the lateral mixing operator 
    69       !!      ========     
     72      !!      ======== 
    7073      !!         zftu = (aht+ahtb0) e2u*e3u/e1u di[ tb ] 
    7174      !!               - aht       e2u*uslp    dk[ mi(mk(tb)) ] 
     
    8992      !! ** Action :   Update pta arrays with the before rotated diffusion 
    9093      !!---------------------------------------------------------------------- 
    91       USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released 
    9294      USE oce     , ONLY:   zftu => ua       , zftv => va            ! (ua,va) used as 3D workspace 
    93       USE wrk_nemo, ONLY:   zdit => wrk_3d_6 , zdjt => wrk_3d_7 , ztfw => wrk_3d_8   ! 3D workspace 
    94       USE wrk_nemo, ONLY:   z2d  => wrk_2d_1                                         ! 2D workspace 
    9595      ! 
    9696      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index 
     97      INTEGER                              , INTENT(in   ) ::   kit000     ! first time step index 
    9798      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
    9899      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers 
    99100      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
    100101      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb        ! before and now tracer fields 
    101       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend  
     102      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend 
    102103      REAL(wp)                             , INTENT(in   ) ::   pahtb0     ! background diffusion coef 
    103104      ! 
     
    108109      REAL(wp) ::  zmskv, zabe2, zcof2, zcoef4   !   -      - 
    109110      REAL(wp) ::  zcoef0, zbtr                  !   -      - 
    110       !REAL(wp), POINTER, DIMENSION(:,:,:) ::   zdkt           ! 2D+1 workspace 
    111111      ! 
    112112      REAL(wp) ::   zslope_skew, zslope_iso, zslope2, zbu, zbv 
     
    116116      REAL(wp) ::   zztmp              ! local scalar 
    117117#endif 
     118      REAL(wp), POINTER, DIMENSION(:,:  ) :: z2d 
     119      REAL(wp), POINTER, DIMENSION(:,:,:) :: zdit, zdjt, ztfw  
     120      REAL(wp), POINTER, DIMENSION(:,:,:) :: zw3d       ! 3D workspace 
    118121      !!---------------------------------------------------------------------- 
    119  
    120       IF( wrk_in_use(3, 6,7,8) .OR. wrk_in_use(2, 1) ) THEN 
    121          CALL ctl_stop('tra_ldf_iso_grif: requested workspace arrays unavailable.')   ;   RETURN 
    122       ENDIF 
    123       ! ARP - line below uses 'bounds re-mapping' which is only defined in 
    124       ! Fortran 2003 and up. We would be OK if code was written to use 
    125       ! zdkt(:,:,1:2) instead as then wouldn't need to re-map bounds. 
    126       ! As it is, we make zdkt a module array and allocate it in _alloc(). 
    127       !zdkt(1:jpi,1:jpj,0:1) => wrk_3d_9(:,:,1:2) 
    128  
    129       IF( kt == nit000 )  THEN 
     122      ! 
     123      IF( nn_timing == 1 )  CALL timing_start('tra_ldf_iso_grif') 
     124      ! 
     125      CALL wrk_alloc( jpi, jpj,      z2d )  
     126      CALL wrk_alloc( jpi, jpj, jpk, zdit, zdjt, ztfw  )  
     127      ! 
     128 
     129      IF( kt == kit000 .AND. .NOT.ALLOCATED(ah_wslp2) )  THEN 
    130130         IF(lwp) WRITE(numout,*) 
    131131         IF(lwp) WRITE(numout,*) 'tra_ldf_iso_grif : rotated laplacian diffusion operator on ', cdtype 
    132          IF(lwp) WRITE(numout,*) '                   WARNING: STILL UNDER TEST, NOT RECOMMENDED. USE AT YOUR OWN PERIL' 
    133132         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
    134          ALLOCATE( ah_wslp2(jpi,jpj,jpk) , zdkt(jpi,jpj,0:1), STAT=ierr ) 
     133         ALLOCATE( ah_wslp2(jpi,jpj,jpk) , zdkt3d(jpi,jpj,0:1), STAT=ierr ) 
    135134         IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
    136135         IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_iso_grif: unable to allocate arrays') 
    137136         IF( ln_traldf_gdia ) THEN 
    138             ALLOCATE( psix_eiv(jpi,jpj,jpk) , psiy_eiv(jpi,jpj,jpk) , STAT=ierr ) 
    139             IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
    140             IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_iso_grif: unable to allocate diagnostics') 
     137            IF (.NOT. ALLOCATED(psix_eiv))THEN 
     138                ALLOCATE( psix_eiv(jpi,jpj,jpk) , psiy_eiv(jpi,jpj,jpk) , STAT=ierr ) 
     139                IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
     140                IF( ierr > 0 )   CALL ctl_stop('STOP', 'tra_ldf_iso_grif: unable to allocate diagnostics') 
     141            ENDIF 
    141142         ENDIF 
    142143      ENDIF 
    143144 
    144145      !!---------------------------------------------------------------------- 
    145       !!   0 - calculate  ah_wslp2, psix_eiv, psiy_eiv  
     146      !!   0 - calculate  ah_wslp2, psix_eiv, psiy_eiv 
    146147      !!---------------------------------------------------------------------- 
    147148 
    148 !!gm Future development: consider using Ah defined at T-points and attached to the 4 t-point triads 
     149      !!gm Future development: consider using Ah defined at T-points and attached to the 4 t-point triads 
    149150 
    150151      ah_wslp2(:,:,:) = 0._wp 
     
    159160               DO jj = 1, jpjm1 
    160161                  DO ji = 1, fs_jpim1 
     162                     ze1ur = 1._wp / e1u(ji,jj) 
    161163                     ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
    162164                     zbu   = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
    163                      zah   = fsahtu(ji,jj,jk)                                  !  fsaht(ji+ip,jj,jk) 
     165                     zah   = fsahtu(ji,jj,jk)                                  ! fsaht(ji+ip,jj,jk) 
    164166                     zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
    165                      zslope2 = zslope_skew - ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * ze1ur * umask(ji,jj,jk+kp) 
     167                     ! Subtract s-coordinate slope at t-points to give slope rel to s surfaces 
     168                     ! (do this by *adding* gradient of depth) 
     169                     zslope2 = zslope_skew + ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * ze1ur * umask(ji,jj,jk+kp) 
    166170                     zslope2 = zslope2 *zslope2 
    167171                     ah_wslp2(ji+ip,jj,jk+kp) = ah_wslp2(ji+ip,jj,jk+kp)    & 
    168172                        &                     + zah * ( zbu * ze3wr / ( e1t(ji+ip,jj) * e2t(ji+ip,jj) ) ) * zslope2 
    169173                     IF( ln_traldf_gdia ) THEN 
    170                         zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew        !fsaeit(ji+ip,jj,jk)*zslope_skew 
     174                        zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew           ! fsaeit(ji+ip,jj,jk)*zslope_skew 
    171175                        psix_eiv(ji,jj,jk+kp) = psix_eiv(ji,jj,jk+kp) + 0.25_wp * zaei_slp 
    172176                     ENDIF 
     
    182186               DO jj = 1, jpjm1 
    183187                  DO ji=1,fs_jpim1 
     188                     ze2vr = 1._wp / e2v(ji,jj) 
    184189                     ze3wr = 1.0_wp / fse3w(ji,jj+jp,jk+kp) 
    185190                     zbv   = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
    186                      zah   = fsahtu(ji,jj,jk)                                       !fsaht(ji,jj+jp,jk) 
     191                     zah   = fsahtv(ji,jj,jk)                                  ! fsaht(ji,jj+jp,jk) 
    187192                     zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
    188                      zslope2 = zslope_skew - ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) * ze2vr * vmask(ji,jj,jk+kp) 
     193                     ! Subtract s-coordinate slope at t-points to give slope rel to s surfaces 
     194                     !    (do this by *adding* gradient of depth) 
     195                     zslope2 = zslope_skew + ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) * ze2vr * vmask(ji,jj,jk+kp) 
    189196                     zslope2 = zslope2 * zslope2 
    190197                     ah_wslp2(ji,jj+jp,jk+kp) = ah_wslp2(ji,jj+jp,jk+kp)   & 
    191198                        &                     + zah * ( zbv * ze3wr / ( e1t(ji,jj+jp) * e2t(ji,jj+jp) ) ) * zslope2 
    192199                     IF( ln_traldf_gdia ) THEN 
    193                         zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew     !fsaeit(ji,jj+jp,jk)*zslope_skew 
     200                        zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew           ! fsaeit(ji,jj+jp,jk)*zslope_skew 
    194201                        psiy_eiv(ji,jj,jk+kp) = psiy_eiv(ji,jj,jk+kp) + 0.25_wp * zaei_slp 
    195202                     ENDIF 
     
    200207      END DO 
    201208      ! 
     209#if defined key_iomput 
     210      IF( ln_traldf_gdia .AND. cdtype == 'TRA' ) THEN 
     211         CALL wrk_alloc( jpi , jpj , jpk  , zw3d ) 
     212         DO jk=1,jpkm1 
     213            zw3d(:,:,jk) = (psix_eiv(:,:,jk+1) - psix_eiv(:,:,jk))/fse3u(:,:,jk)  ! u_eiv = -dpsix/dz 
     214         END DO 
     215         zw3d(:,:,jpk) = 0._wp 
     216         CALL iom_put( "uoce_eiv", zw3d )    ! i-eiv current 
     217 
     218         DO jk=1,jpk-1 
     219            zw3d(:,:,jk) = (psiy_eiv(:,:,jk+1) - psiy_eiv(:,:,jk))/fse3v(:,:,jk)  ! v_eiv = -dpsiy/dz 
     220         END DO 
     221         zw3d(:,:,jpk) = 0._wp 
     222         CALL iom_put( "voce_eiv", zw3d )    ! j-eiv current 
     223 
     224         DO jk=1,jpk-1 
     225            DO jj = 2, jpjm1 
     226               DO ji = fs_2, fs_jpim1  ! vector opt. 
     227                  zw3d(ji,jj,jk) = (psiy_eiv(ji,jj,jk) - psiy_eiv(ji,jj-1,jk))/e2v(ji,jj) + & 
     228                       &    (psix_eiv(ji,jj,jk) - psix_eiv(ji-1,jj,jk))/e1u(ji,jj) ! w_eiv = dpsiy/dy + dpsiy/dx 
     229               END DO 
     230            END DO 
     231         END DO 
     232         zw3d(:,:,jpk) = 0._wp 
     233         CALL iom_put( "woce_eiv", zw3d )    ! vert. eiv current 
     234         CALL wrk_dealloc( jpi , jpj , jpk  , zw3d ) 
     235      ENDIF 
     236#endif 
    202237      !                                                          ! =========== 
    203238      DO jn = 1, kjpt                                            ! tracer loop 
     
    207242         zftu(:,:,:) = 0._wp 
    208243         zftv(:,:,:) = 0._wp 
    209          !                                                
     244         ! 
    210245         DO jk = 1, jpkm1                          !==  before lateral T & S gradients at T-level jk  ==! 
    211246            DO jj = 1, jpjm1 
     
    216251            END DO 
    217252         END DO 
    218          IF( ln_zps ) THEN                               ! partial steps: correction at the last level 
     253         IF( ln_zps.and.l_grad_zps ) THEN              ! partial steps: correction at the last level 
    219254# if defined key_vectopt_loop 
    220255            DO jj = 1, 1 
     
    224259               DO ji = 1, jpim1 
    225260# endif 
    226                   zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)           
    227                   zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)       
     261                  zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn) 
     262                  zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn) 
    228263               END DO 
    229264            END DO 
     
    237272            ! 
    238273            !                    !==  Vertical tracer gradient at level jk and jk+1 
    239             zdkt(:,:,1) = ( ptb(:,:,jk,jn) - ptb(:,:,jk+1,jn) ) * tmask(:,:,jk+1) 
     274            zdkt3d(:,:,1) = ( ptb(:,:,jk,jn) - ptb(:,:,jk+1,jn) ) * tmask(:,:,jk+1) 
    240275            ! 
    241             !                          ! surface boundary condition: zdkt(jk=1)=zdkt(jk=2) 
    242             IF( jk == 1 ) THEN   ;   zdkt(:,:,0) = zdkt(:,:,1) 
    243             ELSE                 ;   zdkt(:,:,0) = ( ptb(:,:,jk-1,jn) - ptb(:,:,jk,jn) ) * tmask(:,:,jk) 
     276            !                    ! surface boundary condition: zdkt3d(jk=0)=zdkt3d(jk=1) 
     277            IF( jk == 1 ) THEN   ;   zdkt3d(:,:,0) = zdkt3d(:,:,1) 
     278            ELSE                 ;   zdkt3d(:,:,0) = ( ptb(:,:,jk-1,jn) - ptb(:,:,jk,jn) ) * tmask(:,:,jk) 
    244279            ENDIF 
    245280 
    246             IF( .NOT. l_triad_iso ) THEN 
    247                triadi = triadi_g 
    248                triadj = triadj_g 
    249             ENDIF 
    250  
    251             DO ip = 0, 1              !==  Horizontal & vertical fluxes 
    252                DO kp = 0, 1 
    253                   DO jj = 1, jpjm1 
    254                      DO ji = 1, fs_jpim1 
    255                         ze1ur = 1._wp / e1u(ji,jj) 
    256                         zdxt = zdit(ji,jj,jk) * ze1ur 
    257                         ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
    258                         zdzt  = zdkt(ji+ip,jj,kp) * ze3wr 
    259                         zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
    260                         zslope_iso  = triadi(ji+ip,jj,jk,1-ip,kp) 
    261  
    262                         zbu = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
    263                         zah = fsahtu(ji,jj,jk)   !*umask(ji,jj,jk+kp)         !fsaht(ji+ip,jj,jk)           ===>>  ???? 
    264                         zah_slp  = zah * zslope_iso 
    265                         zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew    !fsaeit(ji+ip,jj,jk)*zslope_skew 
    266                         zftu(ji,jj,jk) = zftu(ji,jj,jk) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
    267                         ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - (zah_slp + zaei_slp) * zdxt * zbu * ze3wr 
     281 
     282            IF (ln_botmix_grif) THEN 
     283               DO ip = 0, 1              !==  Horizontal & vertical fluxes 
     284                  DO kp = 0, 1 
     285                     DO jj = 1, jpjm1 
     286                        DO ji = 1, fs_jpim1 
     287                           ze1ur = 1._wp / e1u(ji,jj) 
     288                           zdxt  = zdit(ji,jj,jk) * ze1ur 
     289                           ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
     290                           zdzt  = zdkt3d(ji+ip,jj,kp) * ze3wr 
     291                           zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
     292                           zslope_iso  = triadi(ji+ip,jj,jk,1-ip,kp) 
     293 
     294                           zbu = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
     295                           ! ln_botmix_grif is .T. don't mask zah for bottom half cells 
     296                           zah = fsahtu(ji,jj,jk)   !*umask(ji,jj,jk+kp)         !fsaht(ji+ip,jj,jk)           ===>>  ???? 
     297                           zah_slp  = zah * zslope_iso 
     298                           zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew    !fsaeit(ji+ip,jj,jk)*zslope_skew 
     299                           zftu(ji,jj,jk) = zftu(ji,jj,jk) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
     300                           ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - (zah_slp + zaei_slp) * zdxt * zbu * ze3wr 
     301                        END DO 
    268302                     END DO 
    269303                  END DO 
    270304               END DO 
    271             END DO 
    272  
    273             DO jp = 0, 1 
    274                DO kp = 0, 1 
    275                   DO jj = 1, jpjm1 
    276                      DO ji = 1, fs_jpim1 
    277                         ze2vr = 1._wp / e2v(ji,jj) 
    278                         zdyt = zdjt(ji,jj,jk) * ze2vr 
    279                         ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
    280                         zdzt = zdkt(ji,jj+jp,kp) * ze3wr 
    281                         zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
    282                         zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
    283                         zbv = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
    284                         zah = fsahtv(ji,jj,jk)        !*vmask(ji,jj,jk+kp)         !fsaht(ji,jj+jp,jk) 
    285                         zah_slp = zah * zslope_iso 
    286                         zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew    !fsaeit(ji,jj+jp,jk)*zslope_skew 
    287                         zftv(ji,jj,jk) = zftv(ji,jj,jk) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
    288                         ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - (zah_slp + zaei_slp) * zdyt * zbv * ze3wr 
     305 
     306               DO jp = 0, 1 
     307                  DO kp = 0, 1 
     308                     DO jj = 1, jpjm1 
     309                        DO ji = 1, fs_jpim1 
     310                           ze2vr = 1._wp / e2v(ji,jj) 
     311                           zdyt  = zdjt(ji,jj,jk) * ze2vr 
     312                           ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
     313                           zdzt  = zdkt3d(ji,jj+jp,kp) * ze3wr 
     314                           zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
     315                           zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
     316                           zbv = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
     317                           ! ln_botmix_grif is .T. don't mask zah for bottom half cells 
     318                           zah = fsahtv(ji,jj,jk)        !*vmask(ji,jj,jk+kp)  ! fsaht(ji,jj+jp,jk) 
     319                           zah_slp = zah * zslope_iso 
     320                           zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew        ! fsaeit(ji,jj+jp,jk)*zslope_skew 
     321                           zftv(ji,jj,jk) = zftv(ji,jj,jk) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
     322                           ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - (zah_slp + zaei_slp) * zdyt * zbv * ze3wr 
     323                        END DO 
    289324                     END DO 
    290325                  END DO 
    291326               END DO 
    292             END DO 
    293  
    294             !                        !==  divergence and add to the general trend  ==! 
     327            ELSE 
     328               DO ip = 0, 1              !==  Horizontal & vertical fluxes 
     329                  DO kp = 0, 1 
     330                     DO jj = 1, jpjm1 
     331                        DO ji = 1, fs_jpim1 
     332                           ze1ur = 1._wp / e1u(ji,jj) 
     333                           zdxt  = zdit(ji,jj,jk) * ze1ur 
     334                           ze3wr = 1._wp / fse3w(ji+ip,jj,jk+kp) 
     335                           zdzt  = zdkt3d(ji+ip,jj,kp) * ze3wr 
     336                           zslope_skew = triadi_g(ji+ip,jj,jk,1-ip,kp) 
     337                           zslope_iso  = triadi(ji+ip,jj,jk,1-ip,kp) 
     338 
     339                           zbu = 0.25_wp * e1u(ji,jj) * e2u(ji,jj) * fse3u(ji,jj,jk) 
     340                           ! ln_botmix_grif is .F. mask zah for bottom half cells 
     341                           zah = fsahtu(ji,jj,jk) * umask(ji,jj,jk+kp)         ! fsaht(ji+ip,jj,jk)   ===>>  ???? 
     342                           zah_slp  = zah * zslope_iso 
     343                           zaei_slp = fsaeiw(ji+ip,jj,jk) * zslope_skew        ! fsaeit(ji+ip,jj,jk)*zslope_skew 
     344                           zftu(ji,jj,jk) = zftu(ji,jj,jk) - ( zah * zdxt + (zah_slp - zaei_slp) * zdzt ) * zbu * ze1ur 
     345                           ztfw(ji+ip,jj,jk+kp) = ztfw(ji+ip,jj,jk+kp) - (zah_slp + zaei_slp) * zdxt * zbu * ze3wr 
     346                        END DO 
     347                     END DO 
     348                  END DO 
     349               END DO 
     350 
     351               DO jp = 0, 1 
     352                  DO kp = 0, 1 
     353                     DO jj = 1, jpjm1 
     354                        DO ji = 1, fs_jpim1 
     355                           ze2vr = 1._wp / e2v(ji,jj) 
     356                           zdyt  = zdjt(ji,jj,jk) * ze2vr 
     357                           ze3wr = 1._wp / fse3w(ji,jj+jp,jk+kp) 
     358                           zdzt  = zdkt3d(ji,jj+jp,kp) * ze3wr 
     359                           zslope_skew = triadj_g(ji,jj+jp,jk,1-jp,kp) 
     360                           zslope_iso  = triadj(ji,jj+jp,jk,1-jp,kp) 
     361                           zbv = 0.25_wp * e1v(ji,jj) * e2v(ji,jj) * fse3v(ji,jj,jk) 
     362                           ! ln_botmix_grif is .F. mask zah for bottom half cells 
     363                           zah = fsahtv(ji,jj,jk) * vmask(ji,jj,jk+kp)         ! fsaht(ji,jj+jp,jk) 
     364                           zah_slp = zah * zslope_iso 
     365                           zaei_slp = fsaeiw(ji,jj+jp,jk) * zslope_skew        ! fsaeit(ji,jj+jp,jk)*zslope_skew 
     366                           zftv(ji,jj,jk) = zftv(ji,jj,jk) - ( zah * zdyt + (zah_slp - zaei_slp) * zdzt ) * zbv * ze2vr 
     367                           ztfw(ji,jj+jp,jk+kp) = ztfw(ji,jj+jp,jk+kp) - (zah_slp + zaei_slp) * zdyt * zbv * ze3wr 
     368                        END DO 
     369                     END DO 
     370                  END DO 
     371               END DO 
     372            END IF 
     373            !                          !==  divergence and add to the general trend  ==! 
    295374            DO jj = 2 , jpjm1 
    296                DO ji = fs_2, fs_jpim1   ! vector opt. 
     375               DO ji = fs_2, fs_jpim1  ! vector opt. 
    297376                  zbtr = 1._wp / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    298377                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + zbtr * (   zftu(ji-1,jj,jk) - zftu(ji,jj,jk)   & 
     
    303382         END DO 
    304383         ! 
    305          DO jk = 1, jpkm1            !== Divergence of vertical fluxes added to the general tracer trend 
     384         DO jk = 1, jpkm1              !== Divergence of vertical fluxes added to the general tracer trend 
    306385            DO jj = 2, jpjm1 
    307                DO ji = fs_2, fs_jpim1   ! vector opt. 
     386               DO ji = fs_2, fs_jpim1  ! vector opt. 
    308387                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + (  ztfw(ji,jj,jk+1) - ztfw(ji,jj,jk)  )   & 
    309388                     &                                / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
     
    312391         END DO 
    313392         ! 
    314          !                            ! "Poleward" diffusive heat or salt transports (T-S case only) 
     393         !                             ! "Poleward" diffusive heat or salt transports (T-S case only) 
    315394         IF( cdtype == 'TRA' .AND. ln_diaptr .AND. ( MOD( kt, nn_fptr ) == 0 ) ) THEN 
    316395            IF( jn == jp_tem)   htr_ldf(:) = ptr_vj( zftv(:,:,:) )        ! 3.3  names 
     
    320399#if defined key_diaar5 
    321400         IF( cdtype == 'TRA' .AND. jn == jp_tem  ) THEN 
    322             z2d(:,:) = 0._wp  
    323             zztmp = rau0 * rcp  
     401            z2d(:,:) = 0._wp 
     402            zztmp = rau0 * rcp 
    324403            DO jk = 1, jpkm1 
    325404               DO jj = 2, jpjm1 
    326405                  DO ji = fs_2, fs_jpim1   ! vector opt. 
    327                      z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk)  
     406                     z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk) 
    328407                  END DO 
    329408               END DO 
     
    332411            CALL lbc_lnk( z2d, 'U', -1. ) 
    333412            CALL iom_put( "udiff_heattr", z2d )                  ! heat transport in i-direction 
    334             z2d(:,:) = 0._wp  
     413            z2d(:,:) = 0._wp 
    335414            DO jk = 1, jpkm1 
    336415               DO jj = 2, jpjm1 
    337416                  DO ji = fs_2, fs_jpim1   ! vector opt. 
    338                      z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk)  
     417                     z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk) 
    339418                  END DO 
    340419               END DO 
     
    342421            z2d(:,:) = zztmp * z2d(:,:) 
    343422            CALL lbc_lnk( z2d, 'V', -1. ) 
    344             CALL iom_put( "vdiff_heattr", z2d )                  !  heat transport in i-direction 
     423            CALL iom_put( "vdiff_heattr", z2d )                  !  heat transport in j-direction 
    345424         END IF 
    346425#endif 
     
    348427      END DO 
    349428      ! 
    350       IF( wrk_not_released(3, 6,7,8) .OR.   & 
    351           wrk_not_released(2, 1)       )   CALL ctl_stop('tra_ldf_iso_grif: failed to release workspace arrays') 
     429      CALL wrk_dealloc( jpi, jpj,      z2d )  
     430      CALL wrk_dealloc( jpi, jpj, jpk, zdit, zdjt, ztfw  )  
     431      ! 
     432      IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_iso_grif') 
    352433      ! 
    353434  END SUBROUTINE tra_ldf_iso_grif 
     
    357438   !!   default option :   Dummy code   NO rotation of the diffusive tensor 
    358439   !!---------------------------------------------------------------------- 
     440   REAL, PUBLIC, DIMENSION(:,:,:), ALLOCATABLE, SAVE ::   psix_eiv, psiy_eiv   !: eiv stream function (diag only) 
    359441CONTAINS 
    360    SUBROUTINE tra_ldf_iso_grif( kt, cdtype, pgu, pgv, ptb, pta, kjpt, pahtb0 )      ! Empty routine 
     442   SUBROUTINE tra_ldf_iso_grif( kt, kit000, cdtype, pgu, pgv,              & 
     443       &                                   ptb, pta, kjpt, pahtb0 ) 
    361444      CHARACTER(len=3) ::   cdtype 
     445      INTEGER          ::   kit000     ! first time step index 
    362446      REAL, DIMENSION(:,:,:) ::   pgu, pgv   ! tracer gradient at pstep levels 
    363447      REAL, DIMENSION(:,:,:,:) ::   ptb, pta 
Note: See TracChangeset for help on using the changeset viewer.