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 13516 for NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/TRA/traadv_qck.F90 – NEMO

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

Tiling for tra_adv

File:
1 edited

Legend:

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

    r13295 r13516  
    1717   USE oce             ! ocean dynamics and active tracers 
    1818   USE dom_oce         ! ocean space and time domain 
     19   ! TEMP: This change not necessary after trd_tra is tiled 
     20   USE domain, ONLY : dom_tile 
    1921   USE trc_oce         ! share passive tracers/Ocean variables 
    2022   USE trd_oce         ! trends: ocean variables 
     
    9193      INTEGER                                  , INTENT(in   ) ::   kjpt            ! number of tracers 
    9294      REAL(wp)                                 , INTENT(in   ) ::   p2dt            ! tracer time-step 
     95      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
    9396      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pU, pV, pW      ! 3 ocean volume transport components 
    9497      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! tracers and RHS of tracer equation 
    9598      !!---------------------------------------------------------------------- 
    9699      ! 
    97       IF( kt == kit000 )  THEN 
    98          IF(lwp) WRITE(numout,*) 
    99          IF(lwp) WRITE(numout,*) 'tra_adv_qck : 3rd order quickest advection scheme on ', cdtype 
    100          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
    101          IF(lwp) WRITE(numout,*) 
     100      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     101         IF( kt == kit000 )  THEN 
     102            IF(lwp) WRITE(numout,*) 
     103            IF(lwp) WRITE(numout,*) 'tra_adv_qck : 3rd order quickest advection scheme on ', cdtype 
     104            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
     105            IF(lwp) WRITE(numout,*) 
     106         ENDIF 
     107         ! 
     108         l_trd = .FALSE. 
     109         l_ptr = .FALSE. 
     110         IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )   l_trd = .TRUE. 
     111         IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) ) ) l_ptr = .TRUE. 
    102112      ENDIF 
    103       ! 
    104       l_trd = .FALSE. 
    105       l_ptr = .FALSE. 
    106       IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )   l_trd = .TRUE. 
    107       IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) ) ) l_ptr = .TRUE.  
    108113      ! 
    109114      ! 
     
    127132      INTEGER                                  , INTENT(in   ) ::   kjpt       ! number of tracers 
    128133      REAL(wp)                                 , INTENT(in   ) ::   p2dt       ! tracer time-step 
     134      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
    129135      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pU        ! i-velocity components 
    130136      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! active tracers and RHS of tracer equation 
    131137      !! 
     138      ! TEMP: This change not necessary after trd_tra is tiled 
     139      INTEGER  ::   itile 
    132140      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices 
    133141      REAL(wp) ::   ztra, zbtr, zdir, zdx, zmsk   ! local scalars 
    134       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zwx, zfu, zfc, zfd 
     142      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zwx, zfu, zfc, zfd 
     143      ! TEMP: This change not necessary after trd_tra is tiled 
     144      REAL(wp), DIMENSION(:,:,:), SAVE, ALLOCATABLE ::   ztrdx 
    135145      !---------------------------------------------------------------------- 
    136       ! 
     146      ! TEMP: This change not necessary after trd_tra is tiled 
     147      itile = ntile 
     148      ! 
     149      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
     150      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     151         IF( kt == nit000 .AND. l_trd ) THEN 
     152            ALLOCATE( ztrdx(jpi,jpj,jpk) ) 
     153         ENDIF 
     154      ENDIF 
    137155      !                                                          ! =========== 
    138156      DO jn = 1, kjpt                                            ! tracer loop 
     
    200218         END_3D 
    201219         !                                 ! trend diagnostics 
    202          IF( l_trd )   CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, zwx, pU, pt(:,:,:,jn,Kmm) ) 
     220         ! TEMP: These changes not necessary after trd_tra is tiled 
     221         IF( l_trd )  THEN 
     222            DO_3D( 1, 0, 1, 0, 1, jpk ) 
     223               ztrdx(ji,jj,jk) = zwx(ji,jj,jk) 
     224            END_3D 
     225 
     226            IF( ntile == 0 .OR. ntile == nijtile )  THEN                ! Do only for the full domain 
     227               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 )         ! Use full domain 
     228 
     229               ! TODO: TO BE TILED- trd_tra 
     230               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztrdx, pU, pt(:,:,:,jn,Kmm) ) 
     231 
     232               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = itile )     ! Revert to tile domain 
     233            ENDIF 
     234         ENDIF 
    203235         ! 
    204236      END DO 
     
    216248      INTEGER                                  , INTENT(in   ) ::   kjpt       ! number of tracers 
    217249      REAL(wp)                                 , INTENT(in   ) ::   p2dt       ! tracer time-step 
     250      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
    218251      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pV        ! j-velocity components 
    219252      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! active tracers and RHS of tracer equation 
    220253      !! 
     254      ! TEMP: This change not necessary after trd_tra is tiled 
     255      INTEGER  ::   itile 
    221256      INTEGER  :: ji, jj, jk, jn                ! dummy loop indices 
    222257      REAL(wp) :: ztra, zbtr, zdir, zdx, zmsk   ! local scalars 
    223       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zwy, zfu, zfc, zfd   ! 3D workspace 
     258      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zwy, zfu, zfc, zfd   ! 3D workspace 
     259      ! TEMP: This change not necessary after trd_tra is tiled 
     260      REAL(wp), DIMENSION(:,:,:), SAVE, ALLOCATABLE ::   ztrdy 
    224261      !---------------------------------------------------------------------- 
    225       ! 
     262      ! TEMP: This change not necessary after trd_tra is tiled 
     263      itile = ntile 
     264      ! 
     265      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
     266      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     267         IF( kt == nit000 .AND. l_trd ) THEN 
     268            ALLOCATE( ztrdy(jpi,jpj,jpk) ) 
     269         ENDIF 
     270      ENDIF 
    226271      !                                                          ! =========== 
    227272      DO jn = 1, kjpt                                            ! tracer loop 
     
    296341         END_3D 
    297342         !                                 ! trend diagnostics 
    298          IF( l_trd )   CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, zwy, pV, pt(:,:,:,jn,Kmm) ) 
     343         ! TEMP: These changes not necessary after trd_tra is tiled 
     344         IF( l_trd )  THEN 
     345            DO_3D( 1, 0, 1, 0, 1, jpk ) 
     346               ztrdy(ji,jj,jk) = zwy(ji,jj,jk) 
     347            END_3D 
     348 
     349            IF( ntile == 0 .OR. ntile == nijtile )  THEN                ! Do only for the full domain 
     350               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 )         ! Use full domain 
     351 
     352               ! TODO: TO BE TILED- trd_tra 
     353               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztrdy, pV, pt(:,:,:,jn,Kmm) ) 
     354 
     355               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = itile )     ! Revert to tile domain 
     356            ENDIF 
     357         ENDIF 
    299358         !                                 ! "Poleward" heat and salt transports (contribution of upstream fluxes) 
    300359         IF( l_ptr )   CALL dia_ptr_hst( jn, 'adv', zwy(:,:,:) ) 
     
    313372      CHARACTER(len=3)                         , INTENT(in   ) ::   cdtype   ! =TRA or TRC (tracer indicator) 
    314373      INTEGER                                  , INTENT(in   ) ::   kjpt     ! number of tracers 
    315       REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pW      ! vertical velocity  
     374      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
     375      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pW      ! vertical velocity 
    316376      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! active tracers and RHS of tracer equation 
    317377      ! 
     378      ! TEMP: This change not necessary after trd_tra is tiled 
     379      INTEGER  ::   itile 
    318380      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices 
    319       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zwz   ! 3D workspace 
    320       !!---------------------------------------------------------------------- 
    321       ! 
     381      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zwz   ! 3D workspace 
     382      ! TEMP: This change not necessary after trd_tra is tiled 
     383      REAL(wp), DIMENSION(:,:,:), SAVE, ALLOCATABLE ::   ztrdz 
     384      !!---------------------------------------------------------------------- 
     385      ! TEMP: This change not necessary after trd_tra is tiled 
     386      itile = ntile 
     387      ! 
     388      ! TEMP: This can be ST_2D(nn_hls) after trd_tra is tiled 
     389      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     390         IF( kt == nit000 .AND. l_trd ) THEN 
     391            ALLOCATE( ztrdz(jpi,jpj,jpk) ) 
     392         ENDIF 
     393      ENDIF 
     394 
    322395      zwz(:,:, 1 ) = 0._wp       ! surface & bottom values set to zero for all tracers 
    323396      zwz(:,:,jpk) = 0._wp 
     
    331404         END_3D 
    332405         IF( ln_linssh ) THEN                !* top value   (only in linear free surf. as zwz is multiplied by wmask) 
     406            ! TODO: NOT TESTED- requires isf 
    333407            IF( ln_isfcav ) THEN                  ! ice-shelf cavities (top of the ocean) 
    334408               DO_2D( 1, 1, 1, 1 ) 
     
    336410               END_2D 
    337411            ELSE                                   ! no ocean cavities (only ocean surface) 
    338                zwz(:,:,1) = pW(:,:,1) * pt(:,:,1,jn,Kmm) 
     412               DO_2D( 1, 1, 1, 1 ) 
     413                  zwz(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kmm) 
     414               END_2D 
    339415            ENDIF 
    340416         ENDIF 
     
    345421         END_3D 
    346422         !                                 ! Send trends for diagnostic 
    347          IF( l_trd )  CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, zwz, pW, pt(:,:,:,jn,Kmm) ) 
     423         ! TEMP: These changes not necessary after trd_tra is tiled 
     424         IF( l_trd )  THEN 
     425            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     426               ztrdz(ji,jj,jk) = zwz(ji,jj,jk) 
     427            END_3D 
     428 
     429            IF( ntile == 0 .OR. ntile == nijtile )  THEN                ! Do only for the full domain 
     430               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 )         ! Use full domain 
     431 
     432               ! TODO: TO BE TILED- trd_tra 
     433               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, ztrdz, pW, pt(:,:,:,jn,Kmm) ) 
     434 
     435               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = itile )     ! Revert to tile domain 
     436            ENDIF 
     437         ENDIF 
    348438         ! 
    349439      END DO 
     
    359449      !! ** Method :    
    360450      !!---------------------------------------------------------------------- 
    361       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfu   ! second upwind point 
    362       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfd   ! first douwning point 
    363       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfc   ! the central point (or the first upwind point) 
    364       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   puc   ! input as Courant number ; output as flux 
     451      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(in   ) ::   pfu   ! second upwind point 
     452      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(in   ) ::   pfd   ! first douwning point 
     453      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(in   ) ::   pfc   ! the central point (or the first upwind point) 
     454      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   puc   ! input as Courant number ; output as flux 
    365455      !! 
    366456      INTEGER  ::  ji, jj, jk               ! dummy loop indices  
Note: See TracChangeset for help on using the changeset viewer.