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 8544 for branches/UKMO/ROMS_WAD_7832/NEMOGCM/CONFIG/CS15mini/MY_SRC/wet_dry.F90 – NEMO

Ignore:
Timestamp:
2017-09-19T17:41:15+02:00 (7 years ago)
Author:
deazer
Message:

Added reference level (above all potential wet points) to avoid negative depth bathymetry as a work around.
Require reference level to be emedded into the domain configuration file
uses ht_0 instead of ht_wd. This si still in the code but should be removed in time.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/ROMS_WAD_7832/NEMOGCM/CONFIG/CS15mini/MY_SRC/wet_dry.F90

    r8403 r8544  
    7373      !!---------------------------------------------------------------------- 
    7474      NAMELIST/namwad/ ln_wd, ln_rwd, rn_wdmin0, ln_rwd, rn_wdmin1, rn_wdmin2, rn_wdld, nn_wdit, ln_wd_diag, ln_rwd_bc, & 
    75                      & rn_ssh_ref, rn_ht_0, jn_wd_i, jn_wd_j, jn_wd_k,ln_rwd_rmp 
     75                     & rn_ht_0, jn_wd_i, jn_wd_j, jn_wd_k,ln_rwd_rmp 
    7676 
    7777      INTEGER  ::   ios                 ! Local integer output status for namelist read 
     
    104104         WRITE(numout,*) '      T => baroclinic u,v = 0 at dry pts: ln_rwd_bc = ', ln_rwd_bc      
    105105         WRITE(numout,*) '      use a ramp for rwd limiter:      ln_rwd_rmp   = ', ln_rwd_rmp 
    106          WRITE(numout,*) '      height of z = 0 above the geoid: rn_ssh_ref   = ', rn_ssh_ref 
    107106         WRITE(numout,*) '      the height (z) at which ht_0 = 0:rn_ht_0      = ', rn_ht_0   
    108107         WRITE(numout,*) '      i-index for diagnostic point     jn_wd_i      = ', jn_wd_i 
     
    111110      ENDIF 
    112111      ! 
    113 ! CEOD .OR. ln_rwd added to prevent segmentation errors  
    114112      IF(ln_wd .OR. ln_rwd) THEN 
    115113         ALLOCATE( wdmask(jpi,jpj), ht_wd(jpi,jpj),  STAT=ierr ) 
     
    189187 
    190188             IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE   ! we don't care about land cells 
    191              IF( ht_wd(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
     189             IF( ht_0(ji,jj)-rn_ssh_ref > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
    192190 
    193191              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + & 
     
    196194                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp)  
    197195 
    198               zdep2 = ht_wd(ji,jj) + sshb1(ji,jj) - rn_wdmin1 
     196              zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 
    199197              IF(zdep2 .le. 0._wp) THEN  !add more safty, but not necessary 
    200                 sshb1(ji,jj) = rn_wdmin1 - ht_wd(ji,jj) 
     198                sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    201199                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
    202200                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
     
    211209! slwa 
    212210! HPG limiter from jholt 
    213       wdramp(:,:) = min((ht_wd(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp) 
     211      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp) 
    214212!      write(6,*)'wdramp ',wdramp(10,10),wdramp(10,11) 
    215213!jth assume don't need a lbc_lnk here 
     
    238236         
    239237                 IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE  
    240                  IF( ht_wd(ji,jj) > zdepwd )      CYCLE 
     238                 IF( ht_0(ji,jj) > zdepwd )      CYCLE 
    241239         
    242240                 ztmp = e1e2t(ji,jj) 
     
    248246           
    249247                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    250                  zdep2 = ht_wd(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj) 
     248                 zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj) 
    251249           
    252250                 IF( zdep1 > zdep2 ) THEN 
     
    362360 
    363361             IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells 
    364              IF( ht_wd(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
     362             IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry 
    365363 
    366364              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + & 
     
    369367                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp)  
    370368 
    371               zdep2 = ht_wd(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
     369              zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 
    372370              IF(zdep2 .le. 0._wp) THEN  !add more safety, but not necessary 
    373                 sshn_e(ji,jj) = rn_wdmin1 - ht_wd(ji,jj) 
     371                sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj) 
    374372                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp 
    375373                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp 
     
    393391         
    394392                 IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE  
    395                  IF( ht_wd(ji,jj) > zdepwd )      CYCLE 
     393                 IF( ht_0(ji,jj) > zdepwd )      CYCLE 
    396394         
    397395                 ztmp = e1e2t(ji,jj) 
     
    403401           
    404402                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp 
    405                  zdep2 = ht_wd(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
     403                 zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj) 
    406404           
    407405                 IF(zdep1 > zdep2) THEN 
Note: See TracChangeset for help on using the changeset viewer.