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 5239 for branches/UKMO/dev_r5021_nn_etau_revision/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 – NEMO

Ignore:
Timestamp:
2015-04-29T10:35:19+02:00 (9 years ago)
Author:
davestorkey
Message:

Commit changes to UKMO nn_etau revision branch (including clearing
SVN keywords).

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5021_nn_etau_revision/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90

    • Property svn:keywords deleted
    r4990 r5239  
    7676   INTEGER  ::   nn_htau   ! type of tke profile of penetration (=0/1) 
    7777   REAL(wp) ::   rn_efr    ! fraction of TKE surface value which penetrates in the ocean 
     78   REAL(wp) ::   rn_c      ! fraction of TKE added within the mixed layer by nn_etau 
    7879   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not 
    7980   REAL(wp) ::   rn_lc     ! coef to compute vertical velocity of Langmuir cells 
     
    8182   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values) 
    8283   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m] 
     84   REAL(wp) ::   rhtau                     ! coefficient to relate MLD to htau when nn_htau == 2 
    8385   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3) 
    8486   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3) 
    8587 
    8688   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   en             !: now turbulent kinetic energy   [m2/s2] 
     89   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_niw          !: TKE budget- near-inertial waves term 
    8790   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau           ! depth of tke penetration (nn_htau) 
     91   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   efr            ! surface boundary condition for nn_etau = 4 
    8892   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl          ! now mixing lenght of dissipation 
    8993   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   avt_k , avm_k  ! not enhanced Kz 
     
    115119#endif 
    116120         &      en    (jpi,jpj,jpk) , htau  (jpi,jpj)    , dissl(jpi,jpj,jpk) ,     &  
    117          &      avt_k (jpi,jpj,jpk) , avm_k (jpi,jpj,jpk),                          & 
    118          &      avmu_k(jpi,jpj,jpk) , avmv_k(jpi,jpj,jpk), STAT= zdf_tke_alloc      ) 
     121         &      avt_k (jpi,jpj,jpk) , avm_k (jpi,jpj,jpk), efr  (jpi,jpj)     ,     & 
     122         &      avmu_k(jpi,jpj,jpk) , avmv_k(jpi,jpj,jpk), e_niw(jpi,jpj,jpk) ,     & 
     123         &      STAT= zdf_tke_alloc      ) 
    119124         ! 
    120125      IF( lk_mpp             )   CALL mpp_sum ( zdf_tke_alloc ) 
     
    383388      END DO 
    384389 
     390      !                                 ! Save TKE prior to nn_etau addition   
     391      e_niw(:,:,:) = en(:,:,:)   
     392      !   
    385393      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    386394      !                            !  TKE due to surface and internal wave breaking 
    387395      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     396      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale 
     397         DO jj = 2, jpjm1 
     398            DO ji = fs_2, fs_jpim1   ! vector opt. 
     399               htau(ji,jj) = rhtau * hmlp(ji,jj) 
     400            END DO 
     401         END DO 
     402      ENDIF 
     403#if defined key_iomput 
     404      ! 
     405      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time) 
     406#endif 
     407      ! 
    388408      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction) 
    389409         DO jk = 2, jpkm1 
     
    420440            END DO 
    421441         END DO 
     442      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up) 
     443         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau  
     444            DO jj = 2, jpjm1 
     445               DO ji = fs_2, fs_jpim1   ! vector opt. 
     446                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     447               END DO 
     448            END DO 
     449         ENDIF 
     450         DO jk = 2, jpkm1 
     451            DO jj = 2, jpjm1 
     452               DO ji = fs_2, fs_jpim1   ! vector opt. 
     453                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   & 
     454                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk) 
     455               END DO 
     456            END DO 
     457         END DO 
    422458      ENDIF 
    423459      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged) 
     460      ! 
     461      DO jk = 2, jpkm1                             ! TKE budget: near-inertial waves term   
     462         DO jj = 2, jpjm1   
     463            DO ji = fs_2, fs_jpim1   ! vector opt.   
     464               e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk)   
     465            END DO   
     466         END DO   
     467      END DO   
     468      !   
     469      CALL lbc_lnk( e_niw, 'W', 1. )   
    424470      ! 
    425471      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer 
     
    683729         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   & 
    684730         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   & 
    685          &                 nn_etau , nn_htau  , rn_efr    
    686       !!---------------------------------------------------------------------- 
    687       ! 
     731         &                 nn_etau , nn_htau  , rn_efr , rn_c    
     732      !!---------------------------------------------------------------------- 
     733      ! 
     734      ! NB. Default values of namelist parameters should be set in the reference namelist 
     735      !     but setting this one here because we aren't set up to use the reference namelist 
     736      !     properly yet.  
     737      rn_c = 0.8_wp 
     738 
    688739      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy 
    689740      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901) 
     
    718769         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau 
    719770         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr 
     771         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c 
    720772         WRITE(numout,*) 
    721773         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri 
     
    728780      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' ) 
    729781      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' ) 
    730       IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' ) 
     782      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' ) 
    731783      IF( nn_etau == 3 .AND. .NOT. lk_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) 
    732784 
     
    736788      ENDIF 
    737789       
    738       IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln  
     790      IF( nn_etau == 2  .OR. ( nn_etau /= 0 .AND. nn_htau == 2 ) )   CALL zdf_mxl( nit000 - 1 )      ! Initialization of nmln and hmlp 
    739791 
    740792      !                               !* depth of penetration of surface tke 
    741793      IF( nn_etau /= 0 ) THEN       
     794         htau(:,:) = 0._wp 
    742795         SELECT CASE( nn_htau )             ! Choice of the depth of penetration 
    743796         CASE( 0 )                                 ! constant depth penetration (here 10 meters) 
     
    745798         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees 
    746799            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )             
     800         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer 
     801            rhtau = -1._wp / LOG( 1._wp - rn_c ) 
     802         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees 
     803            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   ) 
     804         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south 
     805            DO jj = 2, jpjm1 
     806               DO ji = fs_2, fs_jpim1   ! vector opt. 
     807                  IF( gphit(ji,jj) <= 0._wp ) THEN 
     808                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     809                  ELSE 
     810                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     811                  ENDIF 
     812               END DO 
     813            END DO 
     814         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south, 
     815            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south 
     816               DO ji = fs_2, fs_jpim1   ! vector opt. 
     817                  IF( gphit(ji,jj) <= -30._wp ) THEN 
     818                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   ) 
     819                  ELSE 
     820                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     821                  ENDIF 
     822               END DO 
     823            END DO 
    747824         END SELECT 
     825         ! 
     826         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau 
     827            DO jj = 2, jpjm1 
     828               DO ji = fs_2, fs_jpim1   ! vector opt. 
     829                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     830               END DO 
     831            END DO 
     832         ENDIF 
    748833      ENDIF 
    749834      !                               !* set vertical eddy coef. to the background value 
Note: See TracChangeset for help on using the changeset viewer.