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 13982 for NEMO/trunk/src/OCE/LDF/ldftra.F90 – NEMO

Ignore:
Timestamp:
2020-12-02T11:57:05+01:00 (3 years ago)
Author:
smasson
Message:

trunk: merge dev_r13923_Tiling_Cleanup_MPI3_LoopFusion into the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/LDF/ldftra.F90

    r13558 r13982  
    427427         zaht_min = 0.2_wp * aht0                                       ! minimum value for aht 
    428428         zDaht    = aht0 - zaht_min                                       
    429          DO_2D( 1, 1, 1, 1 ) 
     429         ! NOTE: [tiling-comms-merge] Change needed to preserve results with respect to the trunk 
     430         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 
    430431            !!gm CAUTION : here we assume lat/lon grid in 20deg N/S band (like all ORCA cfg) 
    431432            !!     ==>>>   The Coriolis value is identical for t- & u_points, and for v- and f-points 
     
    725726      !! ** Action  : pu, pv increased by the eiv transport 
    726727      !!---------------------------------------------------------------------- 
    727       INTEGER                         , INTENT(in   ) ::   kt        ! ocean time-step index 
    728       INTEGER                         , INTENT(in   ) ::   kit000    ! first time step index 
    729       INTEGER                         , INTENT(in   ) ::   Kmm, Krhs ! ocean time level indices 
    730       CHARACTER(len=3)                , INTENT(in   ) ::   cdtype    ! =TRA or TRC (tracer indicator) 
    731       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pu      ! in : 3 ocean transport components   [m3/s] 
    732       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pv      ! out: 3 ocean transport components   [m3/s] 
    733       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pw      ! increased by the eiv                [m3/s] 
     728      INTEGER                     , INTENT(in   ) ::   kt        ! ocean time-step index 
     729      INTEGER                     , INTENT(in   ) ::   kit000    ! first time step index 
     730      INTEGER                     , INTENT(in   ) ::   Kmm, Krhs ! ocean time level indices 
     731      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype    ! =TRA or TRC (tracer indicator) 
     732      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pu        ! in : 3 ocean transport components   [m3/s] 
     733      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pv        ! out: 3 ocean transport components   [m3/s] 
     734      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   pw        ! increased by the eiv                [m3/s] 
    734735      !! 
    735736      INTEGER  ::   ji, jj, jk                 ! dummy loop indices 
    736737      REAL(wp) ::   zuwk, zuwk1, zuwi, zuwi1   ! local scalars 
    737738      REAL(wp) ::   zvwk, zvwk1, zvwj, zvwj1   !   -      - 
    738       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpsi_uw, zpsi_vw 
    739       !!---------------------------------------------------------------------- 
    740       ! 
    741       IF( kt == kit000 )  THEN 
    742          IF(lwp) WRITE(numout,*) 
    743          IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :' 
    744          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component' 
     739      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zpsi_uw, zpsi_vw 
     740      !!---------------------------------------------------------------------- 
     741      ! 
     742      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     743         IF( kt == kit000 )  THEN 
     744            IF(lwp) WRITE(numout,*) 
     745            IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :' 
     746            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component' 
     747         ENDIF 
    745748      ENDIF 
    746749 
     
    781784      !! 
    782785      !!---------------------------------------------------------------------- 
    783       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s] 
    784       INTEGER                         , INTENT(in   ) ::   Kmm   ! ocean time level indices 
     786      REAL(wp), DIMENSION(A2D(nn_hls),jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s] 
     787      INTEGER                     , INTENT(in   ) ::   Kmm   ! ocean time level indices 
    785788      ! 
    786789      INTEGER  ::   ji, jj, jk    ! dummy loop indices 
    787790      REAL(wp) ::   zztmp   ! local scalar 
    788       REAL(wp), DIMENSION(jpi,jpj)     ::   zw2d   ! 2D workspace 
    789       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zw3d   ! 3D workspace 
     791      REAL(wp), DIMENSION(A2D(nn_hls))     ::   zw2d   ! 2D workspace 
     792      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zw3d   ! 3D workspace 
    790793      !!---------------------------------------------------------------------- 
    791794      ! 
     
    793796!!gm     to be redesigned....    
    794797      !                                                  !==  eiv stream function: output  ==! 
    795       CALL lbc_lnk_multi( 'ldftra', psi_uw, 'U', -1.0_wp , psi_vw, 'V', -1.0_wp ) 
    796       ! 
    797798!!gm      CALL iom_put( "psi_eiv_uw", psi_uw )                 ! output 
    798799!!gm      CALL iom_put( "psi_eiv_vw", psi_vw ) 
     
    802803      zw3d(:,:,jpk) = 0._wp                                    ! bottom value always 0 
    803804      ! 
    804       DO jk = 1, jpkm1                                         ! e2u e3u u_eiv = -dk[psi_uw] 
    805          zw3d(:,:,jk) = ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) / ( e2u(:,:) * e3u(:,:,jk,Kmm) ) 
    806       END DO 
     805      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                                  ! e2u e3u u_eiv = -dk[psi_uw] 
     806         zw3d(ji,jj,jk) = ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm) ) 
     807      END_3D 
    807808      CALL iom_put( "uoce_eiv", zw3d ) 
    808809      ! 
    809       DO jk = 1, jpkm1                                         ! e1v e3v v_eiv = -dk[psi_vw] 
    810          zw3d(:,:,jk) = ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) / ( e1v(:,:) * e3v(:,:,jk,Kmm) ) 
    811       END DO 
     810      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                                  ! e1v e3v v_eiv = -dk[psi_vw] 
     811         zw3d(ji,jj,jk) = ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm) ) 
     812      END_3D 
    812813      CALL iom_put( "voce_eiv", zw3d ) 
    813814      ! 
     
    816817            &              + psi_uw(ji,jj,jk) - psi_uw(ji-1,jj  ,jk)  ) / e1e2t(ji,jj) 
    817818      END_3D 
    818       CALL lbc_lnk( 'ldftra', zw3d, 'T', 1.0_wp )      ! lateral boundary condition 
    819819      CALL iom_put( "woce_eiv", zw3d ) 
    820820      ! 
    821821      IF( iom_use('weiv_masstr') ) THEN   ! vertical mass transport & its square value 
    822          zw2d(:,:) = rho0 * e1e2t(:,:) 
     822         DO_2D( 0, 0, 0, 0 ) 
     823            zw2d(ji,jj) = rho0 * e1e2t(ji,jj) 
     824         END_2D 
    823825         DO jk = 1, jpk 
    824826            zw3d(:,:,jk) = zw3d(:,:,jk) * zw2d(:,:) 
     
    844846           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 
    845847        END_3D 
    846         CALL lbc_lnk( 'ldftra', zw2d, 'U', -1.0_wp ) 
    847         CALL lbc_lnk( 'ldftra', zw3d, 'U', -1.0_wp ) 
    848848        CALL iom_put( "ueiv_heattr"  , zztmp * zw2d )                  ! heat transport in i-direction 
    849849        CALL iom_put( "ueiv_heattr3d", zztmp * zw3d )                  ! heat transport in i-direction 
     
    865865         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 
    866866      END_3D 
    867       CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp ) 
    868       CALL iom_put( "veiv_heattr", zztmp * zw2d )                  !  heat transport in j-direction 
    869       CALL iom_put( "veiv_heattr", zztmp * zw3d )                  !  heat transport in j-direction 
     867      CALL iom_put( "veiv_heattr"  , zztmp * zw2d )                  !  heat transport in j-direction 
     868      CALL iom_put( "veiv_heattr3d", zztmp * zw3d )                  !  heat transport in j-direction 
    870869      ! 
    871870      IF( iom_use( 'sophteiv' ) )   CALL dia_ptr_hst( jp_tem, 'eiv', 0.5 * zw3d ) 
     
    880879           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 
    881880        END_3D 
    882         CALL lbc_lnk( 'ldftra', zw2d, 'U', -1.0_wp ) 
    883         CALL lbc_lnk( 'ldftra', zw3d, 'U', -1.0_wp ) 
    884881        CALL iom_put( "ueiv_salttr", zztmp * zw2d )                  ! salt transport in i-direction 
    885882        CALL iom_put( "ueiv_salttr3d", zztmp * zw3d )                ! salt transport in i-direction 
     
    892889         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 
    893890      END_3D 
    894       CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp ) 
    895       CALL iom_put( "veiv_salttr", zztmp * zw2d )                  !  salt transport in j-direction 
    896       CALL iom_put( "veiv_salttr", zztmp * zw3d )                  !  salt transport in j-direction 
     891      CALL iom_put( "veiv_salttr"  , zztmp * zw2d )                  !  salt transport in j-direction 
     892      CALL iom_put( "veiv_salttr3d", zztmp * zw3d )                  !  salt transport in j-direction 
    897893      ! 
    898894      IF( iom_use( 'sopsteiv' ) ) CALL dia_ptr_hst( jp_sal, 'eiv', 0.5 * zw3d ) 
Note: See TracChangeset for help on using the changeset viewer.