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_lap_blp.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_lap_blp.F90

    r13295 r14037  
    1313   USE oce            ! ocean dynamics and active tracers 
    1414   USE dom_oce        ! ocean space and time domain 
     15   USE domutl, ONLY : is_tile 
    1516   USE ldftra         ! lateral physics: eddy diffusivity 
    1617   USE traldf_iso     ! iso-neutral lateral diffusion (standard operator)     (tra_ldf_iso   routine) 
     
    4647CONTAINS 
    4748 
    48    SUBROUTINE tra_ldf_lap( kt, Kmm, kit000, cdtype, pahu, pahv  ,               & 
    49       &                                             pgu , pgv   , pgui, pgvi,   & 
    50       &                                             pt  , pt_rhs, kjpt, kpass )  
     49   SUBROUTINE tra_ldf_lap( kt, Kmm, kit000, cdtype, pahu, pahv,             & 
     50      &                                             pgu , pgv , pgui, pgvi, & 
     51      &                                             pt, pt_rhs, kjpt, kpass ) 
     52      !! 
     53      INTEGER                     , INTENT(in   ) ::   kt         ! ocean time-step index 
     54      INTEGER                     , INTENT(in   ) ::   kit000     ! first time step index 
     55      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator) 
     56      INTEGER                     , INTENT(in   ) ::   kjpt       ! number of tracers 
     57      INTEGER                     , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
     58      INTEGER                     , INTENT(in   ) ::   Kmm        ! ocean time level index 
     59      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
     60      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
     61      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
     62      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pt         ! before tracer fields 
     63      REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) ::   pt_rhs     ! tracer trend 
     64      !! 
     65      CALL tra_ldf_lap_t( kt, Kmm, kit000, cdtype, pahu, pahv, is_tile(pahu),                            & 
     66      &                                            pgu , pgv , is_tile(pgu) , pgui, pgvi, is_tile(pgui), & 
     67      &                                            pt, is_tile(pt), pt_rhs, is_tile(pt_rhs), kjpt, kpass ) 
     68   END SUBROUTINE tra_ldf_lap 
     69 
     70 
     71   SUBROUTINE tra_ldf_lap_t( kt, Kmm, kit000, cdtype, pahu, pahv, ktah,                   & 
     72      &                                               pgu , pgv , ktg , pgui, pgvi, ktgi, & 
     73      &                                               pt, ktt, pt_rhs, ktt_rhs, kjpt, kpass ) 
    5174      !!---------------------------------------------------------------------- 
    5275      !!                  ***  ROUTINE tra_ldf_lap  *** 
     
    7295      INTEGER                              , INTENT(in   ) ::   kpass      ! =1/2 first or second passage 
    7396      INTEGER                              , INTENT(in   ) ::   Kmm        ! ocean time level index 
    74       REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
    75       REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
    76       REAL(wp), DIMENSION(jpi,jpj,    kjpt), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
    77       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   pt         ! before tracer fields 
    78       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pt_rhs     ! tracer trend  
    79       ! 
    80       INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices 
    81       REAL(wp) ::   zsign            ! local scalars 
    82       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ztu, ztv, zaheeu, zaheev 
    83       !!---------------------------------------------------------------------- 
    84       ! 
    85       IF( kt == nit000 .AND. lwp )  THEN 
    86          WRITE(numout,*) 
    87          WRITE(numout,*) 'tra_ldf_lap : iso-level laplacian diffusion on ', cdtype, ', pass=', kpass 
    88          WRITE(numout,*) '~~~~~~~~~~~ ' 
    89       ENDIF 
    90       ! 
    91       l_hst = .FALSE. 
    92       l_ptr = .FALSE. 
    93       IF( cdtype == 'TRA' .AND. ( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf' ) ) )     l_ptr = .TRUE. 
    94       IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
    95          &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE. 
     97      INTEGER                              , INTENT(in   ) ::   ktah, ktg, ktgi, ktt, ktt_rhs 
     98      REAL(wp), DIMENSION(A2D_T(ktah),   JPK)     , INTENT(in   ) ::   pahu, pahv ! eddy diffusivity at u- and v-points  [m2/s] 
     99      REAL(wp), DIMENSION(A2D_T(ktg),        KJPT), INTENT(in   ) ::   pgu, pgv   ! tracer gradient at pstep levels 
     100      REAL(wp), DIMENSION(A2D_T(ktgi),       KJPT), INTENT(in   ) ::   pgui, pgvi ! tracer gradient at top   levels 
     101      REAL(wp), DIMENSION(A2D_T(ktt),    JPK,KJPT), INTENT(in   ) ::   pt         ! before tracer fields 
     102      REAL(wp), DIMENSION(A2D_T(ktt_rhs),JPK,KJPT), INTENT(inout) ::   pt_rhs     ! tracer trend 
     103      ! 
     104      INTEGER  ::   ji, jj, jk, jn      ! dummy loop indices 
     105      INTEGER  ::   isi, iei, isj, iej  ! local integers 
     106      REAL(wp) ::   zsign               ! local scalars 
     107      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   ztu, ztv, zaheeu, zaheev 
     108      !!---------------------------------------------------------------------- 
     109      ! 
     110      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     111         IF( kt == nit000 .AND. lwp )  THEN 
     112            WRITE(numout,*) 
     113            WRITE(numout,*) 'tra_ldf_lap : iso-level laplacian diffusion on ', cdtype, ', pass=', kpass 
     114            WRITE(numout,*) '~~~~~~~~~~~ ' 
     115         ENDIF 
     116         ! 
     117         l_hst = .FALSE. 
     118         l_ptr = .FALSE. 
     119         IF( cdtype == 'TRA' .AND. ( iom_use( 'sophtldf' ) .OR. iom_use( 'sopstldf' ) ) )     l_ptr = .TRUE. 
     120         IF( cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR. & 
     121            &                        iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE. 
     122      ENDIF 
    96123      ! 
    97124      !                                !==  Initialization of metric arrays used for all tracers  ==! 
     
    99126      ELSE                    ;   zsign = -1._wp 
    100127      ENDIF 
    101       DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     128 
     129      IF( ntsi == Nis0 ) THEN ; isi = nn_hls - 1 ; ELSE ; isi = 0 ; ENDIF    ! Avoid double-counting when using tiling 
     130      IF( ntsj == Njs0 ) THEN ; isj = nn_hls - 1 ; ELSE ; isj = 0 ; ENDIF 
     131      IF( ntei == Nie0 ) THEN ; iei = nn_hls - 1 ; ELSE ; iei = 0 ; ENDIF 
     132      IF( ntej == Nje0 ) THEN ; iej = nn_hls - 1 ; ELSE ; iej = 0 ; ENDIF 
     133 
     134      DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )            !== First derivative (gradient)  ==! 
    102135         zaheeu(ji,jj,jk) = zsign * pahu(ji,jj,jk) * e2_e1u(ji,jj) * e3u(ji,jj,jk,Kmm)   !!gm   * umask(ji,jj,jk) pah masked! 
    103136         zaheev(ji,jj,jk) = zsign * pahv(ji,jj,jk) * e1_e2v(ji,jj) * e3v(ji,jj,jk,Kmm)   !!gm   * vmask(ji,jj,jk) 
     
    108141         !                          ! =========== !     
    109142         !                                
    110          DO_3D( 1, 0, 1, 0, 1, jpkm1 ) 
     143         DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )            !== First derivative (gradient)  ==! 
    111144            ztu(ji,jj,jk) = zaheeu(ji,jj,jk) * ( pt(ji+1,jj  ,jk,jn) - pt(ji,jj,jk,jn) ) 
    112145            ztv(ji,jj,jk) = zaheev(ji,jj,jk) * ( pt(ji  ,jj+1,jk,jn) - pt(ji,jj,jk,jn) ) 
    113146         END_3D 
    114          IF( ln_zps ) THEN                ! set gradient at bottom/top ocean level 
    115             DO_2D( 1, 0, 1, 0 ) 
     147         IF( ln_zps ) THEN                             ! set gradient at bottom/top ocean level 
     148            DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 )                              ! bottom 
    116149               ztu(ji,jj,mbku(ji,jj)) = zaheeu(ji,jj,mbku(ji,jj)) * pgu(ji,jj,jn) 
    117150               ztv(ji,jj,mbkv(ji,jj)) = zaheev(ji,jj,mbkv(ji,jj)) * pgv(ji,jj,jn) 
    118151            END_2D 
    119             IF( ln_isfcav ) THEN                ! top in ocean cavities only 
    120                DO_2D( 1, 0, 1, 0 ) 
     152            IF( ln_isfcav ) THEN                             ! top in ocean cavities only 
     153               DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 ) 
    121154                  IF( miku(ji,jj) > 1 )   ztu(ji,jj,miku(ji,jj)) = zaheeu(ji,jj,miku(ji,jj)) * pgui(ji,jj,jn)  
    122155                  IF( mikv(ji,jj) > 1 )   ztv(ji,jj,mikv(ji,jj)) = zaheev(ji,jj,mikv(ji,jj)) * pgvi(ji,jj,jn)  
     
    125158         ENDIF 
    126159         ! 
    127          DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     160         ! NOTE: [tiling-comms-merge] Bounds changed to avoid repeating this calculation for overlapping rows when using tiling 
     161         DO_3D( isj, iej, isi, iei, 1, jpkm1 )            !== Second derivative (divergence) added to the general tracer trends  ==! 
    128162            pt_rhs(ji,jj,jk,jn) = pt_rhs(ji,jj,jk,jn) + (  ztu(ji,jj,jk) - ztu(ji-1,jj,jk)     & 
    129163               &                                      +    ztv(ji,jj,jk) - ztv(ji,jj-1,jk) )   & 
     
    142176      !                             ! ================== 
    143177      ! 
    144    END SUBROUTINE tra_ldf_lap 
     178   END SUBROUTINE tra_ldf_lap_t 
    145179    
    146180 
     
    173207      ! 
    174208      INTEGER ::   ji, jj, jk, jn   ! dummy loop indices 
    175       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt) :: zlap         ! laplacian at t-point 
    176       REAL(wp), DIMENSION(jpi,jpj,    kjpt) :: zglu, zglv   ! bottom GRADh of the laplacian (u- and v-points) 
    177       REAL(wp), DIMENSION(jpi,jpj,    kjpt) :: zgui, zgvi   ! top    GRADh of the laplacian (u- and v-points) 
     209      REAL(wp), DIMENSION(A2D(nn_hls),jpk,kjpt) :: zlap         ! laplacian at t-point 
     210      REAL(wp), DIMENSION(A2D(nn_hls),    kjpt) :: zglu, zglv   ! bottom GRADh of the laplacian (u- and v-points) 
     211      REAL(wp), DIMENSION(A2D(nn_hls),    kjpt) :: zgui, zgvi   ! top    GRADh of the laplacian (u- and v-points) 
    178212      !!--------------------------------------------------------------------- 
    179213      ! 
    180       IF( kt == kit000 .AND. lwp )  THEN 
    181          WRITE(numout,*) 
    182          SELECT CASE ( kldf ) 
    183          CASE ( np_blp    )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-level   bilaplacian operator on ', cdtype 
    184          CASE ( np_blp_i  )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-neutral bilaplacian operator on ', cdtype, ' (Standard)' 
    185          CASE ( np_blp_it )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-neutral bilaplacian operator on ', cdtype, ' (triad)' 
    186          END SELECT 
    187          WRITE(numout,*) '~~~~~~~~~~~' 
     214      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     215         IF( kt == kit000 .AND. lwp )  THEN 
     216            WRITE(numout,*) 
     217            SELECT CASE ( kldf ) 
     218            CASE ( np_blp    )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-level   bilaplacian operator on ', cdtype 
     219            CASE ( np_blp_i  )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-neutral bilaplacian operator on ', cdtype, ' (Standard)' 
     220            CASE ( np_blp_it )   ;   WRITE(numout,*) 'tra_ldf_blp : iso-neutral bilaplacian operator on ', cdtype, ' (triad)' 
     221            END SELECT 
     222            WRITE(numout,*) '~~~~~~~~~~~' 
     223         ENDIF 
    188224      ENDIF 
    189225 
     
    200236      END SELECT 
    201237      ! 
     238      ! NOTE: [tiling-comms-merge] Needed for both nn_hls as tra_ldf_iso and tra_ldf_triad have not yet been adjusted to work with nn_hls = 2. In the zps case the lbc_lnk in zps_hde handles this, but in the zco case zlap always needs this lbc_lnk. I did try adjusting the bounds in tra_ldf_iso and tra_ldf_triad so this lbc_lnk was only needed for nn_hls = 1, but this was not correct and I did not have time to figure out why 
    202239      CALL lbc_lnk( 'traldf_lap_blp', zlap(:,:,:,:) , 'T', 1.0_wp )     ! Lateral boundary conditions (unchanged sign) 
    203240      !                                               ! Partial top/bottom cell: GRADh( zlap )   
Note: See TracChangeset for help on using the changeset viewer.