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 13409 for NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/TRA/traadv_ubs.F90 – NEMO

Ignore:
Timestamp:
2020-08-17T15:28:54+02:00 (4 years ago)
Author:
hadcv
Message:

Remaining changes prior to trunk merge

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/TRA/traadv_ubs.F90

    r12377 r13409  
    1414   USE oce            ! ocean dynamics and active tracers 
    1515   USE dom_oce        ! ocean space and time domain 
     16   ! TEMP: This change not necessary after trd_tra is tiled 
     17   USE domain, ONLY : dom_tile 
    1618   USE trc_oce        ! share passive tracers/Ocean variables 
    1719   USE trd_oce        ! trends: ocean variables 
     
    9193      INTEGER                                  , INTENT(in   ) ::   kn_ubs_v        ! number of tracers 
    9294      REAL(wp)                                 , INTENT(in   ) ::   p2dt            ! tracer time-step 
     95      ! TEMP: This can be A2D 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      ! 
     99      ! TEMP: This change not necessary after trd_tra is tiled 
     100      INTEGER  ::   itile 
    96101      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices 
    97102      REAL(wp) ::   ztra, zbtr, zcoef                       ! local scalars 
    98103      REAL(wp) ::   zfp_ui, zfm_ui, zcenut, ztak, zfp_wk, zfm_wk   !   -      - 
    99104      REAL(wp) ::   zfp_vj, zfm_vj, zcenvt, zeeu, zeev, z_hdivn    !   -      - 
    100       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ztu, ztv, zltu, zltv, zti, ztw   ! 3D workspace 
    101       !!---------------------------------------------------------------------- 
    102       ! 
    103       IF( kt == kit000 )  THEN 
    104          IF(lwp) WRITE(numout,*) 
    105          IF(lwp) WRITE(numout,*) 'tra_adv_ubs :  horizontal UBS advection scheme on ', cdtype 
    106          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
     105      REAL(wp), DIMENSION(A2D,jpk) ::   ztu, ztv, zltu, zltv, zti, ztw     ! 3D workspace 
     106      ! TEMP: This change not necessary after trd_tra is tiled 
     107      REAL(wp), DIMENSION(:,:,:), SAVE, ALLOCATABLE ::   ztrdx, ztrdy, ztrdz 
     108      !!---------------------------------------------------------------------- 
     109      ! TEMP: This change not necessary after trd_tra is tiled 
     110      itile = ntile 
     111      ! 
     112      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     113         IF( kt == kit000 )  THEN 
     114            IF(lwp) WRITE(numout,*) 
     115            IF(lwp) WRITE(numout,*) 'tra_adv_ubs :  horizontal UBS advection scheme on ', cdtype 
     116            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
     117         ENDIF 
     118         ! 
     119         l_trd = .FALSE. 
     120         l_hst = .FALSE. 
     121         l_ptr = .FALSE. 
     122         IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )      l_trd = .TRUE. 
     123         IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) )  )   l_ptr = .TRUE. 
     124         IF(   cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
     125            &                          iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) ) l_hst = .TRUE. 
     126 
     127         ! TEMP: This can be A2D after trd_tra is tiled 
     128         IF( kt == kit000 .AND. l_trd ) THEN 
     129            ALLOCATE( ztrdx(jpi,jpj,jpk), ztrdy(jpi,jpj,jpk), ztrdz(jpi,jpj,jpk) ) 
     130         ENDIF 
    107131      ENDIF 
    108       ! 
    109       l_trd = .FALSE. 
    110       l_hst = .FALSE. 
    111       l_ptr = .FALSE. 
    112       IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )      l_trd = .TRUE. 
    113       IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) )  )   l_ptr = .TRUE.  
    114       IF(   cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
    115          &                          iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) ) l_hst = .TRUE. 
    116132      ! 
    117133      ztw (:,:, 1 ) = 0._wp      ! surface & bottom value : set to zero for all tracers 
     
    152168         END_3D 
    153169         ! 
    154          zltu(:,:,:) = pt(:,:,:,jn,Krhs)      ! store the initial trends before its update 
     170         DO_3D_11_11( 1, jpk ) 
     171            zltu(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs)      ! store the initial trends before its update 
     172         END_3D 
    155173         ! 
    156174         DO jk = 1, jpkm1        !==  add the horizontal advective trend  ==! 
     
    163181         END DO 
    164182         ! 
    165          zltu(:,:,:) = pt(:,:,:,jn,Krhs) - zltu(:,:,:)    ! Horizontal advective trend used in vertical 2nd order FCT case 
    166          !                                            ! and/or in trend diagnostic (l_trd=T)  
    167          !                 
     183         DO_3D_11_11( 1, jpk ) 
     184            zltu(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs) - zltu(ji,jj,jk)  ! Horizontal advective trend used in vertical 2nd order FCT case 
     185         END_3D                                                     ! and/or in trend diagnostic (l_trd=T) 
     186         ! 
     187         ! TEMP: These changes not necessary after trd_tra is tiled 
    168188         IF( l_trd ) THEN                  ! trend diagnostics 
    169              CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztu, pU, pt(:,:,:,jn,Kmm) ) 
    170              CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztv, pV, pt(:,:,:,jn,Kmm) ) 
     189            DO_3D_11_11( 1, jpk ) 
     190               ztrdx(ji,jj,jk) = ztu(ji,jj,jk) 
     191               ztrdy(ji,jj,jk) = ztv(ji,jj,jk) 
     192            END_3D 
     193 
     194            IF( ntile == 0 .OR. ntile == nijtile )  THEN                ! Do only for the full domain 
     195               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 )         ! Use full domain 
     196 
     197               ! TODO: TO BE TILED- trd_tra 
     198               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztrdx, pU, pt(:,:,:,jn,Kmm) ) 
     199               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztrdy, pV, pt(:,:,:,jn,Kmm) ) 
     200 
     201               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = itile )     ! Revert to tile domain 
     202            ENDIF 
    171203         END IF 
    172204         !      
     
    183215         CASE(  2  )                   ! 2nd order FCT  
    184216            !          
    185             IF( l_trd )   zltv(:,:,:) = pt(:,:,:,jn,Krhs)          ! store pt(:,:,:,:,Krhs) if trend diag. 
     217            IF( l_trd ) THEN 
     218               DO_3D_11_11( 1, jpk ) 
     219                  zltv(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs)          ! store pt(:,:,:,:,Krhs) if trend diag. 
     220               END_3D 
     221            ENDIF 
    186222            ! 
    187223            !                          !*  upstream advection with initial mass fluxes & intermediate update  ==! 
     
    192228            END_3D 
    193229            IF( ln_linssh ) THEN             ! top ocean value (only in linear free surface as ztw has been w-masked) 
     230               ! TODO: NOT TESTED- requires isf 
    194231               IF( ln_isfcav ) THEN                ! top of the ice-shelf cavities and at the ocean surface 
    195232                  DO_2D_11_11 
     
    197234                  END_2D 
    198235               ELSE                                ! no cavities: only at the ocean surface 
    199                   ztw(:,:,1) = pW(:,:,1) * pt(:,:,1,jn,Kbb) 
     236                  DO_2D_11_11 
     237                     ztw(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kbb) 
     238                  END_2D 
    200239               ENDIF 
    201240            ENDIF 
     
    223262               ztw(ji,jj,jk) = pW(ji,jj,jk) * ztw(ji,jj,jk) * wmask(ji,jj,jk) 
    224263            END_3D 
    225             IF( ln_linssh )   ztw(:,:, 1 ) = pW(:,:,1) * pt(:,:,1,jn,Kmm)     !!gm ISF & 4th COMPACT doesn't work 
     264            IF( ln_linssh ) THEN 
     265               DO_2D_11_11 
     266                  ztw(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kmm)     !!gm ISF & 4th COMPACT doesn't work 
     267               END_2D 
     268            ENDIF 
    226269            ! 
    227270         END SELECT 
     
    231274         END_3D 
    232275         ! 
     276         ! TEMP: These changes not necessary after trd_tra is tiled 
    233277         IF( l_trd )  THEN       ! vertical advective trend diagnostics 
    234278            DO_3D_00_00( 1, jpkm1 ) 
    235                zltv(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs) - zltv(ji,jj,jk)                          & 
     279               ztrdz(ji,jj,jk) = pt(ji,jj,jk,jn,Krhs) - zltv(ji,jj,jk)                          & 
    236280                  &           + pt(ji,jj,jk,jn,Kmm) * (  pW(ji,jj,jk) - pW(ji,jj,jk+1)  )   & 
    237281                  &                              * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
    238282            END_3D 
    239             CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, zltv ) 
     283 
     284            IF( ntile == 0 .OR. ntile == nijtile )  THEN                ! Do only for the full domain 
     285               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 )         ! Use full domain 
     286 
     287               ! TODO: TO BE TILED- trd_tra 
     288               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, ztrdz ) 
     289 
     290               IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = nijtile )   ! Revert to tile domain 
     291            ENDIF 
    240292         ENDIF 
    241293         ! 
     
    258310      !!       in-space based differencing for fluid 
    259311      !!---------------------------------------------------------------------- 
    260       INTEGER , INTENT(in   )                          ::   Kmm    ! time level index 
    261       REAL(wp), INTENT(in   )                          ::   p2dt   ! tracer time-step 
    262       REAL(wp),                DIMENSION (jpi,jpj,jpk) ::   pbef   ! before field 
    263       REAL(wp), INTENT(inout), DIMENSION (jpi,jpj,jpk) ::   paft   ! after field 
    264       REAL(wp), INTENT(inout), DIMENSION (jpi,jpj,jpk) ::   pcc    ! monotonic flux in the k direction 
     312      INTEGER , INTENT(in   )                         ::   Kmm    ! time level index 
     313      REAL(wp), INTENT(in   )                         ::   p2dt   ! tracer time-step 
     314      REAL(wp),                DIMENSION(jpi,jpj,jpk) ::   pbef   ! before field 
     315      REAL(wp), INTENT(inout), DIMENSION(A2D    ,jpk) ::   paft   ! after field 
     316      REAL(wp), INTENT(inout), DIMENSION(A2D    ,jpk) ::   pcc    ! monotonic flux in the k direction 
    265317      ! 
    266318      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    267319      INTEGER  ::   ikm1         ! local integer 
    268320      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn   ! local scalars 
    269       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zbetup, zbetdo     ! 3D workspace 
     321      REAL(wp), DIMENSION(A2D,jpk) ::   zbetup, zbetdo         ! 3D workspace 
    270322      !!---------------------------------------------------------------------- 
    271323      ! 
     
    277329      ! -------------------- 
    278330      !                    ! large negative value (-zbig) inside land 
    279       pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:) - zbig * ( 1.e0 - tmask(:,:,:) ) 
    280       paft(:,:,:) = paft(:,:,:) * tmask(:,:,:) - zbig * ( 1.e0 - tmask(:,:,:) ) 
     331      DO_3D_00_00( 1, jpk ) 
     332         pbef(ji,jj,jk) = pbef(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1.e0 - tmask(ji,jj,jk) ) 
     333         paft(ji,jj,jk) = paft(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1.e0 - tmask(ji,jj,jk) ) 
     334      END_3D 
    281335      ! 
    282336      DO jk = 1, jpkm1     ! search maximum in neighbourhood 
     
    289343      END DO 
    290344      !                    ! large positive value (+zbig) inside land 
    291       pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:) + zbig * ( 1.e0 - tmask(:,:,:) ) 
    292       paft(:,:,:) = paft(:,:,:) * tmask(:,:,:) + zbig * ( 1.e0 - tmask(:,:,:) ) 
     345      DO_3D_00_00( 1, jpk ) 
     346         pbef(ji,jj,jk) = pbef(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1.e0 - tmask(ji,jj,jk) ) 
     347         paft(ji,jj,jk) = paft(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1.e0 - tmask(ji,jj,jk) ) 
     348      END_3D 
    293349      ! 
    294350      DO jk = 1, jpkm1     ! search minimum in neighbourhood 
     
    301357      END DO 
    302358      !                    ! restore masked values to zero 
    303       pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:) 
    304       paft(:,:,:) = paft(:,:,:) * tmask(:,:,:) 
     359      DO_3D_00_00( 1, jpk ) 
     360         pbef(ji,jj,jk) = pbef(ji,jj,jk) * tmask(ji,jj,jk) 
     361         paft(ji,jj,jk) = paft(ji,jj,jk) * tmask(ji,jj,jk) 
     362      END_3D 
    305363      ! 
    306364      ! Positive and negative part of fluxes and beta terms 
Note: See TracChangeset for help on using the changeset viewer.