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 13497 for NEMO/trunk/src/OCE/TRA – NEMO

Ignore:
Timestamp:
2020-09-21T14:37:46+02:00 (4 years ago)
Author:
techene
Message:

re-introduce comments that have been erased by loops transformation see #2525

Location:
NEMO/trunk/src/OCE/TRA
Files:
16 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/TRA/eosbn2.F90

    r13295 r13497  
    873873      IF( ln_timing )   CALL timing_start('bn2') 
    874874      ! 
    875       DO_3D( 1, 1, 1, 1, 2, jpkm1 ) 
     875      DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! interior points only (2=< jk =< jpkm1 ); surface and bottom value set to zero one for all in istate.F90 
    876876         zrw =   ( gdepw(ji,jj,jk  ,Kmm) - gdept(ji,jj,jk,Kmm) )   & 
    877877            &  / ( gdept(ji,jj,jk-1,Kmm) - gdept(ji,jj,jk,Kmm) )  
  • NEMO/trunk/src/OCE/TRA/traadv_cen.F90

    r13457 r13497  
    112112            ztu(:,:,jpk) = 0._wp                   ! Bottom value : flux set to zero 
    113113            ztv(:,:,jpk) = 0._wp 
    114             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     114            DO_3D( 0, 0, 0, 0, 1, jpkm1 )          ! masked gradient 
    115115               ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk) 
    116116               ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk) 
     
    118118            CALL lbc_lnk_multi( 'traadv_cen', ztu, 'U', -1.0_wp , ztv, 'V', -1.0_wp )   ! Lateral boundary cond. 
    119119            ! 
    120             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     120            DO_3D( 0, 0, 0, 0, 1, jpkm1 )           ! Horizontal advective fluxes 
    121121               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! C2 interpolation of T at u- & v-points (x2) 
    122122               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm) 
     
    159159         ENDIF 
    160160         !                
    161          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     161         DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !--  Divergence of advective fluxes  --! 
    162162            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs)    & 
    163163               &             - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )    & 
     
    166166               &                * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
    167167         END_3D 
    168          !                             ! trend diagnostics 
     168         !                               ! trend diagnostics 
    169169         IF( l_trd ) THEN 
    170170            CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, zwx, pU, pt(:,:,:,jn,Kmm) ) 
  • NEMO/trunk/src/OCE/TRA/traadv_fct.F90

    r13295 r13497  
    160160            zwy(ji,jj,jk) = 0.5 * ( zfp_vj * pt(ji,jj,jk,jn,Kbb) + zfm_vj * pt(ji  ,jj+1,jk,jn,Kbb) ) 
    161161         END_3D 
    162          !                    !* upstream tracer flux in the k direction *! 
    163          DO_3D( 1, 1, 1, 1, 2, jpkm1 ) 
     162         !                               !* upstream tracer flux in the k direction *! 
     163         DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! Interior value ( multiplied by wmask) 
    164164            zfp_wk = pW(ji,jj,jk) + ABS( pW(ji,jj,jk) ) 
    165165            zfm_wk = pW(ji,jj,jk) - ABS( pW(ji,jj,jk) ) 
    166166            zwz(ji,jj,jk) = 0.5 * ( zfp_wk * pt(ji,jj,jk,jn,Kbb) + zfm_wk * pt(ji,jj,jk-1,jn,Kbb) ) * wmask(ji,jj,jk) 
    167167         END_3D 
    168          IF( ln_linssh ) THEN    ! top ocean value (only in linear free surface as zwz has been w-masked) 
    169             IF( ln_isfcav ) THEN             ! top of the ice-shelf cavities and at the ocean surface 
     168         IF( ln_linssh ) THEN               ! top ocean value (only in linear free surface as zwz has been w-masked) 
     169            IF( ln_isfcav ) THEN                        ! top of the ice-shelf cavities and at the ocean surface 
    170170               DO_2D( 1, 1, 1, 1 ) 
    171171                  zwz(ji,jj, mikt(ji,jj) ) = pW(ji,jj,mikt(ji,jj)) * pt(ji,jj,mikt(ji,jj),jn,Kbb)   ! linear free surface  
    172172               END_2D 
    173             ELSE                             ! no cavities: only at the ocean surface 
     173            ELSE                                        ! no cavities: only at the ocean surface 
    174174               DO_2D( 1, 1, 1, 1 ) 
    175175                  zwz(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kbb) 
     
    178178         ENDIF 
    179179         !                
    180          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
    181             !                             ! total intermediate advective trends 
     180         DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !* trend and after field with monotonic scheme 
     181            !                               ! total intermediate advective trends 
    182182            ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   & 
    183183               &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   & 
    184184               &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj) 
    185             !                             ! update and guess with monotonic sheme 
     185            !                               ! update and guess with monotonic sheme 
    186186            pt(ji,jj,jk,jn,Krhs) =                   pt(ji,jj,jk,jn,Krhs) +       ztra   & 
    187187               &                                  / e3t(ji,jj,jk,Kmm ) * tmask(ji,jj,jk) 
     
    194194            ! 
    195195            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp ; 
    196             DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     196            DO_3D( 0, 0, 0, 0, 2, jpkm1 )       ! Interior value ( multiplied by wmask) 
    197197               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 
    198198               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 
     
    227227            zltv(:,:,jpk) = 0._wp 
    228228            DO jk = 1, jpkm1                 ! Laplacian 
    229                DO_2D( 1, 0, 1, 0 ) 
     229               DO_2D( 1, 0, 1, 0 )                 ! 1st derivative (gradient) 
    230230                  ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk) 
    231231                  ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk) 
    232232               END_2D 
    233                DO_2D( 0, 0, 0, 0 ) 
     233               DO_2D( 0, 0, 0, 0 )                 ! 2nd derivative * 1/ 6 
    234234                  zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6 
    235235                  zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6 
     
    238238            CALL lbc_lnk_multi( 'traadv_fct', zltu, 'T', 1.0_wp , zltv, 'T', 1.0_wp )   ! Lateral boundary cond. (unchanged sgn) 
    239239            ! 
    240             DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     240            DO_3D( 1, 0, 1, 0, 1, jpkm1 )    ! Horizontal advective fluxes 
    241241               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points 
    242242               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm) 
    243                !                                                  ! C4 minus upstream advective fluxes  
     243               !                                                        ! C4 minus upstream advective fluxes  
    244244               zwx(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * ( zC2t_u + zltu(ji,jj,jk) - zltu(ji+1,jj,jk) ) - zwx(ji,jj,jk) 
    245245               zwy(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * ( zC2t_v + zltv(ji,jj,jk) - zltv(ji,jj+1,jk) ) - zwy(ji,jj,jk) 
     
    249249            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero 
    250250            ztv(:,:,jpk) = 0._wp 
    251             DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     251            DO_3D( 1, 0, 1, 0, 1, jpkm1 )    ! 1st derivative (gradient) 
    252252               ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk) 
    253253               ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk) 
     
    255255            CALL lbc_lnk_multi( 'traadv_fct', ztu, 'U', -1.0_wp , ztv, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn) 
    256256            ! 
    257             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     257            DO_3D( 0, 0, 0, 0, 1, jpkm1 )    ! Horizontal advective fluxes 
    258258               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points (x2) 
    259259               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm) 
     
    288288         !          
    289289         IF ( ll_zAimp ) THEN 
    290             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
    291                !                             ! total intermediate advective trends 
     290            DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !* trend and after field with monotonic scheme 
     291               !                                                ! total intermediate advective trends 
    292292               ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   & 
    293293                  &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   & 
     
    298298            CALL tridia_solver( zwdia, zwsup, zwinf, ztw, ztw , 0 ) 
    299299            ! 
    300             DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     300            DO_3D( 0, 0, 0, 0, 2, jpkm1 )       ! Interior value ( multiplied by wmask) 
    301301               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 
    302302               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 
     
    324324            ! 
    325325            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp 
    326             DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     326            DO_3D( 0, 0, 0, 0, 2, jpkm1 )      ! Interior value ( multiplied by wmask) 
    327327               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) ) 
    328328               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) ) 
     
    454454         pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv ) 
    455455 
    456 ! monotonic flux in the k direction, i.e. pcc 
    457 ! ------------------------------------------- 
     456      ! monotonic flux in the k direction, i.e. pcc 
     457      ! ------------------------------------------- 
    458458         za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) ) 
    459459         zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) ) 
     
    481481      !!---------------------------------------------------------------------- 
    482482       
    483       DO_3D( 1, 1, 1, 1, 3, jpkm1 ) 
     483      DO_3D( 1, 1, 1, 1, 3, jpkm1 )       !==  build the three diagonal matrix  ==! 
    484484         zwd (ji,jj,jk) = 4._wp 
    485485         zwi (ji,jj,jk) = 1._wp 
     
    495495      END_3D 
    496496      ! 
    497       jk = 2                                          ! Switch to second order centered at top 
     497      jk = 2                                    ! Switch to second order centered at top 
    498498      DO_2D( 1, 1, 1, 1 ) 
    499499         zwd (ji,jj,jk) = 1._wp 
     
    504504      ! 
    505505      !                       !==  tridiagonal solve  ==! 
    506       DO_2D( 1, 1, 1, 1 ) 
     506      DO_2D( 1, 1, 1, 1 )           ! first recurrence 
    507507         zwt(ji,jj,2) = zwd(ji,jj,2) 
    508508      END_2D 
     
    511511      END_3D 
    512512      ! 
    513       DO_2D( 1, 1, 1, 1 ) 
     513      DO_2D( 1, 1, 1, 1 )           ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1 
    514514         pt_out(ji,jj,2) = zwrm(ji,jj,2) 
    515515      END_2D 
     
    518518      END_3D 
    519519 
    520       DO_2D( 1, 1, 1, 1 ) 
     520      DO_2D( 1, 1, 1, 1 )           ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk 
    521521         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1) 
    522522      END_2D 
     
    546546      !                      !==  build the three diagonal matrix & the RHS  ==! 
    547547      ! 
    548       DO_3D( 0, 0, 0, 0, 3, jpkm1 ) 
     548      DO_3D( 0, 0, 0, 0, 3, jpkm1 )    ! interior (from jk=3 to jpk-1) 
    549549         zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal 
    550550         zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal 
     
    565565      END IF 
    566566      ! 
    567       DO_2D( 0, 0, 0, 0 ) 
     567      DO_2D( 0, 0, 0, 0 )              ! 2nd order centered at top & bottom 
    568568         ikt = mikt(ji,jj) + 1            ! w-point below the 1st  wet point 
    569569         ikb = MAX(mbkt(ji,jj), 2)        !     -   above the last wet point 
     
    582582      !                       !==  tridiagonal solver  ==! 
    583583      ! 
    584       DO_2D( 0, 0, 0, 0 ) 
     584      DO_2D( 0, 0, 0, 0 )           !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1 
    585585         zwt(ji,jj,2) = zwd(ji,jj,2) 
    586586      END_2D 
     
    589589      END_3D 
    590590      ! 
    591       DO_2D( 0, 0, 0, 0 ) 
     591      DO_2D( 0, 0, 0, 0 )           !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1 
    592592         pt_out(ji,jj,2) = zwrm(ji,jj,2) 
    593593      END_2D 
     
    596596      END_3D 
    597597 
    598       DO_2D( 0, 0, 0, 0 ) 
     598      DO_2D( 0, 0, 0, 0 )           !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk 
    599599         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1) 
    600600      END_2D 
     
    638638      kstart =  1  + klev 
    639639      ! 
    640       DO_2D( 0, 0, 0, 0 ) 
     640      DO_2D( 0, 0, 0, 0 )                         !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1 
    641641         zwt(ji,jj,kstart) = pD(ji,jj,kstart) 
    642642      END_2D 
     
    645645      END_3D 
    646646      ! 
    647       DO_2D( 0, 0, 0, 0 ) 
     647      DO_2D( 0, 0, 0, 0 )                        !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1 
    648648         pt_out(ji,jj,kstart) = pRHS(ji,jj,kstart) 
    649649      END_2D 
     
    652652      END_3D 
    653653 
    654       DO_2D( 0, 0, 0, 0 ) 
     654      DO_2D( 0, 0, 0, 0 )                       !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk 
    655655         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1) 
    656656      END_2D 
  • NEMO/trunk/src/OCE/TRA/traadv_mus.F90

    r13295 r13497  
    148148         END_3D 
    149149         ! 
    150          DO_3D( 0, 1, 0, 1, 1, jpkm1 ) 
     150         DO_3D( 0, 1, 0, 1, 1, jpkm1 )    !-- Slopes limitation 
    151151            zslpx(ji,jj,jk) = SIGN( 1.0_wp, zslpx(ji,jj,jk) ) * MIN(    ABS( zslpx(ji  ,jj,jk) ),   & 
    152152               &                                                     2.*ABS( zwx  (ji-1,jj,jk) ),   & 
     
    157157         END_3D 
    158158         ! 
    159          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     159         DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !-- MUSCL horizontal advective fluxes 
    160160            ! MUSCL fluxes 
    161161            z0u = SIGN( 0.5_wp, pU(ji,jj,jk) ) 
     
    175175         CALL lbc_lnk_multi( 'traadv_mus', zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp )   ! lateral boundary conditions   (changed sign) 
    176176         ! 
    177          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     177         DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !-- Tracer advective trend 
    178178            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )       & 
    179179            &                                     + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  ) )     & 
     
    204204               &            * (  0.25 + SIGN( 0.25_wp, zwx(ji,jj,jk) * zwx(ji,jj,jk+1) )  ) 
    205205         END_3D 
    206          DO_3D( 1, 1, 1, 1, 2, jpkm1 ) 
     206         DO_3D( 1, 1, 1, 1, 2, jpkm1 )    !-- Slopes limitation 
    207207            zslpx(ji,jj,jk) = SIGN( 1.0_wp, zslpx(ji,jj,jk) ) * MIN(    ABS( zslpx(ji,jj,jk  ) ),   & 
    208208               &                                                     2.*ABS( zwx  (ji,jj,jk+1) ),   & 
    209209               &                                                     2.*ABS( zwx  (ji,jj,jk  ) )  ) 
    210210         END_3D 
    211          DO_3D( 0, 0, 0, 0, 1, jpk-2 ) 
     211         DO_3D( 0, 0, 0, 0, 1, jpk-2 )    !-- vertical advective flux 
    212212            z0w = SIGN( 0.5_wp, pW(ji,jj,jk+1) ) 
    213213            zalpha = 0.5 + z0w 
     
    227227         ENDIF 
    228228         ! 
    229          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     229         DO_3D( 0, 0, 0, 0, 1, jpkm1 )     !-- vertical advective trend 
    230230            pt(ji,jj,jk,jn,Krhs) =  pt(ji,jj,jk,jn,Krhs) - ( zwx(ji,jj,jk) - zwx(ji,jj,jk+1) )   & 
    231231               &                                      * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
  • NEMO/trunk/src/OCE/TRA/traadv_qck.F90

    r13295 r13497  
    142142         ! 
    143143!!gm why not using a SHIFT instruction... 
    144          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     144         DO_3D( 0, 0, 0, 0, 1, jpkm1 )     !--- Computation of the ustream and downstream value of the tracer and the mask 
    145145            zfc(ji,jj,jk) = pt(ji-1,jj,jk,jn,Kbb)        ! Upstream   in the x-direction for the tracer 
    146146            zfd(ji,jj,jk) = pt(ji+1,jj,jk,jn,Kbb)        ! Downstream in the x-direction for the tracer 
     
    327327         !                                                       ! =========== 
    328328         ! 
    329          DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     329         DO_3D( 0, 0, 0, 0, 2, jpkm1 )       !* Interior point   (w-masked 2nd order centered flux) 
    330330            zwz(ji,jj,jk) = 0.5 * pW(ji,jj,jk) * ( pt(ji,jj,jk-1,jn,Kmm) + pt(ji,jj,jk,jn,Kmm) ) * wmask(ji,jj,jk) 
    331331         END_3D 
     
    340340         ENDIF 
    341341         ! 
    342          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     342         DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !==  Tracer flux divergence added to the general trend  ==! 
    343343            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( zwz(ji,jj,jk) - zwz(ji,jj,jk+1) )   & 
    344344               &                                * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
  • NEMO/trunk/src/OCE/TRA/traadv_ubs.F90

    r13295 r13497  
    124124         !                                                       ! =========== 
    125125         !                                               
    126          DO jk = 1, jpkm1        !==  horizontal laplacian of before tracer ==! 
    127             DO_2D( 1, 0, 1, 0 ) 
     126         DO jk = 1, jpkm1                !==  horizontal laplacian of before tracer ==! 
     127            DO_2D( 1, 0, 1, 0 )                   ! First derivative (masked gradient) 
    128128               zeeu = e2_e1u(ji,jj) * e3u(ji,jj,jk,Kmm) * umask(ji,jj,jk) 
    129129               zeev = e1_e2v(ji,jj) * e3v(ji,jj,jk,Kmm) * vmask(ji,jj,jk) 
     
    131131               ztv(ji,jj,jk) = zeev * ( pt(ji  ,jj+1,jk,jn,Kbb) - pt(ji,jj,jk,jn,Kbb) ) 
    132132            END_2D 
    133             DO_2D( 0, 0, 0, 0 ) 
     133            DO_2D( 0, 0, 0, 0 )                   ! Second derivative (divergence) 
    134134               zcoef = 1._wp / ( 6._wp * e3t(ji,jj,jk,Kmm) ) 
    135135               zltu(ji,jj,jk) = (  ztu(ji,jj,jk) - ztu(ji-1,jj,jk)  ) * zcoef 
     
    140140         CALL lbc_lnk( 'traadv_ubs', zltu, 'T', 1.0_wp )   ;    CALL lbc_lnk( 'traadv_ubs', zltv, 'T', 1.0_wp )   ! Lateral boundary cond. (unchanged sgn) 
    141141         !     
    142          DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
    143             zfp_ui = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) )      ! upstream transport (x2) 
     142         DO_3D( 1, 0, 1, 0, 1, jpkm1 )   !==  Horizontal advective fluxes  ==!     (UBS) 
     143            zfp_ui = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) )        ! upstream transport (x2) 
    144144            zfm_ui = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) ) 
    145145            zfp_vj = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) ) 
     
    166166         ! 
    167167         zltu(:,:,:) = pt(:,:,:,jn,Krhs) - zltu(:,:,:)    ! Horizontal advective trend used in vertical 2nd order FCT case 
    168          !                                            ! and/or in trend diagnostic (l_trd=T)  
     168         !                                                ! and/or in trend diagnostic (l_trd=T)  
    169169         !                 
    170170         IF( l_trd ) THEN                  ! trend diagnostics 
     
    187187            IF( l_trd )   zltv(:,:,:) = pt(:,:,:,jn,Krhs)          ! store pt(:,:,:,:,Krhs) if trend diag. 
    188188            ! 
    189             !                          !*  upstream advection with initial mass fluxes & intermediate update  ==! 
     189            !                               !*  upstream advection with initial mass fluxes & intermediate update  ==! 
    190190            DO_3D( 1, 1, 1, 1, 2, jpkm1 ) 
    191191               zfp_wk = pW(ji,jj,jk) + ABS( pW(ji,jj,jk) ) 
     
    193193               ztw(ji,jj,jk) = 0.5_wp * (  zfp_wk * pt(ji,jj,jk,jn,Kbb) + zfm_wk * pt(ji,jj,jk-1,jn,Kbb)  ) * wmask(ji,jj,jk) 
    194194            END_3D 
    195             IF( ln_linssh ) THEN             ! top ocean value (only in linear free surface as ztw has been w-masked) 
    196                IF( ln_isfcav ) THEN                ! top of the ice-shelf cavities and at the ocean surface 
     195            IF( ln_linssh ) THEN                ! top ocean value (only in linear free surface as ztw has been w-masked) 
     196               IF( ln_isfcav ) THEN                   ! top of the ice-shelf cavities and at the ocean surface 
    197197                  DO_2D( 1, 1, 1, 1 ) 
    198198                     ztw(ji,jj, mikt(ji,jj) ) = pW(ji,jj,mikt(ji,jj)) * pt(ji,jj,mikt(ji,jj),jn,Kbb)   ! linear free surface  
    199199                  END_2D 
    200                ELSE                                ! no cavities: only at the ocean surface 
     200               ELSE                                   ! no cavities: only at the ocean surface 
    201201                  ztw(:,:,1) = pW(:,:,1) * pt(:,:,1,jn,Kbb) 
    202202               ENDIF 
    203203            ENDIF 
    204204            ! 
    205             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     205            DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !* trend and after field with monotonic scheme 
    206206               ztak = - ( ztw(ji,jj,jk) - ztw(ji,jj,jk+1) )    & 
    207207                  &     * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
     
    230230         END SELECT 
    231231         ! 
    232          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     232         DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !  final trend with corrected fluxes 
    233233            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji,jj,jk+1) )    & 
    234234               &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
    235235         END_3D 
    236236         ! 
    237          IF( l_trd )  THEN       ! vertical advective trend diagnostics 
    238             DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     237         IF( l_trd )  THEN               ! vertical advective trend diagnostics 
     238            DO_3D( 0, 0, 0, 0, 1, jpkm1 )                 ! (compute -w.dk[ptn]= -dk[w.ptn] + ptn.dk[w]) 
    239239               zltv(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs) - zltv(ji,jj,jk)                          & 
    240240                  &           + pt(ji,jj,jk,jn,Kmm) * (  pW(ji,jj,jk) - pW(ji,jj,jk+1)  )   & 
  • NEMO/trunk/src/OCE/TRA/trabbl.F90

    r13295 r13497  
    197197         END_2D 
    198198         !                
    199          DO_2D( 0, 0, 0, 0 ) 
     199         DO_2D( 0, 0, 0, 0 )                               ! Compute the trend 
    200200            ik = mbkt(ji,jj)                            ! bottom T-level index 
    201201            pt_rhs(ji,jj,ik,jn) = pt_rhs(ji,jj,ik,jn)                                                  & 
     
    358358      IF( nn_bbl_ldf == 1 ) THEN          !   diffusive bbl   ! 
    359359         !                                !-------------------! 
    360          DO_2D( 1, 0, 1, 0 ) 
     360         DO_2D( 1, 0, 1, 0 )                   ! (criteria for non zero flux: grad(rho).grad(h) < 0 ) 
    361361            !                                                   ! i-direction 
    362362            za = zab(ji+1,jj,jp_tem) + zab(ji,jj,jp_tem)              ! 2*(alpha,beta) at u-point 
     
    388388         ! 
    389389         CASE( 1 )                                   != use of upper velocity 
    390             DO_2D( 1, 0, 1, 0 ) 
     390            DO_2D( 1, 0, 1, 0 )                              ! criteria: grad(rho).grad(h)<0  and grad(rho).grad(h)<0 
    391391               !                                                  ! i-direction 
    392392               za = zab(ji+1,jj,jp_tem) + zab(ji,jj,jp_tem)               ! 2*(alpha,beta) at u-point 
     
    417417         CASE( 2 )                                 != bbl velocity = F( delta rho ) 
    418418            zgbbl = grav * rn_gambbl 
    419             DO_2D( 1, 0, 1, 0 ) 
     419            DO_2D( 1, 0, 1, 0 )                         ! criteria: rho_up > rho_down 
    420420               !                                                  ! i-direction 
    421421               ! down-slope T-point i/k-index (deep)  &   up-slope T-point i/k-index (shelf) 
     
    509509      ! 
    510510      !                             !* vertical index of  "deep" bottom u- and v-points 
    511       DO_2D( 1, 0, 1, 0 ) 
     511      DO_2D( 1, 0, 1, 0 )                 ! (the "shelf" bottom k-indices are mbku and mbkv) 
    512512         mbku_d(ji,jj) = MAX(  mbkt(ji+1,jj  ) , mbkt(ji,jj)  )   ! >= 1 as mbkt=1 over land 
    513513         mbkv_d(ji,jj) = MAX(  mbkt(ji  ,jj+1) , mbkt(ji,jj)  ) 
     
    530530      END_2D 
    531531      ! 
    532       DO_2D( 1, 0, 1, 0 ) 
     532      DO_2D( 1, 0, 1, 0 )           !* bbl thickness at u- (v-) point; minimum of top & bottom e3u_0 (e3v_0) 
    533533         e3u_bbl_0(ji,jj) = MIN( e3u_0(ji,jj,mbkt(ji+1,jj  )), e3u_0(ji,jj,mbkt(ji,jj)) ) 
    534534         e3v_bbl_0(ji,jj) = MIN( e3v_0(ji,jj,mbkt(ji  ,jj+1)), e3v_0(ji,jj,mbkt(ji,jj)) ) 
  • NEMO/trunk/src/OCE/TRA/traldf_iso.F90

    r13295 r13497  
    205205         END_3D 
    206206         IF( ln_zps ) THEN      ! botton and surface ocean correction of the horizontal gradient 
    207             DO_2D( 1, 0, 1, 0 ) 
     207            DO_2D( 1, 0, 1, 0 )           ! bottom correction (partial bottom cell) 
    208208               zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)           
    209209               zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn) 
     
    229229            ELSE                 ;   zdkt(:,:) = ( pt(:,:,jk-1,jn) - pt(:,:,jk,jn) ) * wmask(:,:,jk) 
    230230            ENDIF 
    231             DO_2D( 1, 0, 1, 0 ) 
     231            DO_2D( 1, 0, 1, 0 )           !==  Horizontal fluxes 
    232232               zabe1 = pahu(ji,jj,jk) * e2_e1u(ji,jj) * e3u(ji,jj,jk,Kmm) 
    233233               zabe2 = pahv(ji,jj,jk) * e1_e2v(ji,jj) * e3v(ji,jj,jk,Kmm) 
     
    250250            END_2D 
    251251            ! 
    252             DO_2D( 0, 0, 0, 0 ) 
     252            DO_2D( 0, 0, 0, 0 )           !== horizontal divergence and add to pta 
    253253               pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)    & 
    254254                  &       + zsign * (  zftu(ji,jj,jk) - zftu(ji-1,jj,jk) + zftv(ji,jj,jk) - zftv(ji,jj-1,jk)  )   & 
     
    266266         ztfw(:,:, 1 ) = 0._wp      ;      ztfw(:,:,jpk) = 0._wp 
    267267          
    268          DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     268         DO_3D( 0, 0, 0, 0, 2, jpkm1 )    ! interior (2=<jk=<jpk-1) 
    269269            ! 
    270270            zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          & 
     
    311311         ENDIF 
    312312         !          
    313          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     313         DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !==  Divergence of vertical fluxes added to pta  ==! 
    314314            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn) + zsign * (  ztfw (ji,jj,jk) - ztfw(ji,jj,jk+1)  ) * r1_e1e2t(ji,jj)   & 
    315315               &                                             / e3t(ji,jj,jk,Kmm) 
  • NEMO/trunk/src/OCE/TRA/traldf_lap_blp.F90

    r13295 r13497  
    108108         !                          ! =========== !     
    109109         !                                
    110          DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     110         DO_3D( 1, 0, 1, 0, 1, jpkm1 )            !== First derivative (gradient)  ==! 
    111111            ztu(ji,jj,jk) = zaheeu(ji,jj,jk) * ( pt(ji+1,jj  ,jk,jn) - pt(ji,jj,jk,jn) ) 
    112112            ztv(ji,jj,jk) = zaheev(ji,jj,jk) * ( pt(ji  ,jj+1,jk,jn) - pt(ji,jj,jk,jn) ) 
    113113         END_3D 
    114          IF( ln_zps ) THEN                ! set gradient at bottom/top ocean level 
    115             DO_2D( 1, 0, 1, 0 ) 
     114         IF( ln_zps ) THEN                             ! set gradient at bottom/top ocean level 
     115            DO_2D( 1, 0, 1, 0 )                              ! bottom 
    116116               ztu(ji,jj,mbku(ji,jj)) = zaheeu(ji,jj,mbku(ji,jj)) * pgu(ji,jj,jn) 
    117117               ztv(ji,jj,mbkv(ji,jj)) = zaheev(ji,jj,mbkv(ji,jj)) * pgv(ji,jj,jn) 
    118118            END_2D 
    119             IF( ln_isfcav ) THEN                ! top in ocean cavities only 
     119            IF( ln_isfcav ) THEN                             ! top in ocean cavities only 
    120120               DO_2D( 1, 0, 1, 0 ) 
    121121                  IF( miku(ji,jj) > 1 )   ztu(ji,jj,miku(ji,jj)) = zaheeu(ji,jj,miku(ji,jj)) * pgui(ji,jj,jn)  
     
    125125         ENDIF 
    126126         ! 
    127          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     127         DO_3D( 0, 0, 0, 0, 1, jpkm1 )            !== Second derivative (divergence) added to the general tracer trends  ==! 
    128128            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn) + (  ztu(ji,jj,jk) - ztu(ji-1,jj,jk)     & 
    129129               &                                      +    ztv(ji,jj,jk) - ztv(ji,jj-1,jk) )   & 
  • NEMO/trunk/src/OCE/TRA/traldf_triad.F90

    r13295 r13497  
    211211         zftv(:,:,:) = 0._wp 
    212212         ! 
    213          DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     213         DO_3D( 1, 0, 1, 0, 1, jpkm1 )    !==  before lateral T & S gradients at T-level jk  ==! 
    214214            zdit(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn) - pt(ji,jj,jk,jn) ) * umask(ji,jj,jk) 
    215215            zdjt(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn) - pt(ji,jj,jk,jn) ) * vmask(ji,jj,jk) 
    216216         END_3D 
    217217         IF( ln_zps .AND. l_grad_zps ) THEN    ! partial steps: correction at top/bottom ocean level 
    218             DO_2D( 1, 0, 1, 0 ) 
     218            DO_2D( 1, 0, 1, 0 )                    ! bottom level 
    219219               zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn) 
    220220               zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn) 
     
    361361         ENDIF 
    362362         ! 
    363          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     363         DO_3D( 0, 0, 0, 0, 1, jpkm1 )      !==  Divergence of vertical fluxes added to pta  ==! 
    364364            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)    & 
    365365            &                                  + zsign * (  ztfw(ji,jj,jk+1) - ztfw(ji,jj,jk)  )   & 
  • NEMO/trunk/src/OCE/TRA/tramle.F90

    r13295 r13497  
    100100      inml_mle(:,:) = mbkt(:,:) + 1                    ! init. to number of ocean w-level (T-level + 1) 
    101101      IF ( nla10 > 0 ) THEN                            ! avoid case where first level is thicker than 10m 
    102          DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 ) 
     102         DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 )        ! from the bottom to nlb10 (10m) 
    103103            IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + rn_rho_c_mle )   inml_mle(ji,jj) = jk      ! Mixed layer 
    104104         END_3D 
     
    110110      zbm (:,:) = 0._wp 
    111111      zn2 (:,:) = 0._wp 
    112       DO_3D( 1, 1, 1, 1, 1, ikmax ) 
     112      DO_3D( 1, 1, 1, 1, 1, ikmax )                    ! MLD and mean buoyancy and N2 over the mixed layer 
    113113         zc = e3t(ji,jj,jk,Kmm) * REAL( MIN( MAX( 0, inml_mle(ji,jj)-jk ) , 1  )  )    ! zc being 0 outside the ML t-points 
    114114         zmld(ji,jj) = zmld(ji,jj) + zc 
     
    182182      zpsi_vw(:,:,:) = 0._wp 
    183183      ! 
    184       DO_3D( 1, 0, 1, 0, 2, ikmax ) 
     184      DO_3D( 1, 0, 1, 0, 2, ikmax )                ! start from 2 : surface value = 0 
    185185         zcuw = 1._wp - ( gdepw(ji+1,jj,jk,Kmm) + gdepw(ji,jj,jk,Kmm) ) * zhu(ji,jj) 
    186186         zcvw = 1._wp - ( gdepw(ji,jj+1,jk,Kmm) + gdepw(ji,jj,jk,Kmm) ) * zhv(ji,jj) 
     
    196196      !                                      !==  transport increased by the MLE induced transport ==! 
    197197      DO jk = 1, ikmax 
    198          DO_2D( 1, 0, 1, 0 ) 
     198         DO_2D( 1, 0, 1, 0 )                      ! CAUTION pu,pv must be defined at row/column i=1 / j=1 
    199199            pu(ji,jj,jk) = pu(ji,jj,jk) + ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) ) 
    200200            pv(ji,jj,jk) = pv(ji,jj,jk) + ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) ) 
     
    283283            IF( ierr /= 0 )   CALL ctl_stop( 'tra_adv_mle_init: failed to allocate arrays' ) 
    284284            z1_t2 = 1._wp / ( rn_time * rn_time ) 
    285             DO_2D( 0, 1, 0, 1 ) 
     285            DO_2D( 0, 1, 0, 1 )                      ! "coriolis+ time^-1" at u- & v-points 
    286286               zfu = ( ff_f(ji,jj) + ff_f(ji,jj-1) ) * 0.5_wp 
    287287               zfv = ( ff_f(ji,jj) + ff_f(ji-1,jj) ) * 0.5_wp 
  • NEMO/trunk/src/OCE/TRA/tranpc.F90

    r13295 r13497  
    103103         inpcc = 0 
    104104         ! 
    105          DO_2D( 0, 0, 0, 0 ) 
     105         DO_2D( 0, 0, 0, 0 )                                ! interior column only 
    106106            ! 
    107107            IF( tmask(ji,jj,2) == 1 ) THEN      ! At least 2 ocean points 
  • NEMO/trunk/src/OCE/TRA/traqsr.F90

    r13333 r13497  
    231231         END_2D 
    232232         ! 
    233          !* interior equi-partition in R-G-B depending on vertical profile of Chl 
     233         !                                    !* interior equi-partition in R-G-B depending on vertical profile of Chl 
    234234         DO_3D( 0, 0, 0, 0, 2, nksr + 1 ) 
    235235            ze3t = e3t(ji,jj,jk-1,Kmm) 
     
    246246         END_3D 
    247247         ! 
    248          DO_3D( 0, 0, 0, 0, 1, nksr ) 
     248         DO_3D( 0, 0, 0, 0, 1, nksr )          !* now qsr induced heat content 
    249249            qsr_hc(ji,jj,jk) = r1_rho0_rcp * ( ztmp3d(ji,jj,jk) - ztmp3d(ji,jj,jk+1) ) 
    250250         END_3D 
     
    256256         zz0 =        rn_abs   * r1_rho0_rcp      ! surface equi-partition in 2-bands 
    257257         zz1 = ( 1. - rn_abs ) * r1_rho0_rcp 
    258          DO_3D( 0, 0, 0, 0, 1, nksr ) 
     258         DO_3D( 0, 0, 0, 0, 1, nksr )             ! solar heat absorbed at T-point in the top 400m  
    259259            zc0 = zz0 * EXP( -gdepw(ji,jj,jk  ,Kmm)*xsi0r ) + zz1 * EXP( -gdepw(ji,jj,jk  ,Kmm)*xsi1r ) 
    260260            zc1 = zz0 * EXP( -gdepw(ji,jj,jk+1,Kmm)*xsi0r ) + zz1 * EXP( -gdepw(ji,jj,jk+1,Kmm)*xsi1r ) 
     
    264264      END SELECT 
    265265      ! 
     266      !                          !-----------------------------! 
     267      !                          !  update to the temp. trend  ! 
    266268      !                          !-----------------------------! 
    267269      DO_3D( 0, 0, 0, 0, 1, nksr ) 
  • NEMO/trunk/src/OCE/TRA/trasbc.F90

    r13295 r13497  
    129129      END_2D 
    130130      IF( ln_linssh ) THEN                !* linear free surface   
    131          DO_2D( 0, 1, 0, 0 ) 
     131         DO_2D( 0, 1, 0, 0 )                    !==>> add concentration/dilution effect due to constant volume cell 
    132132            sbc_tsc(ji,jj,jp_tem) = sbc_tsc(ji,jj,jp_tem) + r1_rho0 * emp(ji,jj) * pts(ji,jj,1,jp_tem,Kmm) 
    133133            sbc_tsc(ji,jj,jp_sal) = sbc_tsc(ji,jj,jp_sal) + r1_rho0 * emp(ji,jj) * pts(ji,jj,1,jp_sal,Kmm) 
    134          END_2D 
     134         END_2D                                 !==>> output c./d. term 
    135135         IF( iom_use('emp_x_sst') )   CALL iom_put( "emp_x_sst", emp (:,:) * pts(:,:,1,jp_tem,Kmm) ) 
    136136         IF( iom_use('emp_x_sss') )   CALL iom_put( "emp_x_sss", emp (:,:) * pts(:,:,1,jp_sal,Kmm) ) 
  • NEMO/trunk/src/OCE/TRA/trazdf.F90

    r13295 r13497  
    208208            !   used as a work space array: its value is modified. 
    209209            ! 
    210             DO_2D( 0, 0, 0, 0 ) 
     210            DO_2D( 0, 0, 0, 0 )      !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1   (increasing k) ! done one for all passive tracers (so included in the IF instruction) 
    211211               zwt(ji,jj,1) = zwd(ji,jj,1) 
    212212            END_2D 
     
    217217         ENDIF  
    218218         !          
    219          DO_2D( 0, 0, 0, 0 ) 
     219         DO_2D( 0, 0, 0, 0 )         !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1 
    220220            pt(ji,jj,1,jn,Kaa) =        e3t(ji,jj,1,Kbb) * pt(ji,jj,1,jn,Kbb)    & 
    221221               &               + p2dt * e3t(ji,jj,1,Kmm) * pt(ji,jj,1,jn,Krhs) 
     
    227227         END_3D 
    228228         ! 
    229          DO_2D( 0, 0, 0, 0 ) 
     229         DO_2D( 0, 0, 0, 0 )         !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk   (result is the after tracer) 
    230230            pt(ji,jj,jpkm1,jn,Kaa) = pt(ji,jj,jpkm1,jn,Kaa) / zwt(ji,jj,jpkm1) * tmask(ji,jj,jpkm1) 
    231231         END_2D 
  • NEMO/trunk/src/OCE/TRA/zpshde.F90

    r13295 r13497  
    167167         CALL eos( ztj, zhj, zrj )        ! at the partial step depth output in  zri, zrj  
    168168         ! 
    169          DO_2D( 1, 0, 1, 0 ) 
     169         DO_2D( 1, 0, 1, 0 )              ! Gradient of density at the last level 
    170170            iku = mbku(ji,jj) 
    171171            ikv = mbkv(ji,jj) 
     
    329329         CALL eos( ztj, zhj, zrj ) 
    330330 
    331          DO_2D( 1, 0, 1, 0 ) 
     331         DO_2D( 1, 0, 1, 0 )            ! Gradient of density at the last level 
    332332            iku = mbku(ji,jj) 
    333333            ikv = mbkv(ji,jj) 
     
    420420         CALL eos( ztj, zhj, zrj )        ! at the partial step depth output in  zri, zrj  
    421421         ! 
    422          DO_2D( 1, 0, 1, 0 ) 
     422         DO_2D( 1, 0, 1, 0 )              ! Gradient of density at the last level 
    423423            iku = miku(ji,jj)  
    424424            ikv = mikv(ji,jj)  
Note: See TracChangeset for help on using the changeset viewer.