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 10456 – NEMO

Changeset 10456


Ignore:
Timestamp:
2019-01-06T17:17:13+01:00 (5 years ago)
Author:
deazer
Message:

Added option to taper sbc fluxes near very shallow water when using WAD
Corrected some IO bugs in dia25h, diatmb for WAD case.
User has control of the tapering. At what depth to start it, and at what fraction to start
the tanh tapering. At the WAD limit SBC is turned off completely.
Dry cells do not have any communication with the atmosphere
To DO: Documentation update.
Although not all sette tests are passed (AGRIF etc.)
it does no worse than the trunk at the revision the branch is made

Location:
NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/cfgs/SHARED/namelist_ref

    r10428 r10456  
    116116   ln_wd_dl    = .false.   !  T/F activation of directional limiter 
    117117   ln_wd_dl_bc = .false.   !  T/F Directional limiteer Baroclinic option 
    118    ln_wd_dl_rmp = .false.   !  T/F Turn on directional limiter ramp 
     118   ln_wd_dl_rmp = .false.  !  T/F Turn on directional limiter ramp 
    119119   rn_wdmin0   =  0.30     !  depth at which WaD starts 
    120120   rn_wdmin1   =  0.2      !  Minimum wet depth on dried cells 
     
    122122   rn_wdld     =  2.5      !  Land elevation below which WaD is allowed 
    123123   nn_wdit     =   20      !  Max iterations for WaD limiter 
     124   rn_wd_sbcdep =  5.0    !  Depth at which to taper sbc fluxes 
     125   rn_wd_sbcfra =  0.999   !  Fraction of SBC fluxes at taper depth (Must be <1) 
     126    
    124127/ 
    125128!----------------------------------------------------------------------- 
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/DIA/dia25h.F90

    r10069 r10456  
    1313   USE in_out_manager  ! I/O units 
    1414   USE iom             ! I/0 library 
     15   USE wet_dry 
    1516 
    1617   IMPLICIT NONE 
     
    211212         CALL iom_put( "salin25h", zw3d  )   ! salinity 
    212213         zw2d(:,:) = sshn_25h(:,:)*tmask(:,:,1) + zmdi*(1.0-tmask(:,:,1)) 
    213          CALL iom_put( "ssh25h", zw2d )   ! sea surface  
     214         IF( ll_wd ) THEN 
     215            CALL iom_put( "ssh25h", zw2d+ssh_ref )   ! sea surface  
     216         ELSE 
     217            CALL iom_put( "ssh25h", zw2d )   ! sea surface 
     218         ENDIF 
    214219         ! Write velocities (instantaneous) 
    215220         zw3d(:,:,:) = un_25h(:,:,:)*umask(:,:,:) + zmdi*(1.0-umask(:,:,:)) 
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/DIA/diatmb.F90

    r10069 r10456  
    1212   USE in_out_manager  ! I/O units 
    1313   USE iom             ! I/0 library 
     14   USE wet_dry 
    1415 
    1516   IMPLICIT NONE 
     
    108109      CALL dia_calctmb( tsn(:,:,:,jp_tem), zwtmb ) 
    109110      !ssh already output but here we output it masked 
    110       CALL iom_put( "sshnmasked", sshn(:,:)*tmask(:,:,1) + zmdi*(1.0 - tmask(:,:,1)) ) 
     111      IF( ll_wd ) THEN 
     112         CALL iom_put( "sshnmasked", (sshn(:,:)+ssh_ref)*tmask(:,:,1) + zmdi*(1.0 - tmask(:,:,1)) ) 
     113      ELSE 
     114         CALL iom_put( "sshnmasked", sshn(:,:)*tmask(:,:,1) + zmdi*(1.0 - tmask(:,:,1)) ) 
     115      ENDIF 
     116 
    111117      CALL iom_put( "top_temp"  , zwtmb(:,:,1) )    ! tmb Temperature 
    112118      CALL iom_put( "mid_temp"  , zwtmb(:,:,2) )    ! tmb Temperature 
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/DOM/istate.F90

    r10068 r10456  
    100100            ! 
    101101            sshb(:,:)   = 0._wp               ! set the ocean at rest 
     102            IF(ll_wd) then 
     103               sshb(:,:) =  -ssh_ref  ! Added in 30 here for bathy that adds 30 as Iterative test CEOD  
     104               ! 
     105               ! Apply minimum wetdepth criterion 
     106               ! 
     107               DO jj = 1,jpj 
     108                  DO ji = 1,jpi 
     109                     IF( ht_0(ji,jj) + sshb(ji,jj)  < rn_wdmin1 ) THEN 
     110                        sshb(ji,jj) = tmask(ji,jj,1)*( rn_wdmin1 - (ht_0(ji,jj)) ) 
     111                     ENDIF 
     112                  END DO 
     113               END DO  
     114            ENDIF  
    102115            ub  (:,:,:) = 0._wp 
    103116            vb  (:,:,:) = 0._wp   
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/DYN/wet_dry.F90

    r10425 r10456  
    4545   REAL(wp), PUBLIC  ::   r_rn_wdmin1 !: 1/minimum water depth on dried cells  
    4646   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells 
     47   REAL(wp), PUBLIC  ::   rn_wd_sbcdep   !: Depth at which to taper sbc fluxes 
     48   REAL(wp), PUBLIC  ::   rn_wd_sbcfra   !: Fraction of SBC at taper depth 
    4749   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying will be considered 
    4850   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter 
     
    7476      !! 
    7577      NAMELIST/namwad/ ln_wd_il, ln_wd_dl   , rn_wdmin0, rn_wdmin1, rn_wdmin2, rn_wdld,   & 
    76          &             nn_wdit , ln_wd_dl_bc, ln_wd_dl_rmp 
     78         &             nn_wdit , ln_wd_dl_bc, ln_wd_dl_rmp, rn_wd_sbcdep,rn_wd_sbcfra 
    7779      !!---------------------------------------------------------------------- 
    7880      ! 
     
    8587      IF(lwm) WRITE ( numond, namwad ) 
    8688      ! 
     89      IF( rn_wd_sbcfra>=1 )   CALL ctl_stop( 'STOP', 'rn_wd_sbcfra >=1 : must be < 1' ) 
    8790      IF(lwp) THEN                  ! control print 
    8891         WRITE(numout,*) 
     
    99102         WRITE(numout,*) '      T => baroclinic u,v=0 at dry pts: ln_wd_dl_bc = ', ln_wd_dl_bc      
    100103         WRITE(numout,*) '      use a ramp for rwd limiter:  ln_wd_dl_rwd_rmp = ', ln_wd_dl_rmp 
     104         WRITE(numout,*) '      cut off depth sbc for wd   rn_wd_sbcdep       = ', rn_wd_sbcdep 
     105         WRITE(numout,*) '      fraction to start sbc wgt rn_wd_sbcfra        = ', rn_wd_sbcfra 
    101106      ENDIF 
    102107      IF( .NOT. ln_read_cfg ) THEN 
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/SBC/sbcmod.F90

    r10425 r10456  
    5858   USE lib_mpp        ! MPP library 
    5959   USE timing         ! Timing 
     60   USE wet_dry 
    6061   USE diurnal_bulk, ONLY:   ln_diurnal_only   ! diurnal SST diagnostic 
    6162 
     
    377378      ! 
    378379      LOGICAL ::   ll_sas, ll_opa   ! local logical 
     380      ! 
     381      REAL(wp) ::     zthscl        ! wd  tanh scale 
     382      REAL(wp), DIMENSION(jpi,jpj) ::  zwdht, zwght  ! wd dep over wd limit, wgt   
     383 
    379384      !!--------------------------------------------------------------------- 
    380385      ! 
     
    461466!!$!clem: it looks like it is necessary for the north fold (in certain circumstances). Don't know why. 
    462467!!$      CALL lbc_lnk( 'sbcmod', emp, 'T', 1. ) 
     468   IF ( ll_wd ) THEN     ! If near WAD point limit the flux for now 
     469    zthscl = atanh(rn_wd_sbcfra)                     ! taper frac default is .999  
     470    zwdht(:,:) = sshn(:,:) + ht_0(:,:) - rn_wdmin1   ! do this calc of water 
     471                                                     ! depth above wd limit once 
     472    WHERE( zwdht(:,:) <= 0.0 ) 
     473            taum(:,:) = 0.0 
     474            utau(:,:) = 0.0 
     475            vtau(:,:) = 0.0 
     476            qns (:,:) = 0.0 
     477            qsr (:,:) = 0.0 
     478            emp (:,:) = min(emp(:,:),0.0) !can allow puddles to grow but not shrink 
     479            sfx (:,:) = 0.0 
     480    END WHERE 
     481    zwght(:,:) = tanh(zthscl*zwdht(:,:)) 
     482    WHERE( zwdht(:,:) > 0.0  .and. zwdht(:,:) < rn_wd_sbcdep ) !  5 m hard limit here is arbitrary 
     483            qsr  (:,:) =  qsr(:,:)  * zwght(:,:) 
     484            qns  (:,:) =  qns(:,:)  * zwght(:,:) 
     485            taum (:,:) =  taum(:,:) * zwght(:,:) 
     486            utau (:,:) =  utau(:,:) * zwght(:,:) 
     487            vtau (:,:) =  vtau(:,:) * zwght(:,:) 
     488            sfx (:,:)  =  sfx(:,:)  * zwght(:,:) 
     489            emp  (:,:) =  emp(:,:)  * zwght(:,:) 
     490    END WHERE 
     491   ENDIF 
    463492      ! 
    464493      IF( kt == nit000 ) THEN                          !   set the forcing field at nit000 - 1    ! 
  • NEMO/branches/UKMO/dev_10448_WAD_SBC_BUGFIX/src/OCE/TRA/trasbc.F90

    r10068 r10456  
    2727   USE trd_oce        ! trends: ocean variables 
    2828   USE trdtra         ! trends manager: tracers  
    29    USE wet_dry,  ONLY : ll_wd, rn_wdmin1, r_rn_wdmin1   ! Wetting and drying 
    3029#if defined key_asminc    
    3130   USE asminc         ! Assimilation increment 
     
    125124      DO jj = 2, jpj 
    126125         DO ji = fs_2, fs_jpim1   ! vector opt. 
    127             IF ( ll_wd ) THEN     ! If near WAD point limit the flux for now 
    128                IF ( sshn(ji,jj) + ht_0(ji,jj) >  2._wp * rn_wdmin1 ) THEN 
    129                   sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)   ! non solar heat flux 
    130                ELSE IF ( sshn(ji,jj) + ht_0(ji,jj) >  rn_wdmin1 ) THEN 
    131                   sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj) & 
    132                        &                * tanh ( 5._wp * ( ( sshn(ji,jj) + ht_0(ji,jj) -  rn_wdmin1 ) * r_rn_wdmin1 ) ) 
    133                ELSE 
    134                   sbc_tsc(ji,jj,jp_tem) = 0._wp 
    135                ENDIF 
    136             ELSE  
    137                sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)   ! non solar heat flux 
    138             ENDIF 
    139  
     126            sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)   ! non solar heat flux 
    140127            sbc_tsc(ji,jj,jp_sal) = r1_rau0     * sfx(ji,jj)   ! salt flux due to freezing/melting 
    141128         END DO 
Note: See TracChangeset for help on using the changeset viewer.