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

Ignore:
Timestamp:
2016-04-21T18:15:17+02:00 (8 years ago)
Author:
davestorkey
Message:

Commit remaining changes to UKMO/r5518_GO6_package branch from following branches:
UKMO/dev_r5021_nn_etau_revision@6238
UKMO/dev_r5107_mld_zint@5534
UKMO/dev_r5107_eorca025_closea@6390
UKMO/restart_datestamp@5539
UKMO/icebergs_latent_heat@5821
UKMO/icebergs_restart_single_file_corrected@6480
UKMO/product_diagnostics@5971
UKMO/antarctic_partial_slip@5961
Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5961 cf. r5958 of /branches/UKMO/antarctic_partial_slip/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r6349 cf. r5962 of /branches/UKMO/product_diagnostics/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r6480 cf. r6479 of /branches/UKMO/icebergs_restart_single_file_corrected/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5986 cf. r5852 of /branches/UKMO/icebergs_restart_single_file/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5821 cf. r5808 of /branches/UKMO/icebergs_latent_heat/NEMOGCM@6490

File:
1 edited

Legend:

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

    r6487 r6491  
    8383   INTEGER  ::   nn_htau   ! type of tke profile of penetration (=0/1) 
    8484   REAL(wp) ::   rn_efr    ! fraction of TKE surface value which penetrates in the ocean 
     85   REAL(wp) ::   rn_c      ! fraction of TKE added within the mixed layer by nn_etau 
    8586   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not 
    8687   REAL(wp) ::   rn_lc     ! coef to compute vertical velocity of Langmuir cells 
     
    8889   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values) 
    8990   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m] 
     91   REAL(wp) ::   rhtau                     ! coefficient to relate MLD to htau when nn_htau == 2 
    9092   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3) 
    9193   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3) 
    9294 
    9395   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau           ! depth of tke penetration (nn_htau) 
     96   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_niw          !: TKE budget- near-inertial waves term 
     97   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   efr            ! surface boundary condition for nn_etau = 4 
    9498   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl          ! now mixing lenght of dissipation 
    9599#if defined key_c1d 
     
    114118      !!---------------------------------------------------------------------- 
    115119      ALLOCATE(                                                                    & 
     120         &      efr  (jpi,jpj)     , e_niw(jpi,jpj,jpk) ,                         &       
    116121#if defined key_c1d 
    117122         &      e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) ,                          & 
     
    420425      END DO 
    421426 
     427      !                                 ! Save TKE prior to nn_etau addition   
     428      e_niw(:,:,:) = en(:,:,:)   
     429      !   
    422430      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    423431      !                            !  TKE due to surface and internal wave breaking 
    424432      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     433      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale 
     434         DO jj = 2, jpjm1 
     435            DO ji = fs_2, fs_jpim1   ! vector opt. 
     436               htau(ji,jj) = rhtau * hmlp(ji,jj) 
     437            END DO 
     438         END DO 
     439      ENDIF 
     440#if defined key_iomput 
     441      ! 
     442      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time) 
     443#endif 
     444      ! 
    425445      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction) 
    426446         DO jk = 2, jpkm1 
     
    457477            END DO 
    458478         END DO 
     479      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up) 
     480         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau  
     481            DO jj = 2, jpjm1 
     482               DO ji = fs_2, fs_jpim1   ! vector opt. 
     483                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     484               END DO 
     485            END DO 
     486         ENDIF 
     487         DO jk = 2, jpkm1 
     488            DO jj = 2, jpjm1 
     489               DO ji = fs_2, fs_jpim1   ! vector opt. 
     490                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   & 
     491                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk) 
     492               END DO 
     493            END DO 
     494         END DO 
    459495      ENDIF 
    460496      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged) 
     497      ! 
     498      DO jk = 2, jpkm1                             ! TKE budget: near-inertial waves term   
     499         DO jj = 2, jpjm1   
     500            DO ji = fs_2, fs_jpim1   ! vector opt.   
     501               e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk)   
     502            END DO   
     503         END DO   
     504      END DO   
     505      !   
     506      CALL lbc_lnk( e_niw, 'W', 1. )   
    461507      ! 
    462508      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer 
     
    722768         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   & 
    723769         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   & 
    724          &                 nn_etau , nn_htau  , rn_efr    
    725       !!---------------------------------------------------------------------- 
    726       ! 
     770         &                 nn_etau , nn_htau  , rn_efr , rn_c    
     771      !!---------------------------------------------------------------------- 
     772 
    727773      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy 
    728774      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901) 
     
    757803         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau 
    758804         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr 
     805         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c 
    759806         WRITE(numout,*) 
    760807         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri 
     
    767814      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' ) 
    768815      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' ) 
    769       IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' ) 
     816      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' ) 
    770817      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) 
    771818 
     
    780827      ENDIF 
    781828 
     829      IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN 
     830          ierr = zdf_mxl_alloc() 
     831          nmln(:,:) = nlb10           ! Initialization of nmln 
     832      ENDIF 
     833 
    782834      !                               !* depth of penetration of surface tke 
    783835      IF( nn_etau /= 0 ) THEN       
     836         htau(:,:) = 0._wp 
    784837         SELECT CASE( nn_htau )             ! Choice of the depth of penetration 
    785838         CASE( 0 )                                 ! constant depth penetration (here 10 meters) 
     
    787840         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees 
    788841            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )             
     842         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer 
     843            rhtau = -1._wp / LOG( 1._wp - rn_c ) 
     844         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees 
     845            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   ) 
     846         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south 
     847            DO jj = 2, jpjm1 
     848               DO ji = fs_2, fs_jpim1   ! vector opt. 
     849                  IF( gphit(ji,jj) <= 0._wp ) THEN 
     850                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     851                  ELSE 
     852                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     853                  ENDIF 
     854               END DO 
     855            END DO 
     856         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south, 
     857            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south 
     858               DO ji = fs_2, fs_jpim1   ! vector opt. 
     859                  IF( gphit(ji,jj) <= -30._wp ) THEN 
     860                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   ) 
     861                  ELSE 
     862                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   ) 
     863                  ENDIF 
     864               END DO 
     865            END DO 
    789866         END SELECT 
     867         ! 
     868         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau 
     869            DO jj = 2, jpjm1 
     870               DO ji = fs_2, fs_jpim1   ! vector opt. 
     871                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 
     872               END DO 
     873            END DO 
     874         ENDIF 
    790875      ENDIF 
    791876      !                               !* set vertical eddy coef. to the background value 
Note: See TracChangeset for help on using the changeset viewer.