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 12928 for NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser/src/OCE/DYN/wet_dry.F90 – NEMO

Ignore:
Timestamp:
2020-05-14T21:46:00+02:00 (4 years ago)
Author:
smueller
Message:

Synchronizing with /NEMO/trunk@12925 (ticket #2170)

Location:
NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser

    • Property svn:externals
      •  

        old new  
        66^/vendors/FCM@HEAD            ext/FCM 
        77^/vendors/IOIPSL@HEAD         ext/IOIPSL 
         8 
         9# SETTE 
         10^/utils/CI/sette@HEAD         sette 
  • NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser/src/OCE/DYN/wet_dry.F90

    r12178 r12928  
    3131   PRIVATE 
    3232 
     33   !! * Substitutions 
     34#  include "do_loop_substitute.h90" 
    3335   !!---------------------------------------------------------------------- 
    3436   !! critical depths,filters, limiters,and masks for  Wetting and Drying 
     
    6163 
    6264   !! * Substitutions 
    63 #  include "vectopt_loop_substitute.h90" 
    6465   !!---------------------------------------------------------------------- 
    6566CONTAINS 
     
    7980      !!---------------------------------------------------------------------- 
    8081      ! 
    81       REWIND( numnam_ref )              ! Namelist namwad in reference namelist : Parameters for Wetting/Drying 
    8282      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905) 
    8383905   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namwad in reference namelist' )  
    84       REWIND( numnam_cfg )              ! Namelist namwad in configuration namelist : Parameters for Wetting/Drying 
    8584      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906) 
    8685906   IF( ios >  0 )   CALL ctl_nam ( ios , 'namwad in configuration namelist' ) 
     
    122121 
    123122 
    124    SUBROUTINE wad_lmt( sshb1, sshemp, z2dt ) 
     123   SUBROUTINE wad_lmt( psshb1, psshemp, z2dt, Kmm, puu, pvv ) 
    125124      !!---------------------------------------------------------------------- 
    126125      !!                  ***  ROUTINE wad_lmt  *** 
     
    132131      !! ** Action  : - calculate flux limiter and W/D flag 
    133132      !!---------------------------------------------------------------------- 
    134       REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1        !!gm DOCTOR names: should start with p ! 
    135       REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sshemp 
    136       REAL(wp)                , INTENT(in   ) ::   z2dt 
     133      REAL(wp), DIMENSION(:,:)            , INTENT(inout) ::   psshb1 
     134      REAL(wp), DIMENSION(:,:)            , INTENT(in   ) ::   psshemp 
     135      REAL(wp)                            , INTENT(in   ) ::   z2dt 
     136      INTEGER                             , INTENT(in   ) ::   Kmm       ! time level index 
     137      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv  ! velocity arrays 
    137138      ! 
    138139      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices 
     
    150151      ! 
    151152      DO jk = 1, jpkm1 
    152          un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:)  
    153          vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:)  
     153         puu(:,:,jk,Kmm) = puu(:,:,jk,Kmm) * zwdlmtu(:,:)  
     154         pvv(:,:,jk,Kmm) = pvv(:,:,jk,Kmm) * zwdlmtv(:,:)  
    154155      END DO 
    155156      jflag  = 0 
     
    165166      ! 
    166167      DO jk = 1, jpkm1     ! Horizontal Flux in u and v direction 
    167          zflxu(:,:) = zflxu(:,:) + e3u_n(:,:,jk) * un(:,:,jk) * umask(:,:,jk) 
    168          zflxv(:,:) = zflxv(:,:) + e3v_n(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk) 
     168         zflxu(:,:) = zflxu(:,:) + e3u(:,:,jk,Kmm) * puu(:,:,jk,Kmm) * umask(:,:,jk) 
     169         zflxv(:,:) = zflxv(:,:) + e3v(:,:,jk,Kmm) * pvv(:,:,jk,Kmm) * vmask(:,:,jk) 
    169170      END DO 
    170171      zflxu(:,:) = zflxu(:,:) * e2u(:,:) 
     
    172173      ! 
    173174      wdmask(:,:) = 1._wp 
    174       DO jj = 2, jpj 
    175          DO ji = 2, jpi  
    176             ! 
    177             IF( tmask(ji,jj,1)        < 0.5_wp )   CYCLE    ! we don't care about land cells 
    178             IF( ht_0(ji,jj) - ssh_ref > zdepwd )   CYCLE    ! and cells which are unlikely to dry 
    179             ! 
    180             zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   & 
    181                &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp )  
    182             zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   & 
    183                &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp )  
    184             ! 
    185             zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 
    186             IF( zdep2 <= 0._wp ) THEN     ! add more safty, but not necessary 
    187                sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    188                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
    189                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
    190                IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp 
    191                IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp  
    192                wdmask(ji,jj) = 0._wp 
    193             END IF 
    194          END DO 
    195       END DO 
     175      DO_2D_01_01 
     176         ! 
     177         IF( tmask(ji,jj,1)        < 0.5_wp )   CYCLE    ! we don't care about land cells 
     178         IF( ht_0(ji,jj) - ssh_ref > zdepwd )   CYCLE    ! and cells which are unlikely to dry 
     179         ! 
     180         zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   & 
     181            &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp )  
     182         zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   & 
     183            &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp )  
     184         ! 
     185         zdep2 = ht_0(ji,jj) + psshb1(ji,jj) - rn_wdmin1 
     186         IF( zdep2 <= 0._wp ) THEN     ! add more safty, but not necessary 
     187            psshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
     188            IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
     189            IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
     190            IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp 
     191            IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp  
     192            wdmask(ji,jj) = 0._wp 
     193         END IF 
     194      END_2D 
    196195      ! 
    197196      !           ! HPG limiter from jholt 
    198       wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp) 
     197      wdramp(:,:) = min((ht_0(:,:) + psshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp) 
    199198      !jth assume don't need a lbc_lnk here 
    200       DO jj = 1, jpjm1 
    201          DO ji = 1, jpim1 
    202             wdrampu(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji+1,jj) ) 
    203             wdrampv(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji,jj+1) ) 
    204          END DO 
    205       END DO 
     199      DO_2D_10_10 
     200         wdrampu(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji+1,jj) ) 
     201         wdrampv(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji,jj+1) ) 
     202      END_2D 
    206203      !           ! end HPG limiter 
    207204      ! 
     
    213210         jflag = 0     ! flag indicating if any further iterations are needed 
    214211         ! 
    215          DO jj = 2, jpj 
    216             DO ji = 2, jpi  
    217                IF( tmask(ji, jj, 1) < 0.5_wp )   CYCLE  
    218                IF( ht_0(ji,jj)      > zdepwd )   CYCLE 
    219                ! 
    220                ztmp = e1e2t(ji,jj) 
    221                ! 
    222                zzflxp = MAX( zflxu1(ji,jj) , 0._wp ) - MIN( zflxu1(ji-1,jj  ) , 0._wp)   & 
    223                   &   + MAX( zflxv1(ji,jj) , 0._wp ) - MIN( zflxv1(ji,  jj-1) , 0._wp)  
    224                zzflxn = MIN( zflxu1(ji,jj) , 0._wp ) - MAX( zflxu1(ji-1,jj  ) , 0._wp)   & 
    225                   &   + MIN( zflxv1(ji,jj) , 0._wp ) - MAX( zflxv1(ji,  jj-1) , 0._wp)  
    226                ! 
    227                zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    228                zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj) 
    229                ! 
    230                IF( zdep1 > zdep2 ) THEN 
    231                   wdmask(ji, jj) = 0._wp 
    232                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt ) 
    233                   !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt ) 
    234                   ! flag if the limiter has been used but stop flagging if the only 
    235                   ! changes have zeroed the coefficient since further iterations will 
    236                   ! not change anything 
    237                   IF( zcoef > 0._wp ) THEN   ;   jflag = 1  
    238                   ELSE                       ;   zcoef = 0._wp 
    239                   ENDIF 
    240                   IF( jk1 > nn_wdit )   zcoef = 0._wp 
    241                   IF( zflxu1(ji  ,jj  ) > 0._wp )   zwdlmtu(ji  ,jj  ) = zcoef 
    242                   IF( zflxu1(ji-1,jj  ) < 0._wp )   zwdlmtu(ji-1,jj  ) = zcoef 
    243                   IF( zflxv1(ji  ,jj  ) > 0._wp )   zwdlmtv(ji  ,jj  ) = zcoef 
    244                   IF( zflxv1(ji  ,jj-1) < 0._wp )   zwdlmtv(ji  ,jj-1) = zcoef 
     212         DO_2D_01_01 
     213            IF( tmask(ji, jj, 1) < 0.5_wp )   CYCLE  
     214            IF( ht_0(ji,jj)      > zdepwd )   CYCLE 
     215            ! 
     216            ztmp = e1e2t(ji,jj) 
     217            ! 
     218            zzflxp = MAX( zflxu1(ji,jj) , 0._wp ) - MIN( zflxu1(ji-1,jj  ) , 0._wp)   & 
     219               &   + MAX( zflxv1(ji,jj) , 0._wp ) - MIN( zflxv1(ji,  jj-1) , 0._wp)  
     220            zzflxn = MIN( zflxu1(ji,jj) , 0._wp ) - MAX( zflxu1(ji-1,jj  ) , 0._wp)   & 
     221               &   + MIN( zflxv1(ji,jj) , 0._wp ) - MAX( zflxv1(ji,  jj-1) , 0._wp)  
     222            ! 
     223            zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
     224            zdep2 = ht_0(ji,jj) + psshb1(ji,jj) - rn_wdmin1 - z2dt * psshemp(ji,jj) 
     225            ! 
     226            IF( zdep1 > zdep2 ) THEN 
     227               wdmask(ji, jj) = 0._wp 
     228               zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt ) 
     229               !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt ) 
     230               ! flag if the limiter has been used but stop flagging if the only 
     231               ! changes have zeroed the coefficient since further iterations will 
     232               ! not change anything 
     233               IF( zcoef > 0._wp ) THEN   ;   jflag = 1  
     234               ELSE                       ;   zcoef = 0._wp 
    245235               ENDIF 
    246             END DO 
    247          END DO 
     236               IF( jk1 > nn_wdit )   zcoef = 0._wp 
     237               IF( zflxu1(ji  ,jj  ) > 0._wp )   zwdlmtu(ji  ,jj  ) = zcoef 
     238               IF( zflxu1(ji-1,jj  ) < 0._wp )   zwdlmtu(ji-1,jj  ) = zcoef 
     239               IF( zflxv1(ji  ,jj  ) > 0._wp )   zwdlmtv(ji  ,jj  ) = zcoef 
     240               IF( zflxv1(ji  ,jj-1) < 0._wp )   zwdlmtv(ji  ,jj-1) = zcoef 
     241            ENDIF 
     242         END_2D 
    248243         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. ) 
    249244         ! 
     
    255250      ! 
    256251      DO jk = 1, jpkm1 
    257          un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:)  
    258          vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:)  
     252         puu(:,:,jk,Kmm) = puu(:,:,jk,Kmm) * zwdlmtu(:,:)  
     253         pvv(:,:,jk,Kmm) = pvv(:,:,jk,Kmm) * zwdlmtv(:,:)  
    259254      END DO 
    260       un_b(:,:) = un_b(:,:) * zwdlmtu(:, :) 
    261       vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :) 
     255      uu_b(:,:,Kmm) = uu_b(:,:,Kmm) * zwdlmtu(:, :) 
     256      vv_b(:,:,Kmm) = vv_b(:,:,Kmm) * zwdlmtv(:, :) 
    262257      ! 
    263258!!gm TO BE SUPPRESSED ?  these lbc_lnk are useless since zwdlmtu and zwdlmtv are defined everywhere ! 
    264       CALL lbc_lnk_multi( 'wet_dry', un  , 'U', -1., vn  , 'V', -1. ) 
    265       CALL lbc_lnk_multi( 'wet_dry', un_b, 'U', -1., vn_b, 'V', -1. ) 
     259      CALL lbc_lnk_multi( 'wet_dry', puu(:,:,:,Kmm)  , 'U', -1., pvv(:,:,:,Kmm)  , 'V', -1. ) 
     260      CALL lbc_lnk_multi( 'wet_dry', uu_b(:,:,Kmm), 'U', -1., vv_b(:,:,Kmm), 'V', -1. ) 
    266261!!gm 
    267262      ! 
    268263      IF(jflag == 1 .AND. lwp)   WRITE(numout,*) 'Need more iterations in wad_lmt!!!' 
    269264      ! 
    270       !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field) 
     265      !IF( ln_rnf      )   CALL sbc_rnf_div( hdiv )          ! runoffs (update hdiv field) 
    271266      ! 
    272267      IF( ln_timing )   CALL timing_stop('wad_lmt')      ! 
     
    275270 
    276271 
    277    SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt ) 
     272   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rDt_e ) 
    278273      !!---------------------------------------------------------------------- 
    279274      !!                  ***  ROUTINE wad_lmt  *** 
     
    285280      !! ** Action  : - calculate flux limiter and W/D flag 
    286281      !!---------------------------------------------------------------------- 
    287       REAL(wp)                , INTENT(in   ) ::   rdtbt    ! ocean time-step index 
     282      REAL(wp)                , INTENT(in   ) ::   rDt_e    ! ocean time-step index 
    288283      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc   
    289284      ! 
     
    304299      zdepwd = 50._wp   ! maximum depth that ocean cells can have W/D processes 
    305300      ! 
    306       z2dt = rdtbt    
     301      z2dt = rDt_e    
    307302      ! 
    308303      zflxp(:,:)   = 0._wp 
     
    311306      zwdlmtv(:,:) = 1._wp 
    312307      ! 
    313       DO jj = 2, jpj      ! Horizontal Flux in u and v direction 
    314          DO ji = 2, jpi  
    315             ! 
    316             IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells 
    317             IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
    318             ! 
    319             zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   & 
    320                &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp )  
    321             zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   & 
    322                &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp )  
    323             ! 
    324             zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
    325             IF( zdep2 <= 0._wp ) THEN  !add more safety, but not necessary 
    326               sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    327               IF( zflxu(ji  ,jj  ) > 0._wp)   zwdlmtu(ji  ,jj  ) = 0._wp 
    328               IF( zflxu(ji-1,jj  ) < 0._wp)   zwdlmtu(ji-1,jj  ) = 0._wp 
    329               IF( zflxv(ji  ,jj  ) > 0._wp)   zwdlmtv(ji  ,jj  ) = 0._wp 
    330               IF( zflxv(ji  ,jj-1) < 0._wp)   zwdlmtv(ji  ,jj-1) = 0._wp  
    331             ENDIF 
    332          END DO 
    333       END DO 
     308      DO_2D_01_01 
     309         ! 
     310         IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells 
     311         IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
     312         ! 
     313         zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   & 
     314            &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp )  
     315         zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   & 
     316            &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp )  
     317         ! 
     318         zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
     319         IF( zdep2 <= 0._wp ) THEN  !add more safety, but not necessary 
     320           sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
     321           IF( zflxu(ji  ,jj  ) > 0._wp)   zwdlmtu(ji  ,jj  ) = 0._wp 
     322           IF( zflxu(ji-1,jj  ) < 0._wp)   zwdlmtu(ji-1,jj  ) = 0._wp 
     323           IF( zflxv(ji  ,jj  ) > 0._wp)   zwdlmtv(ji  ,jj  ) = 0._wp 
     324           IF( zflxv(ji  ,jj-1) < 0._wp)   zwdlmtv(ji  ,jj-1) = 0._wp  
     325         ENDIF 
     326      END_2D 
    334327      ! 
    335328      DO jk1 = 1, nn_wdit + 1      !! start limiter iterations  
     
    339332         jflag = 0     ! flag indicating if any further iterations are needed 
    340333         ! 
    341          DO jj = 2, jpj 
    342             DO ji = 2, jpi  
    343                ! 
    344                IF( tmask(ji, jj, 1 ) < 0.5_wp )   CYCLE  
    345                IF( ht_0(ji,jj)       > zdepwd )   CYCLE 
    346                ! 
    347                ztmp = e1e2t(ji,jj) 
    348                ! 
    349                zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp)   & 
    350                   &   + max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp)  
    351                zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp)   & 
    352                   &   + min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp)  
    353            
    354                zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    355                zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
    356            
    357                IF(zdep1 > zdep2) THEN 
    358                  zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt ) 
    359                  !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt ) 
    360                  ! flag if the limiter has been used but stop flagging if the only 
    361                  ! changes have zeroed the coefficient since further iterations will 
    362                  ! not change anything 
    363                  IF( zcoef > 0._wp ) THEN 
    364                     jflag = 1  
    365                  ELSE 
    366                     zcoef = 0._wp 
    367                  ENDIF 
    368                  IF(jk1 > nn_wdit) zcoef = 0._wp 
    369                  IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef 
    370                  IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef 
    371                  IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef 
    372                  IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef 
    373                END IF 
    374             END DO ! ji loop 
    375          END DO  ! jj loop 
     334         DO_2D_01_01 
     335            ! 
     336            IF( tmask(ji, jj, 1 ) < 0.5_wp )   CYCLE  
     337            IF( ht_0(ji,jj)       > zdepwd )   CYCLE 
     338            ! 
     339            ztmp = e1e2t(ji,jj) 
     340            ! 
     341            zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp)   & 
     342               &   + max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp)  
     343            zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp)   & 
     344               &   + min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp)  
     345        
     346            zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
     347            zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
     348        
     349            IF(zdep1 > zdep2) THEN 
     350              zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt ) 
     351              !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt ) 
     352              ! flag if the limiter has been used but stop flagging if the only 
     353              ! changes have zeroed the coefficient since further iterations will 
     354              ! not change anything 
     355              IF( zcoef > 0._wp ) THEN 
     356                 jflag = 1  
     357              ELSE 
     358                 zcoef = 0._wp 
     359              ENDIF 
     360              IF(jk1 > nn_wdit) zcoef = 0._wp 
     361              IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef 
     362              IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef 
     363              IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef 
     364              IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef 
     365            END IF 
     366         END_2D 
    376367         ! 
    377368         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. ) 
     
    392383      IF( jflag == 1 .AND. lwp )   WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!' 
    393384      ! 
    394       !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field) 
     385      !IF( ln_rnf      )   CALL sbc_rnf_div( hdiv )          ! runoffs (update hdiv field) 
    395386      ! 
    396387      IF( ln_timing )   CALL timing_stop('wad_lmt_bt')      ! 
Note: See TracChangeset for help on using the changeset viewer.