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 13515 for NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/LDF/ldftra.F90 – NEMO

Ignore:
Timestamp:
2020-09-24T20:32:14+02:00 (4 years ago)
Author:
hadcv
Message:

Tiling for tra_ldf

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/LDF/ldftra.F90

    r13295 r13515  
    726726      !! ** Action  : pu, pv increased by the eiv transport 
    727727      !!---------------------------------------------------------------------- 
    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(jpi,jpj,jpk), INTENT(inout) ::   pu      ! in : 3 ocean transport components   [m3/s] 
    733       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pv      ! out: 3 ocean transport components   [m3/s] 
    734       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(ST_2D(nn_hls),jpk), INTENT(inout) ::   pu        ! in : 3 ocean transport components   [m3/s] 
     733      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   pv        ! out: 3 ocean transport components   [m3/s] 
     734      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   pw        ! increased by the eiv                [m3/s] 
    735735      !! 
    736736      INTEGER  ::   ji, jj, jk                 ! dummy loop indices 
    737737      REAL(wp) ::   zuwk, zuwk1, zuwi, zuwi1   ! local scalars 
    738738      REAL(wp) ::   zvwk, zvwk1, zvwj, zvwj1   !   -      - 
    739       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpsi_uw, zpsi_vw 
    740       !!---------------------------------------------------------------------- 
    741       ! 
    742       IF( kt == kit000 )  THEN 
    743          IF(lwp) WRITE(numout,*) 
    744          IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :' 
    745          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component' 
     739      REAL(wp), DIMENSION(ST_2D(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 
    746748      ENDIF 
    747749 
     
    782784      !! 
    783785      !!---------------------------------------------------------------------- 
    784       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s] 
    785       INTEGER                         , INTENT(in   ) ::   Kmm   ! ocean time level indices 
     786      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s] 
     787      INTEGER                     , INTENT(in   ) ::   Kmm   ! ocean time level indices 
    786788      ! 
    787789      INTEGER  ::   ji, jj, jk    ! dummy loop indices 
    788790      REAL(wp) ::   zztmp   ! local scalar 
    789       REAL(wp), DIMENSION(jpi,jpj)     ::   zw2d   ! 2D workspace 
    790       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zw3d   ! 3D workspace 
     791      REAL(wp), DIMENSION(ST_2D(nn_hls))     ::   zw2d   ! 2D workspace 
     792      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zw3d   ! 3D workspace 
    791793      !!---------------------------------------------------------------------- 
    792794      ! 
     
    803805      zw3d(:,:,jpk) = 0._wp                                    ! bottom value always 0 
    804806      ! 
    805       DO jk = 1, jpkm1                                         ! e2u e3u u_eiv = -dk[psi_uw] 
    806          zw3d(:,:,jk) = ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) / ( e2u(:,:) * e3u(:,:,jk,Kmm) ) 
    807       END DO 
     807      DO_3D( 1, 1, 1, 1, 1, jpkm1 )                                  ! e2u e3u u_eiv = -dk[psi_uw] 
     808         zw3d(ji,jj,jk) = ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm) ) 
     809      END_3D 
    808810      CALL iom_put( "uoce_eiv", zw3d ) 
    809811      ! 
    810       DO jk = 1, jpkm1                                         ! e1v e3v v_eiv = -dk[psi_vw] 
    811          zw3d(:,:,jk) = ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) / ( e1v(:,:) * e3v(:,:,jk,Kmm) ) 
    812       END DO 
     812      DO_3D( 1, 1, 1, 1, 1, jpkm1 )                                  ! e1v e3v v_eiv = -dk[psi_vw] 
     813         zw3d(ji,jj,jk) = ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm) ) 
     814      END_3D 
    813815      CALL iom_put( "voce_eiv", zw3d ) 
    814816      ! 
     
    821823      ! 
    822824      IF( iom_use('weiv_masstr') ) THEN   ! vertical mass transport & its square value 
    823          zw2d(:,:) = rho0 * e1e2t(:,:) 
     825         DO_2D( 1, 1, 1, 1 ) 
     826            zw2d(ji,jj) = rho0 * e1e2t(ji,jj) 
     827         END_2D 
    824828         DO jk = 1, jpk 
    825829            zw3d(:,:,jk) = zw3d(:,:,jk) * zw2d(:,:) 
     
    867871      END_3D 
    868872      CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp ) 
    869       CALL iom_put( "veiv_heattr", zztmp * zw2d )                  !  heat transport in j-direction 
    870       CALL iom_put( "veiv_heattr", zztmp * zw3d )                  !  heat transport in j-direction 
     873      CALL lbc_lnk( 'ldftra', zw3d, 'V', -1. ) 
     874      CALL iom_put( "veiv_heattr"  , zztmp * zw2d )                  !  heat transport in j-direction 
     875      CALL iom_put( "veiv_heattr3d", zztmp * zw3d )                  !  heat transport in j-direction 
    871876      ! 
    872877      IF( iom_use( 'sophteiv' ) )   CALL dia_ptr_hst( jp_tem, 'eiv', 0.5 * zw3d ) 
     
    894899      END_3D 
    895900      CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp ) 
    896       CALL iom_put( "veiv_salttr", zztmp * zw2d )                  !  salt transport in j-direction 
    897       CALL iom_put( "veiv_salttr", zztmp * zw3d )                  !  salt transport in j-direction 
     901      CALL lbc_lnk( 'ldftra', zw3d, 'V', -1. ) 
     902      CALL iom_put( "veiv_salttr"  , zztmp * zw2d )                  !  salt transport in j-direction 
     903      CALL iom_put( "veiv_salttr3d", zztmp * zw3d )                  !  salt transport in j-direction 
    898904      ! 
    899905      IF( iom_use( 'sopsteiv' ) ) CALL dia_ptr_hst( jp_sal, 'eiv', 0.5 * zw3d ) 
Note: See TracChangeset for help on using the changeset viewer.