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 14037 for NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE/TRA/traldf_iso.F90 – NEMO

Ignore:
Timestamp:
2020-12-03T12:20:38+01:00 (3 years ago)
Author:
ayoung
Message:

Updated to trunk at 14020. Sette tests passed with change of results for configurations with non-linear ssh. Ticket #2506.

Location:
NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG

    • Property svn:externals
      •  

        old new  
        88 
        99# SETTE 
        10 ^/utils/CI/sette@13292        sette 
         10^/utils/CI/sette_wave@13990         sette 
  • NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE/TRA/traldf_iso.F90

    r13295 r14037  
    1919   USE oce            ! ocean dynamics and active tracers 
    2020   USE dom_oce        ! ocean space and time domain 
     21   USE domutl, ONLY : is_tile 
    2122   USE trc_oce        ! share passive tracers/Ocean variables 
    2223   USE zdf_oce        ! ocean vertical physics 
     
    4950CONTAINS 
    5051 
    51   SUBROUTINE tra_ldf_iso( kt, Kmm, kit000, cdtype, pahu, pahv,                    & 
    52       &                                            pgu , pgv    ,   pgui, pgvi,   & 
    53       &                                       pt , pt2 , pt_rhs , kjpt  , kpass ) 
     52   SUBROUTINE tra_ldf_iso( kt, Kmm, kit000, cdtype, pahu, pahv,             & 
     53      &                                             pgu , pgv , pgui, pgvi, & 
     54      &                                             pt, pt2, pt_rhs, kjpt, kpass ) 
     55      !! 
     56      INTEGER                     , INTENT(in   ) ::   kt         ! ocean time-step index 
     57      INTEGER                     , INTENT(in   ) ::   kit000     ! first time step index 
     58      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
     59      INTEGER                     , INTENT(in   ) ::   kjpt       ! number of tracers 
     60      INTEGER                     , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
     61      INTEGER                     , INTENT(in   ) ::   Kmm        ! ocean time level index 
     62      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
     63      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
     64      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
     65      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2) 
     66      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2) 
     67      REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) ::   pt_rhs     ! tracer trend 
     68      !! 
     69      CALL tra_ldf_iso_t( kt, Kmm, kit000, cdtype, pahu, pahv, is_tile(pahu),                             & 
     70         &                                         pgu , pgv , is_tile(pgu) , pgui, pgvi, is_tile(pgui),  & 
     71         &                                         pt, is_tile(pt), pt2, is_tile(pt2), pt_rhs, is_tile(pt_rhs), kjpt, kpass ) 
     72   END SUBROUTINE tra_ldf_iso 
     73 
     74 
     75  SUBROUTINE tra_ldf_iso_t( kt, Kmm, kit000, cdtype, pahu, pahv, ktah,                    & 
     76      &                                              pgu , pgv , ktg , pgui, pgvi, ktgi,  & 
     77      &                                              pt, ktt, pt2, ktt2, pt_rhs, ktt_rhs, kjpt, kpass ) 
    5478      !!---------------------------------------------------------------------- 
    5579      !!                  ***  ROUTINE tra_ldf_iso  *** 
     
    92116      !! ** Action :   Update pt_rhs arrays with the before rotated diffusion 
    93117      !!---------------------------------------------------------------------- 
    94       INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index 
    95       INTEGER                              , INTENT(in   ) ::   kit000     ! first time step index 
    96       CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
    97       INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers 
    98       INTEGER                              , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
    99       INTEGER                              , INTENT(in   ) ::   Kmm        ! ocean time level index 
    100       REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
    101       REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
    102       REAL(wp), DIMENSION(jpi,jpj,    kjpt), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
    103       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2) 
    104       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2) 
    105       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pt_rhs     ! tracer trend 
     118      INTEGER                                   , INTENT(in   ) ::   kt         ! ocean time-step index 
     119      INTEGER                                   , INTENT(in   ) ::   kit000     ! first time step index 
     120      CHARACTER(len=3)                          , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
     121      INTEGER                                   , INTENT(in   ) ::   kjpt       ! number of tracers 
     122      INTEGER                                   , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
     123      INTEGER                                   , INTENT(in   ) ::   Kmm        ! ocean time level index 
     124      INTEGER                                   , INTENT(in   ) ::   ktah, ktg, ktgi, ktt, ktt2, ktt_rhs 
     125      REAL(wp), DIMENSION(A2D_T(ktah)   ,JPK)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
     126      REAL(wp), DIMENSION(A2D_T(ktg)        ,KJPT), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
     127      REAL(wp), DIMENSION(A2D_T(ktgi)       ,KJPT), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
     128      REAL(wp), DIMENSION(A2D_T(ktt)    ,JPK,KJPT), INTENT(in   ) ::   pt         ! tracer (kpass=1) or laplacian of tracer (kpass=2) 
     129      REAL(wp), DIMENSION(A2D_T(ktt2)   ,JPK,KJPT), INTENT(in   ) ::   pt2        ! tracer (only used in kpass=2) 
     130      REAL(wp), DIMENSION(A2D_T(ktt_rhs),JPK,KJPT), INTENT(inout) ::   pt_rhs     ! tracer trend 
    106131      ! 
    107132      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices 
     
    111136      REAL(wp) ::  zmskv, zahv_w, zabe2, zcof2, zcoef4   !   -      - 
    112137      REAL(wp) ::  zcoef0, ze3w_2, zsign                 !   -      - 
    113       REAL(wp), DIMENSION(jpi,jpj)     ::   zdkt, zdk1t, z2d 
    114       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zdit, zdjt, zftu, zftv, ztfw  
     138      REAL(wp), DIMENSION(A2D(nn_hls))     ::   zdkt, zdk1t, z2d 
     139      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   zdit, zdjt, zftu, zftv, ztfw 
    115140      !!---------------------------------------------------------------------- 
    116141      ! 
    117142      IF( kpass == 1 .AND. kt == kit000 )  THEN 
    118          IF(lwp) WRITE(numout,*) 
    119          IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype 
    120          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
    121          ! 
    122          akz     (:,:,:) = 0._wp       
    123          ah_wslp2(:,:,:) = 0._wp 
     143         IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     144            IF(lwp) WRITE(numout,*) 
     145            IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype 
     146            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
     147         ENDIF 
     148         ! 
     149         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     150            akz     (ji,jj,jk) = 0._wp 
     151            ah_wslp2(ji,jj,jk) = 0._wp 
     152         END_3D 
    124153      ENDIF 
    125       !    
    126       l_hst = .FALSE. 
    127       l_ptr = .FALSE. 
    128       IF( cdtype == 'TRA' .AND. ( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf' ) ) )     l_ptr = .TRUE.  
    129       IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
    130          &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )   l_hst = .TRUE. 
     154      ! 
     155      IF( ntile == 0 .OR. ntile == 1 )  THEN                           ! Do only on the first tile 
     156         l_hst = .FALSE. 
     157         l_ptr = .FALSE. 
     158         IF( cdtype == 'TRA' .AND. ( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf' ) ) )     l_ptr = .TRUE. 
     159         IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
     160            &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )   l_hst = .TRUE. 
     161      ENDIF 
    131162      ! 
    132163      ! 
     
    167198            ! 
    168199            IF( ln_traldf_blp ) THEN                ! bilaplacian operator 
    169                DO_3D( 1, 0, 1, 0, 2, jpkm1 ) 
     200               DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
    170201                  akz(ji,jj,jk) = 16._wp   & 
    171202                     &   * ah_wslp2   (ji,jj,jk)   & 
     
    175206               END_3D 
    176207            ELSEIF( ln_traldf_lap ) THEN              ! laplacian operator 
    177                DO_3D( 1, 0, 1, 0, 2, jpkm1 ) 
     208               DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
    178209                  ze3w_2 = e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm) 
    179210                  zcoef0 = rDt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
     
    183214           ! 
    184215         ELSE                                    ! 33 flux set to zero with akz=ah_wslp2 ==>> computed in full implicit 
    185             akz(:,:,:) = ah_wslp2(:,:,:)       
     216            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     217               akz(ji,jj,jk) = ah_wslp2(ji,jj,jk) 
     218            END_3D 
    186219         ENDIF 
    187220      ENDIF 
     
    195228         !!---------------------------------------------------------------------- 
    196229!!gm : bug.... why (x,:,:)?   (1,jpj,:) and (jpi,1,:) should be sufficient.... 
    197          zdit (1,:,:) = 0._wp     ;     zdit (jpi,:,:) = 0._wp 
    198          zdjt (1,:,:) = 0._wp     ;     zdjt (jpi,:,:) = 0._wp 
     230         zdit (ntsi-nn_hls,:,:) = 0._wp     ;     zdit (ntei+nn_hls,:,:) = 0._wp 
     231         zdjt (ntsi-nn_hls,:,:) = 0._wp     ;     zdjt (ntei+nn_hls,:,:) = 0._wp 
    199232         !!end 
    200233 
     
    205238         END_3D 
    206239         IF( ln_zps ) THEN      ! botton and surface ocean correction of the horizontal gradient 
    207             DO_2D( 1, 0, 1, 0 ) 
     240            DO_2D( 1, 0, 1, 0 )           ! bottom correction (partial bottom cell) 
    208241               zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)           
    209242               zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn) 
     
    223256         DO jk = 1, jpkm1                                 ! Horizontal slab 
    224257            ! 
    225             !                             !== Vertical tracer gradient 
    226             zdk1t(:,:) = ( pt(:,:,jk,jn) - pt(:,:,jk+1,jn) ) * wmask(:,:,jk+1)     ! level jk+1 
    227             ! 
    228             IF( jk == 1 ) THEN   ;   zdkt(:,:) = zdk1t(:,:)                          ! surface: zdkt(jk=1)=zdkt(jk=2) 
    229             ELSE                 ;   zdkt(:,:) = ( pt(:,:,jk-1,jn) - pt(:,:,jk,jn) ) * wmask(:,:,jk) 
    230             ENDIF 
    231             DO_2D( 1, 0, 1, 0 ) 
     258            DO_2D( 1, 1, 1, 1 ) 
     259               !                             !== Vertical tracer gradient 
     260               zdk1t(ji,jj) = ( pt(ji,jj,jk,jn) - pt(ji,jj,jk+1,jn) ) * wmask(ji,jj,jk+1)     ! level jk+1 
     261               ! 
     262               IF( jk == 1 ) THEN   ;   zdkt(ji,jj) = zdk1t(ji,jj)                            ! surface: zdkt(jk=1)=zdkt(jk=2) 
     263               ELSE                 ;   zdkt(ji,jj) = ( pt(ji,jj,jk-1,jn) - pt(ji,jj,jk,jn) ) * wmask(ji,jj,jk) 
     264               ENDIF 
     265            END_2D 
     266            ! 
     267            DO_2D( 1, 0, 1, 0 )           !==  Horizontal fluxes 
    232268               zabe1 = pahu(ji,jj,jk) * e2_e1u(ji,jj) * e3u(ji,jj,jk,Kmm) 
    233269               zabe2 = pahv(ji,jj,jk) * e1_e2v(ji,jj) * e3v(ji,jj,jk,Kmm) 
     
    250286            END_2D 
    251287            ! 
    252             DO_2D( 0, 0, 0, 0 ) 
     288            DO_2D( 0, 0, 0, 0 )           !== horizontal divergence and add to pta 
    253289               pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn)    & 
    254290                  &       + zsign * (  zftu(ji,jj,jk) - zftu(ji-1,jj,jk) + zftv(ji,jj,jk) - zftv(ji,jj-1,jk)  )   & 
     
    266302         ztfw(:,:, 1 ) = 0._wp      ;      ztfw(:,:,jpk) = 0._wp 
    267303          
    268          DO_3D( 0, 0, 0, 0, 2, jpkm1 ) 
     304         DO_3D( 0, 0, 0, 0, 2, jpkm1 )    ! interior (2=<jk=<jpk-1) 
    269305            ! 
    270306            zmsku = wmask(ji,jj,jk) / MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)          & 
     
    311347         ENDIF 
    312348         !          
    313          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     349         DO_3D( 0, 0, 0, 0, 1, jpkm1 )    !==  Divergence of vertical fluxes added to pta  ==! 
    314350            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn) + zsign * (  ztfw (ji,jj,jk) - ztfw(ji,jj,jk+1)  ) * r1_e1e2t(ji,jj)   & 
    315351               &                                             / e3t(ji,jj,jk,Kmm) 
     
    330366      END DO                                                      ! end tracer loop 
    331367      ! 
    332    END SUBROUTINE tra_ldf_iso 
     368   END SUBROUTINE tra_ldf_iso_t 
    333369 
    334370   !!============================================================================== 
Note: See TracChangeset for help on using the changeset viewer.