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 13469 for NEMO/branches/2020/temporary_r4_trunk/tests/CANAL/MY_SRC/diawri.F90 – NEMO

Ignore:
Timestamp:
2020-09-15T12:49:18+02:00 (4 years ago)
Author:
smasson
Message:

r4_trunk: first change of DO loops for routines to be merged, see #2523

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/temporary_r4_trunk/tests/CANAL/MY_SRC/diawri.F90

    r13466 r13469  
    150150      CALL iom_put(  "sst", tsn(:,:,1,jp_tem) )    ! surface temperature 
    151151      IF ( iom_use("sbt") ) THEN 
    152          DO jj = 1, jpj 
    153             DO ji = 1, jpi 
    154                ikbot = mbkt(ji,jj) 
    155                z2d(ji,jj) = tsn(ji,jj,ikbot,jp_tem) 
    156             END DO 
    157          END DO 
     152         DO_2D_11_11 
     153            ikbot = mbkt(ji,jj) 
     154            z2d(ji,jj) = tsn(ji,jj,ikbot,jp_tem) 
     155         END_2D 
    158156         CALL iom_put( "sbt", z2d )                ! bottom temperature 
    159157      ENDIF 
     
    162160      CALL iom_put(  "sss", tsn(:,:,1,jp_sal) )    ! surface salinity 
    163161      IF ( iom_use("sbs") ) THEN 
    164          DO jj = 1, jpj 
    165             DO ji = 1, jpi 
    166                ikbot = mbkt(ji,jj) 
    167                z2d(ji,jj) = tsn(ji,jj,ikbot,jp_sal) 
    168             END DO 
    169          END DO 
     162         DO_2D_11_11 
     163            ikbot = mbkt(ji,jj) 
     164            z2d(ji,jj) = tsn(ji,jj,ikbot,jp_sal) 
     165         END_2D 
    170166         CALL iom_put( "sbs", z2d )                ! bottom salinity 
    171167      ENDIF 
     
    174170         zztmp = rau0 * 0.25 
    175171         z2d(:,:) = 0._wp 
    176          DO jj = 2, jpjm1 
    177             DO ji = fs_2, fs_jpim1   ! vector opt. 
    178                zztmp2 = (  ( rCdU_bot(ji+1,jj)+rCdU_bot(ji  ,jj) ) * un(ji  ,jj,mbku(ji  ,jj))  )**2   & 
    179                   &   + (  ( rCdU_bot(ji  ,jj)+rCdU_bot(ji-1,jj) ) * un(ji-1,jj,mbku(ji-1,jj))  )**2   & 
    180                   &   + (  ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj  ) ) * vn(ji,jj  ,mbkv(ji,jj  ))  )**2   & 
    181                   &   + (  ( rCdU_bot(ji,jj  )+rCdU_bot(ji,jj-1) ) * vn(ji,jj-1,mbkv(ji,jj-1))  )**2 
    182                z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1)  
    183                ! 
    184             END DO 
    185          END DO 
     172         DO_2D_00_00 
     173            zztmp2 = (  ( rCdU_bot(ji+1,jj)+rCdU_bot(ji  ,jj) ) * uu(ji  ,jj,mbku(ji  ,jj),Nii)  )**2   & 
     174               &   + (  ( rCdU_bot(ji  ,jj)+rCdU_bot(ji-1,jj) ) * uu(ji-1,jj,mbku(ji-1,jj),Nii)  )**2   & 
     175               &   + (  ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj  ) ) * vv(ji,jj  ,mbkv(ji,jj  ),Nii)  )**2   & 
     176               &   + (  ( rCdU_bot(ji,jj  )+rCdU_bot(ji,jj-1) ) * vv(ji,jj-1,mbkv(ji,jj-1),Nii)  )**2 
     177            z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1)  
     178            ! 
     179         END_2D 
    186180         CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 
    187181         CALL iom_put( "taubot", z2d )            
    188182      ENDIF 
    189183          
    190       CALL iom_put( "uoce", un(:,:,:) )            ! 3D i-current 
    191       CALL iom_put(  "ssu", un(:,:,1) )            ! surface i-current 
     184      CALL iom_put( "uoce", uu(:,:,:,Nii) )            ! 3D i-current 
     185      CALL iom_put(  "ssu", uu(:,:,1,Nii) )            ! surface i-current 
    192186      IF ( iom_use("sbu") ) THEN 
    193          DO jj = 1, jpj 
    194             DO ji = 1, jpi 
    195                ikbot = mbku(ji,jj) 
    196                z2d(ji,jj) = un(ji,jj,ikbot) 
    197             END DO 
    198          END DO 
     187         DO_2D_11_11 
     188            ikbot = mbku(ji,jj) 
     189            z2d(ji,jj) = uu(ji,jj,ikbot,Nii) 
     190         END_2D 
    199191         CALL iom_put( "sbu", z2d )                ! bottom i-current 
    200192      ENDIF 
    201193       
    202       CALL iom_put( "voce", vn(:,:,:) )            ! 3D j-current 
    203       CALL iom_put(  "ssv", vn(:,:,1) )            ! surface j-current 
     194      CALL iom_put( "voce", vv(:,:,:,Nii) )            ! 3D j-current 
     195      CALL iom_put(  "ssv", vv(:,:,1,Nii) )            ! surface j-current 
    204196      IF ( iom_use("sbv") ) THEN 
    205          DO jj = 1, jpj 
    206             DO ji = 1, jpi 
    207                ikbot = mbkv(ji,jj) 
    208                z2d(ji,jj) = vn(ji,jj,ikbot) 
    209             END DO 
    210          END DO 
     197         DO_2D_11_11 
     198            ikbot = mbkv(ji,jj) 
     199            z2d(ji,jj) = vv(ji,jj,ikbot,Nii) 
     200         END_2D 
    211201         CALL iom_put( "sbv", z2d )                ! bottom j-current 
    212202      ENDIF 
     
    217207         z2d(:,:) = rau0 * e1e2t(:,:) 
    218208         DO jk = 1, jpk 
    219             z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     209            z3d(:,:,jk) = ww(:,:,jk,Nii) * z2d(:,:) 
    220210         END DO 
    221211         CALL iom_put( "w_masstr" , z3d )   
     
    232222      IF ( iom_use("socegrad") .OR. iom_use("socegrad2") ) THEN 
    233223         z3d(:,:,jpk) = 0. 
    234          DO jk = 1, jpkm1 
    235             DO jj = 2, jpjm1                                    ! sal gradient 
    236                DO ji = fs_2, fs_jpim1   ! vector opt. 
    237                   zztmp  = tsn(ji,jj,jk,jp_sal) 
    238                   zztmpx = ( tsn(ji+1,jj,jk,jp_sal) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,jk,jp_sal) ) * r1_e1u(ji-1,jj) 
    239                   zztmpy = ( tsn(ji,jj+1,jk,jp_sal) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,jk,jp_sal) ) * r1_e2v(ji,jj-1) 
    240                   z3d(ji,jj,jk) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
    241                      &                 * umask(ji,jj,jk) * umask(ji-1,jj,jk) * vmask(ji,jj,jk) * umask(ji,jj-1,jk) 
    242                END DO 
    243             END DO 
    244          END DO 
     224         DO_3D_00_00( 1, jpkm1 ) 
     225            zztmp  = tsn(ji,jj,jk,jp_sal) 
     226            zztmpx = ( tsn(ji+1,jj,jk,jp_sal) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,jk,jp_sal) ) * r1_e1u(ji-1,jj) 
     227            zztmpy = ( tsn(ji,jj+1,jk,jp_sal) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,jk,jp_sal) ) * r1_e2v(ji,jj-1) 
     228            z3d(ji,jj,jk) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
     229               &                 * umask(ji,jj,jk) * umask(ji-1,jj,jk) * vmask(ji,jj,jk) * umask(ji,jj-1,jk) 
     230         END_3D 
    245231         CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) 
    246232         CALL iom_put( "socegrad2",  z3d )          ! square of module of sal gradient 
     
    250236          
    251237      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
    252          DO jj = 2, jpjm1                                    ! sst gradient 
    253             DO ji = fs_2, fs_jpim1   ! vector opt. 
    254                zztmp  = tsn(ji,jj,1,jp_tem) 
    255                zztmpx = ( tsn(ji+1,jj,1,jp_tem) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) * r1_e1u(ji-1,jj) 
    256                zztmpy = ( tsn(ji,jj+1,1,jp_tem) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) * r1_e2v(ji,jj-1) 
    257                z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
    258                   &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 
    259             END DO 
    260          END DO 
     238         DO_2D_00_00 
     239            zztmp  = tsn(ji,jj,1,jp_tem) 
     240            zztmpx = ( tsn(ji+1,jj,1,jp_tem) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) * r1_e1u(ji-1,jj) 
     241            zztmpy = ( tsn(ji,jj+1,1,jp_tem) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) * r1_e2v(ji,jj-1) 
     242            z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
     243               &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 
     244         END_2D 
    261245         CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 
    262246         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient 
     
    268252      IF( iom_use("heatc") ) THEN 
    269253         z2d(:,:)  = 0._wp  
    270          DO jk = 1, jpkm1 
    271             DO jj = 1, jpj 
    272                DO ji = 1, jpi 
    273                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
    274                END DO 
    275             END DO 
    276          END DO 
     254         DO_3D_11_11( 1, jpkm1 ) 
     255            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
     256         END_3D 
    277257         CALL iom_put( "heatc", rau0_rcp * z2d )   ! vertically integrated heat content (J/m2) 
    278258      ENDIF 
     
    280260      IF( iom_use("saltc") ) THEN 
    281261         z2d(:,:)  = 0._wp  
    282          DO jk = 1, jpkm1 
    283             DO jj = 1, jpj 
    284                DO ji = 1, jpi 
    285                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
    286                END DO 
    287             END DO 
    288          END DO 
     262         DO_3D_11_11( 1, jpkm1 ) 
     263            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
     264         END_3D 
    289265         CALL iom_put( "saltc", rau0 * z2d )          ! vertically integrated salt content (PSU*kg/m2) 
    290266      ENDIF 
     
    292268      IF( iom_use("salt2c") ) THEN 
    293269         z2d(:,:)  = 0._wp  
    294          DO jk = 1, jpkm1 
    295             DO jj = 1, jpj 
    296                DO ji = 1, jpi 
    297                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
    298                END DO 
    299             END DO 
    300          END DO 
     270         DO_3D_11_11( 1, jpkm1 ) 
     271            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
     272         END_3D 
    301273         CALL iom_put( "salt2c", rau0 * z2d )          ! vertically integrated squared salt content (PSU*kg/m2) 
    302274      ENDIF 
     
    304276      IF ( iom_use("eken") .OR. iom_use("eken_int") ) THEN 
    305277         z3d(:,:,jpk) = 0._wp  
    306          DO jk = 1, jpkm1 
    307             DO jj = 2, jpjm1 
    308                DO ji = 2, jpim1 
    309                   zztmpx = 0.5 * ( un(ji-1,jj  ,jk) + un(ji,jj,jk) ) 
    310                   zztmpy = 0.5 * ( vn(ji  ,jj-1,jk) + vn(ji,jj,jk) ) 
    311                   z3d(ji,jj,jk) = 0.5 * ( zztmpx*zztmpx + zztmpy*zztmpy ) 
    312                END DO 
    313             END DO 
    314          END DO 
     278         DO_3D_00_00( 1, jpkm1 ) 
     279            zztmpx = 0.5 * ( uu(ji-1,jj  ,jk,Nii) + uu(ji,jj,jk,Nii) ) 
     280            zztmpy = 0.5 * ( vv(ji  ,jj-1,jk,Nii) + vv(ji,jj,jk,Nii) ) 
     281            z3d(ji,jj,jk) = 0.5 * ( zztmpx*zztmpx + zztmpy*zztmpy ) 
     282         END_3D 
    315283         CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) 
    316284         CALL iom_put( "eken", z3d )                 ! kinetic energy 
    317285 
    318286         z2d(:,:)  = 0._wp  
    319          DO jk = 1, jpkm1 
    320             DO jj = 1, jpj 
    321                DO ji = 1, jpi 
    322                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * z3d(ji,jj,jk) * e1e2t(ji,jj) * tmask(ji,jj,jk) 
    323                END DO 
    324             END DO 
    325          END DO 
     287         DO_3D_11_11( 1, jpkm1 ) 
     288            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * z3d(ji,jj,jk) * e1e2t(ji,jj) * tmask(ji,jj,jk) 
     289         END_3D 
    326290         CALL iom_put( "eken_int", z2d )   ! vertically integrated kinetic energy 
    327291      ENDIF 
     
    332296          
    333297         z3d(:,:,jpk) = 0._wp  
    334          DO jk = 1, jpkm1 
    335             DO jj = 1, jpjm1 
    336                DO ji = 1, fs_jpim1   ! vector opt. 
    337                   z3d(ji,jj,jk) = (  e2v(ji+1,jj  ) * vn(ji+1,jj  ,jk) - e2v(ji,jj) * vn(ji,jj,jk)    & 
    338                      &             - e1u(ji  ,jj+1) * un(ji  ,jj+1,jk) + e1u(ji,jj) * un(ji,jj,jk)  ) * r1_e1e2f(ji,jj) 
    339                END DO 
    340             END DO 
    341          END DO 
     298         DO_3D_10_10( 1, jpkm1 ) 
     299            z3d(ji,jj,jk) = (  e2v(ji+1,jj  ) * vv(ji+1,jj  ,jk,Nii) - e2v(ji,jj) * vv(ji,jj,jk,Nii)    & 
     300               &             - e1u(ji  ,jj+1) * uu(ji  ,jj+1,jk,Nii) + e1u(ji,jj) * uu(ji,jj,jk,Nii)  ) * r1_e1e2f(ji,jj) 
     301         END_3D 
    342302         CALL lbc_lnk( 'diawri', z3d, 'F', 1. ) 
    343303         CALL iom_put( "relvor", z3d )                  ! relative vorticity 
    344304 
    345          DO jk = 1, jpkm1 
    346             DO jj = 1, jpj 
    347                DO ji = 1, jpi 
    348                   z3d(ji,jj,jk) = ff_f(ji,jj) + z3d(ji,jj,jk)  
    349                END DO 
    350             END DO 
    351          END DO 
     305         DO_3D_11_11( 1, jpkm1 ) 
     306            z3d(ji,jj,jk) = ff_f(ji,jj) + z3d(ji,jj,jk)  
     307         END_3D 
    352308         CALL iom_put( "absvor", z3d )                  ! absolute vorticity 
    353309 
    354          DO jk = 1, jpkm1 
    355             DO jj = 1, jpjm1 
    356                DO ji = 1, fs_jpim1   ! vector opt. 
    357                   ze3  = (  e3t_n(ji,jj+1,jk)*tmask(ji,jj+1,jk) + e3t_n(ji+1,jj+1,jk)*tmask(ji+1,jj+1,jk)   & 
    358                      &    + e3t_n(ji,jj  ,jk)*tmask(ji,jj  ,jk) + e3t_n(ji+1,jj  ,jk)*tmask(ji+1,jj  ,jk)  ) 
    359                   IF( ze3 /= 0._wp ) THEN   ;   ze3 = 4._wp / ze3 
    360                   ELSE                      ;   ze3 = 0._wp 
    361                   ENDIF 
    362                   z3d(ji,jj,jk) = ze3 * z3d(ji,jj,jk)  
    363                END DO 
    364             END DO 
    365          END DO 
     310         DO_3D_10_10( 1, jpkm1 ) 
     311            ze3  = (  e3t_n(ji,jj+1,jk)*tmask(ji,jj+1,jk) + e3t_n(ji+1,jj+1,jk)*tmask(ji+1,jj+1,jk)   & 
     312               &    + e3t_n(ji,jj  ,jk)*tmask(ji,jj  ,jk) + e3t_n(ji+1,jj  ,jk)*tmask(ji+1,jj  ,jk)  ) 
     313            IF( ze3 /= 0._wp ) THEN   ;   ze3 = 4._wp / ze3 
     314            ELSE                      ;   ze3 = 0._wp 
     315            ENDIF 
     316            z3d(ji,jj,jk) = ze3 * z3d(ji,jj,jk)  
     317         END_3D 
    366318         CALL lbc_lnk( 'diawri', z3d, 'F', 1. ) 
    367319         CALL iom_put( "potvor", z3d )                  ! potential vorticity 
     
    374326         z2d(:,:) = 0.e0 
    375327         DO jk = 1, jpkm1 
    376             z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * e3u_n(:,:,jk) * umask(:,:,jk) 
     328            z3d(:,:,jk) = rau0 * uu(:,:,jk,Nii) * e2u(:,:) * e3u_n(:,:,jk) * umask(:,:,jk) 
    377329            z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 
    378330         END DO 
     
    383335      IF( iom_use("u_heattr") ) THEN 
    384336         z2d(:,:) = 0._wp  
    385          DO jk = 1, jpkm1 
    386             DO jj = 2, jpjm1 
    387                DO ji = fs_2, fs_jpim1   ! vector opt. 
    388                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) 
    389                END DO 
    390             END DO 
    391          END DO 
     337         DO_3D_00_00( 1, jpkm1 ) 
     338            z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) 
     339         END_3D 
    392340         CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 
    393341         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction 
     
    396344      IF( iom_use("u_salttr") ) THEN 
    397345         z2d(:,:) = 0.e0  
    398          DO jk = 1, jpkm1 
    399             DO jj = 2, jpjm1 
    400                DO ji = fs_2, fs_jpim1   ! vector opt. 
    401                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) 
    402                END DO 
    403             END DO 
    404          END DO 
     346         DO_3D_00_00( 1, jpkm1 ) 
     347            z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) 
     348         END_3D 
    405349         CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 
    406350         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction 
     
    411355         z3d(:,:,jpk) = 0.e0 
    412356         DO jk = 1, jpkm1 
    413             z3d(:,:,jk) = rau0 * vn(:,:,jk) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk) 
     357            z3d(:,:,jk) = rau0 * vv(:,:,jk,Nii) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk) 
    414358         END DO 
    415359         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction 
     
    418362      IF( iom_use("v_heattr") ) THEN 
    419363         z2d(:,:) = 0.e0  
    420          DO jk = 1, jpkm1 
    421             DO jj = 2, jpjm1 
    422                DO ji = fs_2, fs_jpim1   ! vector opt. 
    423                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) 
    424                END DO 
    425             END DO 
    426          END DO 
     364         DO_3D_00_00( 1, jpkm1 ) 
     365            z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) 
     366         END_3D 
    427367         CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 
    428368         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction 
     
    431371      IF( iom_use("v_salttr") ) THEN 
    432372         z2d(:,:) = 0._wp  
    433          DO jk = 1, jpkm1 
    434             DO jj = 2, jpjm1 
    435                DO ji = fs_2, fs_jpim1   ! vector opt. 
    436                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) 
    437                END DO 
    438             END DO 
    439          END DO 
     373         DO_3D_00_00( 1, jpkm1 ) 
     374            z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) 
     375         END_3D 
    440376         CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 
    441377         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction 
     
    444380      IF( iom_use("tosmint") ) THEN 
    445381         z2d(:,:) = 0._wp 
    446          DO jk = 1, jpkm1 
    447             DO jj = 2, jpjm1 
    448                DO ji = fs_2, fs_jpim1   ! vector opt. 
    449                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
    450                END DO 
    451             END DO 
    452          END DO 
     382         DO_3D_00_00( 1, jpkm1 ) 
     383            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
     384         END_3D 
    453385         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    454386         CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature 
     
    456388      IF( iom_use("somint") ) THEN 
    457389         z2d(:,:)=0._wp 
    458          DO jk = 1, jpkm1 
    459             DO jj = 2, jpjm1 
    460                DO ji = fs_2, fs_jpim1   ! vector opt. 
    461                   z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
    462                END DO 
    463             END DO 
    464          END DO 
     390         DO_3D_00_00( 1, jpkm1 ) 
     391            z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
     392         END_3D 
    465393         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    466394         CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity 
Note: See TracChangeset for help on using the changeset viewer.