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

Ignore:
Timestamp:
2016-05-03T14:28:12+02:00 (8 years ago)
Author:
timgraham
Message:

First attempt at merging in science changes from GO6 package branch at v3.6 stable (Note-namelists not yet dealt with)

File:
1 edited

Legend:

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

    r6503 r6507  
    8181   INTEGER  ::   nn_htau   ! type of tke profile of penetration (=0/1) 
    8282   REAL(wp) ::   rn_efr    ! fraction of TKE surface value which penetrates in the ocean 
     83   REAL(wp) ::   rn_c      ! fraction of TKE added within the mixed layer by nn_etau 
    8384   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not 
    8485   REAL(wp) ::   rn_lc     ! coef to compute vertical velocity of Langmuir cells 
     
    8687   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values) 
    8788   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m] 
     89   REAL(wp) ::   rhtau                     ! coefficient to relate MLD to htau when nn_htau == 2 
    8890   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3) 
    8991   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3) 
    9092 
    9193   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau           ! depth of tke penetration (nn_htau) 
     94   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   efr            ! surface boundary condition for nn_etau = 4 
    9295   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl          ! now mixing lenght of dissipation 
    9396   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr          ! now mixing lenght of dissipation 
     
    112115      !!---------------------------------------------------------------------- 
    113116      ALLOCATE(                                                                    & 
     117         &      efr  (jpi,jpj)     ,                                               &       
    114118#if defined key_c1d 
    115119         &      e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) ,                          & 
     
    446450      !                            !  TKE due to surface and internal wave breaking 
    447451      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     452      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale 
     453         DO jj = 2, jpjm1 
     454            DO ji = fs_2, fs_jpim1   ! vector opt. 
     455               htau(ji,jj) = rhtau * hmlp(ji,jj) 
     456            END DO 
     457         END DO 
     458      ENDIF 
     459#if defined key_iomput 
     460      ! 
     461      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time) 
     462#endif 
     463      ! 
    448464!!gm BUG : in the exp  remove the depth of ssh !!! 
    449        
    450465       
    451466      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction) 
     
    477492                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw_n(ji,jj,jk) / htau(ji,jj) )   & 
    478493                     &                        * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) 
     494               END DO 
     495            END DO 
     496         END DO 
     497      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up) 
     498         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau  
     499            DO jj = 2, jpjm1 
     500               DO ji = fs_2, fs_jpim1   ! vector opt. 
     501                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     502               END DO 
     503            END DO 
     504         ENDIF 
     505         DO jk = 2, jpkm1 
     506            DO jj = 2, jpjm1 
     507               DO ji = fs_2, fs_jpim1   ! vector opt. 
     508                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -gdepw_n(ji,jj,jk) / htau(ji,jj) )   & 
     509                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk) 
    479510               END DO 
    480511            END DO 
     
    723754      !!---------------------------------------------------------------------- 
    724755      INTEGER ::   ji, jj, jk   ! dummy loop indices 
    725       INTEGER ::   ios 
     756      INTEGER ::   ios, ierr 
    726757      !! 
    727758      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   & 
    728759         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   & 
    729760         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   & 
    730          &                 nn_etau , nn_htau  , rn_efr    
     761         &                 nn_etau , nn_htau  , rn_efr , rn_c    
    731762      !!---------------------------------------------------------------------- 
    732763      ! 
     
    774805         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau 
    775806         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr 
     807         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c 
    776808         WRITE(numout,*) 
    777809         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri 
     
    784816      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' ) 
    785817      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' ) 
    786       IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' ) 
     818      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' ) 
    787819      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) 
    788820 
     
    792824      ENDIF 
    793825       
    794       IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln  
     826      IF( nn_etau == 2  ) THEN 
     827          ierr = zdf_mxl_alloc() 
     828          nmln(:,:) = nlb10           ! Initialization of nmln 
     829      ENDIF 
     830 
     831      IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN 
     832          ierr = zdf_mxl_alloc() 
     833          nmln(:,:) = nlb10           ! Initialization of nmln 
     834      ENDIF 
    795835 
    796836      !                               !* depth of penetration of surface tke 
    797837      IF( nn_etau /= 0 ) THEN       
     838         htau(:,:) = 0._wp 
    798839         SELECT CASE( nn_htau )             ! Choice of the depth of penetration 
    799840         CASE( 0 )                                 ! constant depth penetration (here 10 meters) 
     
    801842         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees 
    802843            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )             
     844         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer 
     845            rhtau = -1._wp / LOG( 1._wp - rn_c ) 
     846         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees 
     847            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   ) 
     848         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south 
     849            DO jj = 2, jpjm1 
     850               DO ji = fs_2, fs_jpim1   ! vector opt. 
     851                  IF( gphit(ji,jj) <= 0._wp ) THEN 
     852                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     853                  ELSE 
     854                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     855                  ENDIF 
     856               END DO 
     857            END DO 
     858         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south, 
     859            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south 
     860               DO ji = fs_2, fs_jpim1   ! vector opt. 
     861                  IF( gphit(ji,jj) <= -30._wp ) THEN 
     862                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   ) 
     863                  ELSE 
     864                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     865                  ENDIF 
     866               END DO 
     867            END DO 
    803868         END SELECT 
     869         ! 
     870         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau 
     871            DO jj = 2, jpjm1 
     872               DO ji = fs_2, fs_jpim1   ! vector opt. 
     873                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     874               END DO 
     875            END DO 
     876         ENDIF 
    804877      ENDIF 
    805878      !                               !* set vertical eddy coef. to the background value 
Note: See TracChangeset for help on using the changeset viewer.