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 6497 for trunk/NEMOGCM/NEMO – NEMO

Changeset 6497 for trunk/NEMOGCM/NEMO


Ignore:
Timestamp:
2016-04-27T09:33:46+02:00 (8 years ago)
Author:
gm
Message:

#1720 - trunk: add Casimir tidal parameterization

Location:
trunk/NEMOGCM/NEMO/OPA_SRC/ZDF
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfddm.F90

    r6140 r6497  
    174174                  &                             +  0.15 * zrau(ji,jj)          * zmskd2(ji,jj)  ) 
    175175               ! add to the eddy viscosity coef. previously computed 
     176# if defined key_zdftmx_new 
     177               ! key_zdftmx_new: New internal wave-driven param: use avs value computed by zdftmx 
     178               avs (ji,jj,jk) = avs(ji,jj,jk) + zavfs + zavds 
     179# else 
    176180               avs (ji,jj,jk) = avt(ji,jj,jk) + zavfs + zavds 
     181# endif 
    177182               avt (ji,jj,jk) = avt(ji,jj,jk) + zavft + zavdt 
    178183               avm (ji,jj,jk) = avm(ji,jj,jk) + MAX( zavft + zavdt, zavfs + zavds ) 
  • trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90

    r6140 r6497  
    323323                  zwlc = zind * rn_lc * zus * SIN( rpi * gdepw_n(ji,jj,jk) / zhlc(ji,jj) ) 
    324324                  !                                           ! TKE Langmuir circulation source term 
    325                   en(ji,jj,jk) = en(ji,jj,jk) + rdt * (1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1) 
     325                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * (1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc )   & 
     326                     &                              / zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1) 
    326327               END DO 
    327328            END DO 
     
    375376            DO ji = fs_2, fs_jpim1   ! vector opt. 
    376377               zcof   = zfact1 * tmask(ji,jj,jk) 
     378# if defined key_zdftmx_new 
     379               ! key_zdftmx_new: New internal wave-driven param: set a minimum value for Kz on TKE (ensure numerical stability) 
     380               zzd_up = zcof * MAX( avm(ji,jj,jk+1) + avm(ji,jj,jk), 2.e-5_wp )   &  ! upper diagonal 
     381                  &          / (  e3t_n(ji,jj,jk  ) * e3w_n(ji,jj,jk  )  ) 
     382               zzd_lw = zcof * MAX( avm(ji,jj,jk) + avm(ji,jj,jk-1), 2.e-5_wp )   &  ! lower diagonal 
     383                  &          / (  e3t_n(ji,jj,jk-1) * e3w_n(ji,jj,jk  )  ) 
     384# else 
    377385               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal 
    378386                  &          / ( e3t_n(ji,jj,jk  ) * e3w_n(ji,jj,jk  ) ) 
    379387               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal 
    380388                  &          / ( e3t_n(ji,jj,jk-1) * e3w_n(ji,jj,jk  ) ) 
     389# endif 
    381390               !                                   ! shear prod. at w-point weightened by mask 
    382391               zesh2  =  ( z3du(ji-1,jj,jk) + z3du(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   & 
     
    732741      ! 
    733742      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number 
     743# if defined key_zdftmx_new 
     744      ! key_zdftmx_new: New internal wave-driven param: specified value of rn_emin & rmxl_min are used 
     745      rn_emin  = 1.e-10_wp 
     746      rmxl_min = 1.e-03_wp 
     747      IF(lwp) THEN                  ! Control print 
     748         WRITE(numout,*) 
     749         WRITE(numout,*) 'zdf_tke_init :  New tidal mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3 ' 
     750         WRITE(numout,*) '~~~~~~~~~~~~' 
     751      ENDIF 
     752# else 
    734753      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity 
     754# endif 
    735755      ! 
    736756      IF(lwp) THEN                    !* Control print 
  • trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftmx.F90

    r6140 r6497  
    541541   END SUBROUTINE zdf_tmx_init 
    542542 
     543#elif defined key_zdftmx_new 
     544   !!---------------------------------------------------------------------- 
     545   !!   'key_zdftmx_new'               Internal wave-driven vertical mixing 
     546   !!---------------------------------------------------------------------- 
     547   !!   zdf_tmx       : global     momentum & tracer Kz with wave induced Kz 
     548   !!   zdf_tmx_init  : global     momentum & tracer Kz with wave induced Kz 
     549   !!---------------------------------------------------------------------- 
     550   USE oce            ! ocean dynamics and tracers variables 
     551   USE dom_oce        ! ocean space and time domain variables 
     552   USE zdf_oce        ! ocean vertical physics variables 
     553   USE zdfddm         ! ocean vertical physics: double diffusive mixing 
     554   USE lbclnk         ! ocean lateral boundary conditions (or mpp link) 
     555   USE eosbn2         ! ocean equation of state 
     556   USE phycst         ! physical constants 
     557   USE prtctl         ! Print control 
     558   USE in_out_manager ! I/O manager 
     559   USE iom            ! I/O Manager 
     560   USE lib_mpp        ! MPP library 
     561   USE wrk_nemo       ! work arrays 
     562   USE timing         ! Timing 
     563   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)   
     564 
     565   IMPLICIT NONE 
     566   PRIVATE 
     567 
     568   PUBLIC   zdf_tmx         ! called in step module  
     569   PUBLIC   zdf_tmx_init    ! called in nemogcm module  
     570   PUBLIC   zdf_tmx_alloc   ! called in nemogcm module 
     571 
     572   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftmx = .TRUE.    !: wave-driven mixing flag 
     573 
     574   !                       !!* Namelist  namzdf_tmx : internal wave-driven mixing * 
     575   INTEGER  ::  nn_zpyc     ! pycnocline-intensified mixing energy proportional to N (=1) or N^2 (=2) 
     576   LOGICAL  ::  ln_mevar    ! variable (=T) or constant (=F) mixing efficiency 
     577   LOGICAL  ::  ln_tsdiff   ! account for differential T/S wave-driven mixing (=T) or not (=F) 
     578 
     579   REAL(wp) ::  r1_6 = 1._wp / 6._wp 
     580 
     581   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ebot_tmx     ! power available from high-mode wave breaking (W/m2) 
     582   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   epyc_tmx     ! power available from low-mode, pycnocline-intensified wave breaking (W/m2) 
     583   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ecri_tmx     ! power available from low-mode, critical slope wave breaking (W/m2) 
     584   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hbot_tmx     ! WKB decay scale for high-mode energy dissipation (m) 
     585   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hcri_tmx     ! decay scale for low-mode critical slope dissipation (m) 
     586   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   emix_tmx     ! local energy density available for mixing (W/kg) 
     587   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   bflx_tmx     ! buoyancy flux Kz * N^2 (W/kg) 
     588   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   pcmap_tmx    ! vertically integrated buoyancy flux (W/m2) 
     589   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   zav_ratio    ! S/T diffusivity ratio (only for ln_tsdiff=T) 
     590   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   zav_wave     ! Internal wave-induced diffusivity 
     591 
     592   !! * Substitutions 
     593#  include "zdfddm_substitute.h90" 
     594#  include "vectopt_loop_substitute.h90" 
     595   !!---------------------------------------------------------------------- 
     596   !! NEMO/OPA 4.0 , NEMO Consortium (2016) 
     597   !! $Id$ 
     598   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     599   !!---------------------------------------------------------------------- 
     600CONTAINS 
     601 
     602   INTEGER FUNCTION zdf_tmx_alloc() 
     603      !!---------------------------------------------------------------------- 
     604      !!                ***  FUNCTION zdf_tmx_alloc  *** 
     605      !!---------------------------------------------------------------------- 
     606      ALLOCATE(     ebot_tmx(jpi,jpj),  epyc_tmx(jpi,jpj),  ecri_tmx(jpi,jpj)    ,   & 
     607      &             hbot_tmx(jpi,jpj),  hcri_tmx(jpi,jpj),  emix_tmx(jpi,jpj,jpk),   & 
     608      &         bflx_tmx(jpi,jpj,jpk), pcmap_tmx(jpi,jpj), zav_ratio(jpi,jpj,jpk),   &  
     609      &         zav_wave(jpi,jpj,jpk), STAT=zdf_tmx_alloc     ) 
     610      ! 
     611      IF( lk_mpp             )   CALL mpp_sum ( zdf_tmx_alloc ) 
     612      IF( zdf_tmx_alloc /= 0 )   CALL ctl_warn('zdf_tmx_alloc: failed to allocate arrays') 
     613   END FUNCTION zdf_tmx_alloc 
     614 
     615 
     616   SUBROUTINE zdf_tmx( kt ) 
     617      !!---------------------------------------------------------------------- 
     618      !!                  ***  ROUTINE zdf_tmx  *** 
     619      !!                    
     620      !! ** Purpose :   add to the vertical mixing coefficients the effect of 
     621      !!              breaking internal waves. 
     622      !! 
     623      !! ** Method  : - internal wave-driven vertical mixing is given by: 
     624      !!                  Kz_wave = min(  100 cm2/s, f(  Reb = emix_tmx /( Nu * N^2 )  ) 
     625      !!              where emix_tmx is the 3D space distribution of the wave-breaking  
     626      !!              energy and Nu the molecular kinematic viscosity. 
     627      !!              The function f(Reb) is linear (constant mixing efficiency) 
     628      !!              if the namelist parameter ln_mevar = F and nonlinear if ln_mevar = T. 
     629      !! 
     630      !!              - Compute emix_tmx, the 3D power density that allows to compute 
     631      !!              Reb and therefrom the wave-induced vertical diffusivity. 
     632      !!              This is divided into three components: 
     633      !!                 1. Bottom-intensified low-mode dissipation at critical slopes 
     634      !!                     emix_tmx(z) = ( ecri_tmx / rau0 ) * EXP( -(H-z)/hcri_tmx ) 
     635      !!                                   / ( 1. - EXP( - H/hcri_tmx ) ) * hcri_tmx 
     636      !!              where hcri_tmx is the characteristic length scale of the bottom  
     637      !!              intensification, ecri_tmx a map of available power, and H the ocean depth. 
     638      !!                 2. Pycnocline-intensified low-mode dissipation 
     639      !!                     emix_tmx(z) = ( epyc_tmx / rau0 ) * ( sqrt(rn2(z))^nn_zpyc ) 
     640      !!                                   / SUM( sqrt(rn2(z))^nn_zpyc * e3w(z) ) 
     641      !!              where epyc_tmx is a map of available power, and nn_zpyc 
     642      !!              is the chosen stratification-dependence of the internal wave 
     643      !!              energy dissipation. 
     644      !!                 3. WKB-height dependent high mode dissipation 
     645      !!                     emix_tmx(z) = ( ebot_tmx / rau0 ) * rn2(z) * EXP(-z_wkb(z)/hbot_tmx) 
     646      !!                                   / SUM( rn2(z) * EXP(-z_wkb(z)/hbot_tmx) * e3w(z) ) 
     647      !!              where hbot_tmx is the characteristic length scale of the WKB bottom  
     648      !!              intensification, ebot_tmx is a map of available power, and z_wkb is the 
     649      !!              WKB-stretched height above bottom defined as 
     650      !!                    z_wkb(z) = H * SUM( sqrt(rn2(z'>=z)) * e3w(z'>=z) ) 
     651      !!                                 / SUM( sqrt(rn2(z'))    * e3w(z')    ) 
     652      !! 
     653      !!              - update the model vertical eddy viscosity and diffusivity:  
     654      !!                     avt  = avt  +    av_wave 
     655      !!                     avm  = avm  +    av_wave 
     656      !!                     avmu = avmu + mi(av_wave) 
     657      !!                     avmv = avmv + mj(av_wave) 
     658      !! 
     659      !!              - if namelist parameter ln_tsdiff = T, account for differential mixing: 
     660      !!                     avs  = avt  +    av_wave * diffusivity_ratio(Reb) 
     661      !! 
     662      !! ** Action  : - Define emix_tmx used to compute internal wave-induced mixing 
     663      !!              - avt, avs, avm, avmu, avmv increased by internal wave-driven mixing     
     664      !! 
     665      !! References :  de Lavergne et al. 2015, JPO; 2016, in prep. 
     666      !!---------------------------------------------------------------------- 
     667      INTEGER, INTENT(in) ::   kt   ! ocean time-step  
     668      ! 
     669      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     670      REAL(wp) ::   ztpc         ! scalar workspace 
     671      REAL(wp), DIMENSION(:,:)  , POINTER ::  zfact     ! Used for vertical structure 
     672      REAL(wp), DIMENSION(:,:)  , POINTER ::  zhdep     ! Ocean depth 
     673      REAL(wp), DIMENSION(:,:,:), POINTER ::  zwkb      ! WKB-stretched height above bottom 
     674      REAL(wp), DIMENSION(:,:,:), POINTER ::  zweight   ! Weight for high mode vertical distribution 
     675      REAL(wp), DIMENSION(:,:,:), POINTER ::  znu_t     ! Molecular kinematic viscosity (T grid) 
     676      REAL(wp), DIMENSION(:,:,:), POINTER ::  znu_w     ! Molecular kinematic viscosity (W grid) 
     677      REAL(wp), DIMENSION(:,:,:), POINTER ::  zReb      ! Turbulence intensity parameter 
     678      !!---------------------------------------------------------------------- 
     679      ! 
     680      IF( nn_timing == 1 )   CALL timing_start('zdf_tmx') 
     681      ! 
     682      CALL wrk_alloc( jpi,jpj,       zfact, zhdep ) 
     683      CALL wrk_alloc( jpi,jpj,jpk,   zwkb, zweight, znu_t, znu_w, zReb ) 
     684 
     685      !                          ! ----------------------------- ! 
     686      !                          !  Internal wave-driven mixing  !  (compute zav_wave) 
     687      !                          ! ----------------------------- ! 
     688      !                              
     689      !                        !* Critical slope mixing: distribute energy over the time-varying ocean depth, 
     690      !                                                 using an exponential decay from the seafloor. 
     691      DO jj = 1, jpj                ! part independent of the level 
     692         DO ji = 1, jpi 
     693            zhdep(ji,jj) = gdepw_0(ji,jj,mbkt(ji,jj)+1)       ! depth of the ocean 
     694            zfact(ji,jj) = rau0 * (  1._wp - EXP( -zhdep(ji,jj) / hcri_tmx(ji,jj) )  ) 
     695            IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = ecri_tmx(ji,jj) / zfact(ji,jj) 
     696         END DO 
     697      END DO 
     698 
     699      DO jk = 2, jpkm1              ! complete with the level-dependent part 
     700         emix_tmx(:,:,jk) = zfact(:,:) * (  EXP( ( gde3w_n(:,:,jk  ) - zhdep(:,:) ) / hcri_tmx(:,:) )                      & 
     701            &                             - EXP( ( gde3w_n(:,:,jk-1) - zhdep(:,:) ) / hcri_tmx(:,:) )  ) * wmask(:,:,jk)   & 
     702            &                          / ( gde3w_n(:,:,jk) - gde3w_n(:,:,jk-1) ) 
     703      END DO 
     704 
     705      !                        !* Pycnocline-intensified mixing: distribute energy over the time-varying  
     706      !                        !* ocean depth as proportional to sqrt(rn2)^nn_zpyc 
     707 
     708      SELECT CASE ( nn_zpyc ) 
     709 
     710      CASE ( 1 )               ! Dissipation scales as N (recommended) 
     711 
     712         zfact(:,:) = 0._wp 
     713         DO jk = 2, jpkm1              ! part independent of the level 
     714            zfact(:,:) = zfact(:,:) + e3w_n(:,:,jk) * SQRT(  MAX( 0._wp, rn2(:,:,jk) )  ) * wmask(:,:,jk) 
     715         END DO 
     716 
     717         DO jj = 1, jpj 
     718            DO ji = 1, jpi 
     719               IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = epyc_tmx(ji,jj) / ( rau0 * zfact(ji,jj) ) 
     720            END DO 
     721         END DO 
     722 
     723         DO jk = 2, jpkm1              ! complete with the level-dependent part 
     724            emix_tmx(:,:,jk) = emix_tmx(:,:,jk) + zfact(:,:) * SQRT(  MAX( 0._wp, rn2(:,:,jk) )  ) * wmask(:,:,jk) 
     725         END DO 
     726 
     727      CASE ( 2 )               ! Dissipation scales as N^2 
     728 
     729         zfact(:,:) = 0._wp 
     730         DO jk = 2, jpkm1              ! part independent of the level 
     731            zfact(:,:) = zfact(:,:) + e3w_n(:,:,jk) * MAX( 0._wp, rn2(:,:,jk) ) * wmask(:,:,jk) 
     732         END DO 
     733 
     734         DO jj= 1, jpj 
     735            DO ji = 1, jpi 
     736               IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = epyc_tmx(ji,jj) / ( rau0 * zfact(ji,jj) ) 
     737            END DO 
     738         END DO 
     739 
     740         DO jk = 2, jpkm1              ! complete with the level-dependent part 
     741            emix_tmx(:,:,jk) = emix_tmx(:,:,jk) + zfact(:,:) * MAX( 0._wp, rn2(:,:,jk) ) * wmask(:,:,jk) 
     742         END DO 
     743 
     744      END SELECT 
     745 
     746      !                        !* WKB-height dependent mixing: distribute energy over the time-varying  
     747      !                        !* ocean depth as proportional to rn2 * exp(-z_wkb/rn_hbot) 
     748       
     749      zwkb(:,:,:) = 0._wp 
     750      zfact(:,:) = 0._wp 
     751      DO jk = 2, jpkm1 
     752         zfact(:,:) = zfact(:,:) + e3w_n(:,:,jk) * SQRT(  MAX( 0._wp, rn2(:,:,jk) )  ) * wmask(:,:,jk) 
     753         zwkb(:,:,jk) = zfact(:,:) 
     754      END DO 
     755 
     756      DO jk = 2, jpkm1 
     757         DO jj = 1, jpj 
     758            DO ji = 1, jpi 
     759               IF( zfact(ji,jj) /= 0 )   zwkb(ji,jj,jk) = zhdep(ji,jj) * ( zfact(ji,jj) - zwkb(ji,jj,jk) )   & 
     760                                            &           * tmask(ji,jj,jk) / zfact(ji,jj) 
     761            END DO 
     762         END DO 
     763      END DO 
     764      zwkb(:,:,1) = zhdep(:,:) * tmask(:,:,1) 
     765 
     766      zweight(:,:,:) = 0._wp 
     767      DO jk = 2, jpkm1 
     768         zweight(:,:,jk) = MAX( 0._wp, rn2(:,:,jk) ) * hbot_tmx(:,:) * wmask(:,:,jk)                    & 
     769            &   * (  EXP( -zwkb(:,:,jk) / hbot_tmx(:,:) ) - EXP( -zwkb(:,:,jk-1) / hbot_tmx(:,:) )  ) 
     770      END DO 
     771 
     772      zfact(:,:) = 0._wp 
     773      DO jk = 2, jpkm1              ! part independent of the level 
     774         zfact(:,:) = zfact(:,:) + zweight(:,:,jk) 
     775      END DO 
     776 
     777      DO jj = 1, jpj 
     778         DO ji = 1, jpi 
     779            IF( zfact(ji,jj) /= 0 )   zfact(ji,jj) = ebot_tmx(ji,jj) / ( rau0 * zfact(ji,jj) ) 
     780         END DO 
     781      END DO 
     782 
     783      DO jk = 2, jpkm1              ! complete with the level-dependent part 
     784         emix_tmx(:,:,jk) = emix_tmx(:,:,jk) + zweight(:,:,jk) * zfact(:,:) * wmask(:,:,jk)   & 
     785            &                                / ( gde3w_n(:,:,jk) - gde3w_n(:,:,jk-1) ) 
     786      END DO 
     787 
     788 
     789      ! Calculate molecular kinematic viscosity 
     790      znu_t(:,:,:) = 1.e-4_wp * (  17.91_wp - 0.53810_wp * tsn(:,:,:,jp_tem) + 0.00694_wp * tsn(:,:,:,jp_tem) * tsn(:,:,:,jp_tem)  & 
     791         &                                  + 0.02305_wp * tsn(:,:,:,jp_sal)  ) * tmask(:,:,:) * r1_rau0 
     792      DO jk = 2, jpkm1 
     793         znu_w(:,:,jk) = 0.5_wp * ( znu_t(:,:,jk-1) + znu_t(:,:,jk) ) * wmask(:,:,jk) 
     794      END DO 
     795 
     796      ! Calculate turbulence intensity parameter Reb 
     797      DO jk = 2, jpkm1 
     798         zReb(:,:,jk) = emix_tmx(:,:,jk) / MAX( 1.e-20_wp, znu_w(:,:,jk) * rn2(:,:,jk) ) 
     799      END DO 
     800 
     801      ! Define internal wave-induced diffusivity 
     802      DO jk = 2, jpkm1 
     803         zav_wave(:,:,jk) = znu_w(:,:,jk) * zReb(:,:,jk) * r1_6   ! This corresponds to a constant mixing efficiency of 1/6 
     804      END DO 
     805 
     806      IF( ln_mevar ) THEN              ! Variable mixing efficiency case : modify zav_wave in the 
     807         DO jk = 2, jpkm1              ! energetic (Reb > 480) and buoyancy-controlled (Reb <10.224 ) regimes 
     808            DO jj = 1, jpj 
     809               DO ji = 1, jpi 
     810                  IF( zReb(ji,jj,jk) > 480.00_wp ) THEN 
     811                     zav_wave(ji,jj,jk) = 3.6515_wp * znu_w(ji,jj,jk) * SQRT( zReb(ji,jj,jk) ) 
     812                  ELSEIF( zReb(ji,jj,jk) < 10.224_wp ) THEN 
     813                     zav_wave(ji,jj,jk) = 0.052125_wp * znu_w(ji,jj,jk) * zReb(ji,jj,jk) * SQRT( zReb(ji,jj,jk) ) 
     814                  ENDIF 
     815               END DO 
     816            END DO 
     817         END DO 
     818      ENDIF 
     819 
     820      DO jk = 2, jpkm1                 ! Bound diffusivity by molecular value and 100 cm2/s 
     821         zav_wave(:,:,jk) = MIN(  MAX( 1.4e-7_wp, zav_wave(:,:,jk) ), 1.e-2_wp  ) * wmask(:,:,jk) 
     822      END DO 
     823 
     824      IF( kt == nit000 ) THEN        !* Control print at first time-step: diagnose the energy consumed by zav_wave 
     825         ztpc = 0._wp 
     826         DO jk = 2, jpkm1 
     827            DO jj = 1, jpj 
     828               DO ji = 1, jpi 
     829                  ztpc = ztpc + e3w_n(ji,jj,jk) * e1e2t(ji,jj)   & 
     830                     &         * MAX( 0._wp, rn2(ji,jj,jk) ) * zav_wave(ji,jj,jk) * wmask(ji,jj,jk) * tmask_i(ji,jj) 
     831               END DO 
     832            END DO 
     833         END DO 
     834         IF( lk_mpp )   CALL mpp_sum( ztpc ) 
     835         ztpc = rau0 * ztpc ! Global integral of rauo * Kz * N^2 = power contributing to mixing  
     836  
     837         IF(lwp) THEN 
     838            WRITE(numout,*) 
     839            WRITE(numout,*) 'zdf_tmx : Internal wave-driven mixing (tmx)' 
     840            WRITE(numout,*) '~~~~~~~ ' 
     841            WRITE(numout,*) 
     842            WRITE(numout,*) '      Total power consumption by av_wave: ztpc =  ', ztpc * 1.e-12_wp, 'TW' 
     843         ENDIF 
     844      ENDIF 
     845 
     846      !                          ! ----------------------- ! 
     847      !                          !   Update  mixing coefs  !                           
     848      !                          ! ----------------------- ! 
     849      !       
     850      IF( ln_tsdiff ) THEN          !* Option for differential mixing of salinity and temperature 
     851         DO jk = 2, jpkm1              ! Calculate S/T diffusivity ratio as a function of Reb 
     852            DO jj = 1, jpj 
     853               DO ji = 1, jpi 
     854                  zav_ratio(ji,jj,jk) = ( 0.505_wp + 0.495_wp *                                                                  & 
     855                      &   TANH(    0.92_wp * (   LOG10(  MAX( 1.e-20_wp, zReb(ji,jj,jk) * 5._wp * r1_6 )  ) - 0.60_wp   )    )   & 
     856                      &                 ) * wmask(ji,jj,jk) 
     857               END DO 
     858            END DO 
     859         END DO 
     860         CALL iom_put( "av_ratio", zav_ratio ) 
     861         DO jk = 2, jpkm1           !* update momentum & tracer diffusivity with wave-driven mixing 
     862            fsavs(:,:,jk) = avt(:,:,jk) + zav_wave(:,:,jk) * zav_ratio(:,:,jk) 
     863            avt  (:,:,jk) = avt(:,:,jk) + zav_wave(:,:,jk) 
     864            avm  (:,:,jk) = avm(:,:,jk) + zav_wave(:,:,jk) 
     865         END DO 
     866         ! 
     867      ELSE                          !* update momentum & tracer diffusivity with wave-driven mixing 
     868         DO jk = 2, jpkm1 
     869            fsavs(:,:,jk) = avt(:,:,jk) + zav_wave(:,:,jk) 
     870            avt  (:,:,jk) = avt(:,:,jk) + zav_wave(:,:,jk) 
     871            avm  (:,:,jk) = avm(:,:,jk) + zav_wave(:,:,jk) 
     872         END DO 
     873      ENDIF 
     874 
     875      DO jk = 2, jpkm1              !* update momentum diffusivity at wu and wv points 
     876         DO jj = 2, jpjm1 
     877            DO ji = fs_2, fs_jpim1  ! vector opt. 
     878               avmu(ji,jj,jk) = avmu(ji,jj,jk) + 0.5_wp * ( zav_wave(ji,jj,jk) + zav_wave(ji+1,jj  ,jk) ) * wumask(ji,jj,jk) 
     879               avmv(ji,jj,jk) = avmv(ji,jj,jk) + 0.5_wp * ( zav_wave(ji,jj,jk) + zav_wave(ji  ,jj+1,jk) ) * wvmask(ji,jj,jk) 
     880            END DO 
     881         END DO 
     882      END DO 
     883      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! lateral boundary condition 
     884 
     885      !                             !* output internal wave-driven mixing coefficient 
     886      CALL iom_put( "av_wave", zav_wave ) 
     887                                    !* output useful diagnostics: N^2, Kz * N^2 (bflx_tmx),  
     888                                    !  vertical integral of rau0 * Kz * N^2 (pcmap_tmx), energy density (emix_tmx) 
     889      IF( iom_use("bflx_tmx") .OR. iom_use("pcmap_tmx") ) THEN 
     890         bflx_tmx(:,:,:) = MAX( 0._wp, rn2(:,:,:) ) * zav_wave(:,:,:) 
     891         pcmap_tmx(:,:) = 0._wp 
     892         DO jk = 2, jpkm1 
     893            pcmap_tmx(:,:) = pcmap_tmx(:,:) + e3w_n(:,:,jk) * bflx_tmx(:,:,jk) * wmask(:,:,jk) 
     894         END DO 
     895         pcmap_tmx(:,:) = rau0 * pcmap_tmx(:,:) 
     896         CALL iom_put( "bflx_tmx", bflx_tmx ) 
     897         CALL iom_put( "pcmap_tmx", pcmap_tmx ) 
     898      ENDIF 
     899      CALL iom_put( "bn2", rn2 ) 
     900      CALL iom_put( "emix_tmx", emix_tmx ) 
     901       
     902      CALL wrk_dealloc( jpi,jpj,       zfact, zhdep ) 
     903      CALL wrk_dealloc( jpi,jpj,jpk,   zwkb, zweight, znu_t, znu_w, zReb ) 
     904 
     905      IF(ln_ctl)   CALL prt_ctl(tab3d_1=zav_wave , clinfo1=' tmx - av_wave: ', tab3d_2=avt, clinfo2=' avt: ', ovlap=1, kdim=jpk) 
     906      ! 
     907      IF( nn_timing == 1 )   CALL timing_stop('zdf_tmx') 
     908      ! 
     909   END SUBROUTINE zdf_tmx 
     910 
     911 
     912   SUBROUTINE zdf_tmx_init 
     913      !!---------------------------------------------------------------------- 
     914      !!                  ***  ROUTINE zdf_tmx_init  *** 
     915      !!                      
     916      !! ** Purpose :   Initialization of the wave-driven vertical mixing, reading 
     917      !!              of input power maps and decay length scales in netcdf files. 
     918      !! 
     919      !! ** Method  : - Read the namzdf_tmx namelist and check the parameters 
     920      !! 
     921      !!              - Read the input data in NetCDF files : 
     922      !!              power available from high-mode wave breaking (mixing_power_bot.nc) 
     923      !!              power available from pycnocline-intensified wave-breaking (mixing_power_pyc.nc) 
     924      !!              power available from critical slope wave-breaking (mixing_power_cri.nc) 
     925      !!              WKB decay scale for high-mode wave-breaking (decay_scale_bot.nc) 
     926      !!              decay scale for critical slope wave-breaking (decay_scale_cri.nc) 
     927      !! 
     928      !! ** input   : - Namlist namzdf_tmx 
     929      !!              - NetCDF files : mixing_power_bot.nc, mixing_power_pyc.nc, mixing_power_cri.nc, 
     930      !!              decay_scale_bot.nc decay_scale_cri.nc 
     931      !! 
     932      !! ** Action  : - Increase by 1 the nstop flag is setting problem encounter 
     933      !!              - Define ebot_tmx, epyc_tmx, ecri_tmx, hbot_tmx, hcri_tmx 
     934      !! 
     935      !! References : de Lavergne et al. 2015, JPO; 2016, in prep. 
     936      !!          
     937      !!---------------------------------------------------------------------- 
     938      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     939      INTEGER  ::   inum         ! local integer 
     940      INTEGER  ::   ios 
     941      REAL(wp) ::   zbot, zpyc, zcri   ! local scalars 
     942      !! 
     943      NAMELIST/namzdf_tmx_new/ nn_zpyc, ln_mevar, ln_tsdiff 
     944      !!---------------------------------------------------------------------- 
     945      ! 
     946      IF( nn_timing == 1 )  CALL timing_start('zdf_tmx_init') 
     947      ! 
     948      REWIND( numnam_ref )              ! Namelist namzdf_tmx in reference namelist : Wave-driven mixing 
     949      READ  ( numnam_ref, namzdf_tmx_new, IOSTAT = ios, ERR = 901) 
     950901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tmx in reference namelist', lwp ) 
     951      ! 
     952      REWIND( numnam_cfg )              ! Namelist namzdf_tmx in configuration namelist : Wave-driven mixing 
     953      READ  ( numnam_cfg, namzdf_tmx_new, IOSTAT = ios, ERR = 902 ) 
     954902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tmx in configuration namelist', lwp ) 
     955      IF(lwm) WRITE ( numond, namzdf_tmx_new ) 
     956      ! 
     957      IF(lwp) THEN                  ! Control print 
     958         WRITE(numout,*) 
     959         WRITE(numout,*) 'zdf_tmx_init : internal wave-driven mixing' 
     960         WRITE(numout,*) '~~~~~~~~~~~~' 
     961         WRITE(numout,*) '   Namelist namzdf_tmx_new : set wave-driven mixing parameters' 
     962         WRITE(numout,*) '      Pycnocline-intensified diss. scales as N (=1) or N^2 (=2) = ', nn_zpyc 
     963         WRITE(numout,*) '      Variable (T) or constant (F) mixing efficiency            = ', ln_mevar 
     964         WRITE(numout,*) '      Differential internal wave-driven mixing (T) or not (F)   = ', ln_tsdiff 
     965      ENDIF 
     966       
     967      ! The new wave-driven mixing parameterization elevates avt and avm in the interior, and 
     968      ! ensures that avt remains larger than its molecular value (=1.4e-7). Therefore, avtb should  
     969      ! be set here to a very small value, and avmb to its (uniform) molecular value (=1.4e-6). 
     970      avmb(:) = 1.4e-6_wp        ! viscous molecular value 
     971      avtb(:) = 1.e-10_wp        ! very small diffusive minimum (background avt is specified in zdf_tmx)     
     972      avtb_2d(:,:) = 1.e0_wp     ! uniform  
     973      IF(lwp) THEN                  ! Control print 
     974         WRITE(numout,*) 
     975         WRITE(numout,*) '   Force the background value applied to avm & avt in TKE to be everywhere ',   & 
     976            &               'the viscous molecular value & a very small diffusive value, resp.' 
     977      ENDIF 
     978       
     979      IF( .NOT.lk_zdfddm )   CALL ctl_stop( 'STOP', 'zdf_tmx_init_new : key_zdftmx_new requires key_zdfddm' ) 
     980       
     981      !                             ! allocate tmx arrays 
     982      IF( zdf_tmx_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tmx_init : unable to allocate tmx arrays' ) 
     983      ! 
     984      !                             ! read necessary fields 
     985      CALL iom_open('mixing_power_bot',inum)       ! energy flux for high-mode wave breaking [W/m2] 
     986      CALL iom_get  (inum, jpdom_data, 'field', ebot_tmx, 1 )  
     987      CALL iom_close(inum) 
     988      ! 
     989      CALL iom_open('mixing_power_pyc',inum)       ! energy flux for pynocline-intensified wave breaking [W/m2] 
     990      CALL iom_get  (inum, jpdom_data, 'field', epyc_tmx, 1 ) 
     991      CALL iom_close(inum) 
     992      ! 
     993      CALL iom_open('mixing_power_cri',inum)       ! energy flux for critical slope wave breaking [W/m2] 
     994      CALL iom_get  (inum, jpdom_data, 'field', ecri_tmx, 1 ) 
     995      CALL iom_close(inum) 
     996      ! 
     997      CALL iom_open('decay_scale_bot',inum)        ! spatially variable decay scale for high-mode wave breaking [m] 
     998      CALL iom_get  (inum, jpdom_data, 'field', hbot_tmx, 1 ) 
     999      CALL iom_close(inum) 
     1000      ! 
     1001      CALL iom_open('decay_scale_cri',inum)        ! spatially variable decay scale for critical slope wave breaking [m] 
     1002      CALL iom_get  (inum, jpdom_data, 'field', hcri_tmx, 1 ) 
     1003      CALL iom_close(inum) 
     1004 
     1005      ebot_tmx(:,:) = ebot_tmx(:,:) * ssmask(:,:) 
     1006      epyc_tmx(:,:) = epyc_tmx(:,:) * ssmask(:,:) 
     1007      ecri_tmx(:,:) = ecri_tmx(:,:) * ssmask(:,:) 
     1008 
     1009      ! Set once for all to zero the first and last vertical levels of appropriate variables 
     1010      emix_tmx (:,:, 1 ) = 0._wp 
     1011      emix_tmx (:,:,jpk) = 0._wp 
     1012      zav_ratio(:,:, 1 ) = 0._wp 
     1013      zav_ratio(:,:,jpk) = 0._wp 
     1014      zav_wave (:,:, 1 ) = 0._wp 
     1015      zav_wave (:,:,jpk) = 0._wp 
     1016 
     1017      zbot = glob_sum( e1e2t(:,:) * ebot_tmx(:,:) ) 
     1018      zpyc = glob_sum( e1e2t(:,:) * epyc_tmx(:,:) ) 
     1019      zcri = glob_sum( e1e2t(:,:) * ecri_tmx(:,:) ) 
     1020      IF(lwp) THEN 
     1021         WRITE(numout,*) '      High-mode wave-breaking energy:             ', zbot * 1.e-12_wp, 'TW' 
     1022         WRITE(numout,*) '      Pycnocline-intensifed wave-breaking energy: ', zpyc * 1.e-12_wp, 'TW' 
     1023         WRITE(numout,*) '      Critical slope wave-breaking energy:        ', zcri * 1.e-12_wp, 'TW' 
     1024      ENDIF 
     1025      ! 
     1026      IF( nn_timing == 1 )  CALL timing_stop('zdf_tmx_init') 
     1027      ! 
     1028   END SUBROUTINE zdf_tmx_init 
     1029 
    5431030#else 
    5441031   !!---------------------------------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.