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 8865 for branches/UKMO/ROMS_WAD_7832/NEMOGCM/NEMO/OPA_SRC/DYN/wet_dry.F90 – NEMO

Ignore:
Timestamp:
2017-12-01T05:41:32+01:00 (6 years ago)
Author:
deazer
Message:

Moving Changes from CS15mini config into NEMO main src
Updating TEST configs to run wit this version of the code
all sette tests pass at this revision other than AGRIF
Includes updates to dynnxt and tranxt required for 3D rives in WAD case to be conservative.

Next commit will update naming conventions and tidy the code.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/ROMS_WAD_7832/NEMOGCM/NEMO/OPA_SRC/DYN/wet_dry.F90

    r7646 r8865  
    11MODULE wet_dry 
     2 
     3   !! includes updates to namelist namwad for diagnostic outputs of ROMS wetting and drying 
     4 
    25   !!============================================================================== 
    36   !!                       ***  MODULE  wet_dry  *** 
     
    3235 
    3336   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdmask         !: u- and v- limiter  
    34    REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   ht_wd          !: wetting and drying t-pt depths 
    3537                                                                     !  (can include negative depths) 
     38   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdramp, wdrampu, wdrampv ! for hpg limiting 
    3639 
    3740   LOGICAL,  PUBLIC  ::   ln_wd       !: Wetting/drying activation switch (T:on,F:off) 
     41   REAL(wp), PUBLIC  ::   rn_wdmin0   !: depth at which wetting/drying starts 
     42   LOGICAL,  PUBLIC  ::   ln_rwd      !: ROMS Wetting/drying activation switch (T:on,F:off) 
    3843   REAL(wp), PUBLIC  ::   rn_wdmin1   !: minimum water depth on dried cells 
    39    REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerrance of minimum water depth on dried cells 
     44   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells 
    4045   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying  
    4146                                           !: will be considered 
    4247   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter 
     48   LOGICAL,  PUBLIC  ::   ln_rwd_bc   !: ROMS scheme: True implies 3D velocities are set to the barotropic values at points  
     49                                      !: where the flow is from wet points on less than half the barotropic sub-steps   
     50   LOGICAL,  PUBLIC  ::   ln_rwd_rmp  !: use a ramp for the rwd flux limiter between 2 rn_wdmin1 and rn_wdmin1 (rather than a cut-off at rn_wdmin1)       
     51   REAL(wp), PUBLIC  ::   rn_ssh_ref  !: height of z=0 with respect to the geoid;  
     52   REAL(wp), PUBLIC  ::   rn_ht_0     !: the height at which ht_0 = 0    
     53 
     54   LOGICAL,  PUBLIC  ::   ln_wd_diag  ! True implies wad diagnostic at chosen point are printed out 
     55   INTEGER , PUBLIC  ::   jn_wd_i, jn_wd_j, jn_wd_k   ! indices at which diagnostic outputs are generated  
    4356 
    4457   PUBLIC   wad_init                  ! initialisation routine called by step.F90 
     
    5871      !! ** input   : - namwad namelist 
    5972      !!---------------------------------------------------------------------- 
    60       NAMELIST/namwad/ ln_wd, rn_wdmin1, rn_wdmin2, rn_wdld, nn_wdit 
     73      NAMELIST/namwad/ ln_wd, ln_rwd, rn_wdmin0, ln_rwd, rn_wdmin1, rn_wdmin2, rn_wdld, nn_wdit, ln_wd_diag, ln_rwd_bc, & 
     74                     &  rn_ht_0, jn_wd_i, jn_wd_j, jn_wd_k,ln_rwd_rmp 
     75 
    6176      INTEGER  ::   ios                 ! Local integer output status for namelist read 
    6277      INTEGER  ::   ierr                ! Local integer status array allocation  
     
    7893         WRITE(numout,*) '~~~~~~~~' 
    7994         WRITE(numout,*) '   Namelist namwad' 
    80          WRITE(numout,*) '      Logical activation                 ln_wd      = ', ln_wd 
     95         WRITE(numout,*) '      Logical for NOC  wd scheme         ln_wd      = ', ln_wd 
     96         WRITE(numout,*) '      Logical for ROMS wd scheme         ln_rwd     = ', ln_rwd 
     97         WRITE(numout,*) '      Depth at which wet/drying starts rn_wdmin0    = ', rn_wdmin0 
    8198         WRITE(numout,*) '      Minimum wet depth on dried cells rn_wdmin1    = ', rn_wdmin1 
    8299         WRITE(numout,*) '      Tolerance of min wet depth     rn_wdmin2      = ', rn_wdmin2 
    83100         WRITE(numout,*) '      land elevation threshold         rn_wdld      = ', rn_wdld 
    84101         WRITE(numout,*) '      Max iteration for W/D limiter    nn_wdit      = ', nn_wdit 
     102         WRITE(numout,*) '      Logical for WAD diagnostics      ln_wd_diag   = ', ln_wd_diag 
     103         WRITE(numout,*) '      T => baroclinic u,v = 0 at dry pts: ln_rwd_bc = ', ln_rwd_bc      
     104         WRITE(numout,*) '      use a ramp for rwd limiter:      ln_rwd_rmp   = ', ln_rwd_rmp 
     105         WRITE(numout,*) '      the height (z) at which ht_0 = 0:rn_ht_0      = ', rn_ht_0   
     106         WRITE(numout,*) '      i-index for diagnostic point     jn_wd_i      = ', jn_wd_i 
     107         WRITE(numout,*) '      j-index for diagnostic point     jn_wd_j      = ', jn_wd_j 
     108         WRITE(numout,*) '      k-index for diagnostic point     jn_wd_k      = ', jn_wd_k 
    85109      ENDIF 
    86110      ! 
    87       IF(ln_wd) THEN 
    88          ALLOCATE( wdmask(jpi,jpj), ht_wd(jpi,jpj),  STAT=ierr ) 
     111      IF(ln_wd .OR. ln_rwd) THEN 
     112         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr ) 
     113         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr )  
    89114         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error') 
    90115      ENDIF 
     
    161186 
    162187             IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE   ! we don't care about land cells 
    163              IF( ht_wd(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
     188             IF( ht_0(ji,jj)-rn_ssh_ref > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
    164189 
    165190              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + & 
     
    168193                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp)  
    169194 
    170               zdep2 = ht_wd(ji,jj) + sshb1(ji,jj) - rn_wdmin1 
     195              zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 
    171196              IF(zdep2 .le. 0._wp) THEN  !add more safty, but not necessary 
    172                 sshb1(ji,jj) = rn_wdmin1 - ht_wd(ji,jj) 
     197                sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    173198                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
    174199                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
     
    180205        END DO 
    181206 
     207 
     208! slwa 
     209! HPG limiter from jholt 
     210      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp) 
     211!      write(6,*)'wdramp ',wdramp(10,10),wdramp(10,11) 
     212!jth assume don't need a lbc_lnk here 
     213       DO jj = 1, jpjm1 
     214         DO ji = 1, jpim1 
     215           wdrampu(ji,jj) = min(wdramp(ji,jj),wdramp(ji+1,jj)) 
     216           wdrampv(ji,jj) = min(wdramp(ji,jj),wdramp(ji,jj+1)) 
     217         ENDDO 
     218       ENDDO 
     219        !wdrampu(:,:)=1.0_wp 
     220        !wdrampv(:,:)=1.0_wp 
     221! end HPG limiter 
     222 
     223 
    182224       
    183225        !! start limiter iterations  
     
    193235         
    194236                 IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE  
    195                  IF( ht_wd(ji,jj) > zdepwd )      CYCLE 
     237                 IF( ht_0(ji,jj) > zdepwd )      CYCLE 
    196238         
    197239                 ztmp = e1e2t(ji,jj) 
     
    203245           
    204246                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    205                  zdep2 = ht_wd(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj) 
     247                 zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj) 
    206248           
    207249                 IF( zdep1 > zdep2 ) THEN 
     
    317359 
    318360             IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells 
    319              IF( ht_wd(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
     361             IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
    320362 
    321363              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + & 
     
    324366                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp)  
    325367 
    326               zdep2 = ht_wd(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
     368              zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
    327369              IF(zdep2 .le. 0._wp) THEN  !add more safety, but not necessary 
    328                 sshn_e(ji,jj) = rn_wdmin1 - ht_wd(ji,jj) 
     370                sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    329371                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
    330372                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
     
    348390         
    349391                 IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE  
    350                  IF( ht_wd(ji,jj) > zdepwd )      CYCLE 
     392                 IF( ht_0(ji,jj) > zdepwd )      CYCLE 
    351393         
    352394                 ztmp = e1e2t(ji,jj) 
     
    358400           
    359401                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    360                  zdep2 = ht_wd(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
     402                 zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
    361403           
    362404                 IF(zdep1 > zdep2) THEN 
Note: See TracChangeset for help on using the changeset viewer.