Changeset 9987 for branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90
- Timestamp:
- 2018-07-23T11:33:03+02:00 (6 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90
r9188 r9987 53 53 USE timing ! Timing 54 54 USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 55 #if defined key_agrif 56 USE agrif_opa_interp 57 USE agrif_opa_update 58 #endif 59 60 55 61 56 62 IMPLICIT NONE … … 77 83 INTEGER :: nn_htau ! type of tke profile of penetration (=0/1) 78 84 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 79 86 LOGICAL :: ln_lc ! Langmuir cells (LC) as a source term of TKE or not 80 87 REAL(wp) :: rn_lc ! coef to compute vertical velocity of Langmuir cells … … 82 89 REAL(wp) :: ri_cri ! critic Richardson number (deduced from rn_ediff and rn_ediss values) 83 90 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 84 92 REAL(wp) :: rhftau_add = 1.e-3_wp ! add offset applied to HF part of taum (nn_etau=3) 85 93 REAL(wp) :: rhftau_scl = 1.0_wp ! scale factor applied to HF part of taum (nn_etau=3) 86 94 87 95 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 88 98 REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: dissl ! now mixing lenght of dissipation 89 99 #if defined key_c1d … … 108 118 !!---------------------------------------------------------------------- 109 119 ALLOCATE( & 120 & efr (jpi,jpj) , e_niw(jpi,jpj,jpk) , & 110 121 #if defined key_c1d 111 122 & e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) , & … … 184 195 avmv_k(:,:,:) = avmv(:,:,:) 185 196 ! 197 #if defined key_agrif 198 ! Update child grid f => parent grid 199 IF( .NOT.Agrif_Root() ) CALL Agrif_Update_Tke( kt ) ! children only 200 #endif 201 ! 186 202 END SUBROUTINE zdf_tke 187 203 … … 312 328 zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) ) 313 329 ! ! TKE Langmuir circulation source term 314 en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1) 330 en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( 1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) / & 331 & zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1) 315 332 END DO 316 333 END DO … … 345 362 DO ji = fs_2, fs_jpim1 ! vector opt. 346 363 zcof = zfact1 * tmask(ji,jj,jk) 364 # if defined key_zdftmx_new 365 ! key_zdftmx_new: New internal wave-driven param: set a minimum value for Kz on TKE (ensure numerical stability) 366 zzd_up = zcof * ( MAX( avm(ji,jj,jk+1) + avm(ji,jj,jk), 2.e-5_wp ) ) & ! upper diagonal 367 & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk ) ) 368 zzd_lw = zcof * ( MAX( avm(ji,jj,jk) + avm(ji,jj,jk-1), 2.e-5_wp ) ) & ! lower diagonal 369 & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) 370 # else 347 371 zzd_up = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) & ! upper diagonal 348 372 & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk ) ) 349 373 zzd_lw = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) & ! lower diagonal 350 374 & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) 375 # endif 351 376 ! ! shear prod. at w-point weightened by mask 352 377 zesh2 = ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) & … … 408 433 END DO 409 434 435 ! ! Save TKE prior to nn_etau addition 436 e_niw(:,:,:) = en(:,:,:) 437 ! 410 438 ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 411 439 ! ! TKE due to surface and internal wave breaking 412 440 ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 441 IF( nn_htau == 2 ) THEN !* mixed-layer depth dependant length scale 442 DO jj = 2, jpjm1 443 DO ji = fs_2, fs_jpim1 ! vector opt. 444 htau(ji,jj) = rhtau * hmlp(ji,jj) 445 END DO 446 END DO 447 ENDIF 448 #if defined key_iomput 449 ! 450 CALL iom_put( "htau", htau(:,:) ) ! Check htau (even if constant in time) 451 #endif 452 ! 413 453 IF( nn_etau == 1 ) THEN !* penetration below the mixed layer (rn_efr fraction) 414 454 DO jk = 2, jpkm1 … … 445 485 END DO 446 486 END DO 487 ELSEIF( nn_etau == 4 ) THEN !* column integral independant of htau (rn_efr must be scaled up) 488 IF( nn_htau == 2 ) THEN ! efr dependant on time-varying htau 489 DO jj = 2, jpjm1 490 DO ji = fs_2, fs_jpim1 ! vector opt. 491 efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 492 END DO 493 END DO 494 ENDIF 495 DO jk = 2, jpkm1 496 DO jj = 2, jpjm1 497 DO ji = fs_2, fs_jpim1 ! vector opt. 498 en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) ) & 499 & * ( 1._wp - fr_i(ji,jj) ) * tmask(ji,jj,jk) 500 END DO 501 END DO 502 END DO 447 503 ENDIF 448 504 CALL lbc_lnk( en, 'W', 1. ) ! Lateral boundary conditions (sign unchanged) 505 ! 506 DO jk = 2, jpkm1 ! TKE budget: near-inertial waves term 507 DO jj = 2, jpjm1 508 DO ji = fs_2, fs_jpim1 ! vector opt. 509 e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk) 510 END DO 511 END DO 512 END DO 513 ! 514 CALL lbc_lnk( e_niw, 'W', 1. ) 449 515 ! 450 516 CALL wrk_dealloc( jpi,jpj, imlc ) ! integer … … 705 771 !!---------------------------------------------------------------------- 706 772 INTEGER :: ji, jj, jk ! dummy loop indices 707 INTEGER :: ios 773 INTEGER :: ios, ierr 708 774 !! 709 775 NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin , & 710 776 & rn_emin0, rn_bshear, nn_mxl , ln_mxl0 , & 711 777 & rn_mxl0 , nn_pdl , ln_lc , rn_lc , & 712 & nn_etau , nn_htau , rn_efr 778 & nn_etau , nn_htau , rn_efr , rn_c 713 779 !!---------------------------------------------------------------------- 714 ! 780 715 781 REWIND( numnam_ref ) ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy 716 782 READ ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901) … … 723 789 ! 724 790 ri_cri = 2._wp / ( 2._wp + rn_ediss / rn_ediff ) ! resulting critical Richardson number 791 # if defined key_zdftmx_new 792 ! key_zdftmx_new: New internal wave-driven param: specified value of rn_emin & rmxl_min are used 793 rn_emin = 1.e-10_wp 794 rmxl_min = 1.e-03_wp 795 IF(lwp) THEN ! Control print 796 WRITE(numout,*) 797 WRITE(numout,*) 'zdf_tke_init : New tidal mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3 ' 798 WRITE(numout,*) '~~~~~~~~~~~~' 799 ENDIF 800 # else 725 801 rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) ) ! resulting minimum length to recover molecular viscosity 802 # endif 726 803 ! 727 804 IF(lwp) THEN !* Control print … … 745 822 WRITE(numout,*) ' flag for computation of exp. tke profile nn_htau = ', nn_htau 746 823 WRITE(numout,*) ' fraction of en which pene. the thermocline rn_efr = ', rn_efr 824 WRITE(numout,*) ' fraction of TKE added within the mixed layer by nn_etau rn_c = ', rn_c 747 825 WRITE(numout,*) 748 826 WRITE(numout,*) ' critical Richardson nb with your parameters ri_cri = ', ri_cri … … 755 833 IF( nn_mxl < 0 .OR. nn_mxl > 3 ) CALL ctl_stop( 'bad flag: nn_mxl is 0, 1 or 2 ' ) 756 834 IF( nn_pdl < 0 .OR. nn_pdl > 1 ) CALL ctl_stop( 'bad flag: nn_pdl is 0 or 1 ' ) 757 IF( nn_htau < 0 .OR. nn_htau > 1 ) CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2' )835 IF( nn_htau < 0 .OR. nn_htau > 5 ) CALL ctl_stop( 'bad flag: nn_htau is 0 to 5 ' ) 758 836 IF( nn_etau == 3 .AND. .NOT. ln_cpl ) CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) 759 837 … … 763 841 ENDIF 764 842 765 IF( nn_etau == 2 ) CALL zdf_mxl( nit000 ) ! Initialization of nmln 843 IF( nn_etau == 2 ) THEN 844 ierr = zdf_mxl_alloc() 845 nmln(:,:) = nlb10 ! Initialization of nmln 846 ENDIF 847 848 IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN 849 ierr = zdf_mxl_alloc() 850 nmln(:,:) = nlb10 ! Initialization of nmln 851 ENDIF 766 852 767 853 ! !* depth of penetration of surface tke 768 854 IF( nn_etau /= 0 ) THEN 855 htau(:,:) = 0._wp 769 856 SELECT CASE( nn_htau ) ! Choice of the depth of penetration 770 857 CASE( 0 ) ! constant depth penetration (here 10 meters) … … 772 859 CASE( 1 ) ! F(latitude) : 0.5m to 30m poleward of 40 degrees 773 860 htau(:,:) = MAX( 0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) ) ) 861 CASE( 2 ) ! fraction of depth-integrated TKE within mixed-layer 862 rhtau = -1._wp / LOG( 1._wp - rn_c ) 863 CASE( 3 ) ! F(latitude) : 0.5m to 15m poleward of 20 degrees 864 htau(:,:) = MAX( 0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) ) ) 865 CASE( 4 ) ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south 866 DO jj = 2, jpjm1 867 DO ji = fs_2, fs_jpim1 ! vector opt. 868 IF( gphit(ji,jj) <= 0._wp ) THEN 869 htau(ji,jj) = MAX( 0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) ) ) 870 ELSE 871 htau(ji,jj) = MAX( 0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) ) ) 872 ENDIF 873 END DO 874 END DO 875 CASE ( 5 ) ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south, 876 DO jj = 2, jpjm1 ! 10m to 30m between 30/45 degrees south 877 DO ji = fs_2, fs_jpim1 ! vector opt. 878 IF( gphit(ji,jj) <= -30._wp ) THEN 879 htau(ji,jj) = MAX( 10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) ) ) 880 ELSE 881 htau(ji,jj) = MAX( 0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) ) ) 882 ENDIF 883 END DO 884 END DO 774 885 END SELECT 886 ! 887 IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN ! efr dependant on constant htau 888 DO jj = 2, jpjm1 889 DO ji = fs_2, fs_jpim1 ! vector opt. 890 efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) ) 891 END DO 892 END DO 893 ENDIF 775 894 ENDIF 776 895 ! !* set vertical eddy coef. to the background value
Note: See TracChangeset
for help on using the changeset viewer.