Changeset 3959
- Timestamp:
- 2013-07-08T18:46:13+02:00 (11 years ago)
- Location:
- branches/2013/dev_r3948_NOC_FK/NEMOGCM
- Files:
-
- 1 added
- 7 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2013/dev_r3948_NOC_FK/NEMOGCM/CONFIG/GYRE/EXP00/iodef.xml
r3940 r3959 47 47 <field field_ref="mldr10_1" name="somxl010" /> 48 48 <field field_ref="mldkz5" name="somixhgt" /> 49 </file> 49 <!-- variables available with MLE --> 50 <field field_ref="Lf_NHpf" name="Lf_NHpf" long_name="MLE:_Lf=NH/f" /> 51 </file> 50 52 51 53 <file id="file2" name_suffix="_grid_U" description="ocean U grid variables" > 52 54 <field field_ref="uoce" name="vozocrtx" /> 53 55 <field field_ref="utau" name="sozotaux" /> 56 <!-- variables available with MLE --> 57 <field field_ref="psiu_mle" name="psiu_mle" long_name="MLE_streamfunction_along_i-axis" /> 54 58 </file> 55 59 … … 57 61 <field field_ref="voce" name="vomecrty" /> 58 62 <field field_ref="vtau" name="sometauy" /> 63 <!-- variables available with MLE --> 64 <field field_ref="psiv_mle" name="psiv_mle" long_name="MLE_streamfunction_along_j-axis" /> 59 65 </file> 60 66 -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/CONFIG/GYRE/EXP00/namelist
r3795 r3959 586 586 ln_traadv_qck = .false. ! QUICKEST scheme 587 587 ln_traadv_msc_ups= .false. ! use upstream scheme within muscl 588 / 589 !----------------------------------------------------------------------- 590 &namtra_adv_mle ! mixed layer eddy parametrisation (Fox-Kemper param) 591 !----------------------------------------------------------------------- 592 ln_mle = .true. ! (T) use the Mixed Layer Eddy (MLE) parameterisation 593 rn_ce = 0.06 ! magnitude of the MLE (typical value: 0.06 to 0.08) 594 nn_mle = 1 ! MLE type: =0 standard Fox-Kemper ; =1 new formulation 595 rn_lf = 5.e+3 ! typical scale of mixed layer front (meters) (case rn_mle=0) 596 rn_time = 172800. ! time scale for mixing momentum across the mixed layer (seconds) (case rn_mle=0) 597 rn_lat = 20. ! reference latitude (degrees) of MLE coef. (case rn_mle=1) 598 nn_mld_uv = 0 ! space interpolation of MLD at u- & v-pts (0=min,1=averaged,2=max) 599 nn_conv = 0 ! =1 no MLE in case of convection ; =0 always MLE 600 rn_rho_c_mle = 0.01 ! delta rho criterion used to calculate MLD for FK 588 601 / 589 602 !---------------------------------------------------------------------------------- -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/CONFIG/ORCA2_LIM/EXP00/iodef_default.xml
r3940 r3959 64 64 <field field_ref="qsr" name="rsntds" long_name="surface_net_downward_shortwave_flux" /> 65 65 <field field_ref="qt" name="tohfls" long_name="surface_net_downward_total_heat_flux" /> 66 <!-- variables available with MLE --> 67 <field field_ref="Lf_NHpf" name="Lf_NHpf" long_name="MLE:_Lf=NH/f" /> 66 68 <field field_ref="taum" /> 67 69 <field field_ref="mldkz5" /> … … 73 75 <field field_ref="suoce" name="uos" long_name="sea_surface_x_velocity" /> 74 76 <field field_ref="utau" name="tauuo" long_name="surface_downward_x_stress" /> 77 <!-- variables available with MLE --> 78 <field field_ref="psiu_mle" name="psiu_mle" long_name="MLE_streamfunction_along_i-axis" /> 75 79 </file> 76 80 … … 79 83 <field field_ref="svoce" name="vos" long_name="sea_surface_y_velocity" /> 80 84 <field field_ref="vtau" name="tauvo" long_name="surface_downward_y_stress" /> 85 <!-- variables available with MLE --> 86 <field field_ref="psiv_mle" name="psiv_mle" long_name="MLE_streamfunction_along_j-axis" /> 81 87 </file> 82 88 -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/CONFIG/ORCA2_LIM/EXP00/namelist
r3795 r3959 581 581 ln_traadv_qck = .false. ! QUICKEST scheme 582 582 ln_traadv_msc_ups= .false. ! use upstream scheme within muscl 583 / 584 !----------------------------------------------------------------------- 585 &namtra_adv_mle ! mixed layer eddy parametrisation (Fox-Kemper param) 586 !----------------------------------------------------------------------- 587 ln_mle = .true. ! (T) use the Mixed Layer Eddy (MLE) parameterisation 588 rn_ce = 0.06 ! magnitude of the MLE (typical value: 0.06 to 0.08) 589 nn_mle = 1 ! MLE type: =0 standard Fox-Kemper ; =1 new formulation 590 rn_lf = 5.e+3 ! typical scale of mixed layer front (meters) (case rn_mle=0) 591 rn_time = 172800. ! time scale for mixing momentum across the mixed layer (seconds) (case rn_mle=0) 592 rn_lat = 20. ! reference latitude (degrees) of MLE coef. (case rn_mle=1) 593 nn_mld_uv = 0 ! space interpolation of MLD at u- & v-pts (0=min,1=averaged,2=max) 594 nn_conv = 0 ! =1 no MLE in case of convection ; =0 always MLE 595 rn_rho_c_mle = 0.01 ! delta rho criterion used to calculate MLD for FK 583 596 / 584 597 !---------------------------------------------------------------------------------- -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/CONFIG/SHARED/field_def.xml
r3905 r3959 29 29 <field id="mldr10_1" long_name="Mixed Layer Depth 0.01 ref.10m" unit="m" /> 30 30 <field id="rhop" long_name="potential density (sigma0)" unit="kg/m3" grid_ref="grid_T_3D"/> 31 <!-- variables available with MLE --> 32 <field id="Lf_NHpf" long_name="MLE: Lf = N H / f" unit="m" /> 31 33 <!-- next variables available with key_diahth --> 32 34 <field id="mlddzt" long_name="Thermocline Depth (max dT/dz)" unit="m" /> … … 140 142 <field id="uoce" long_name="ocean current along i-axis" unit="m/s" grid_ref="grid_U_3D" /> 141 143 <field id="uocetr_eff" long_name="Effective ocean transport along i-axis" unit="m3/s" grid_ref="grid_U_3D" /> 144 <!-- variables available with MLE --> 145 <field id="psiu_mle" long_name="MLE streamfunction along i-axis" unit="m3/s" grid_ref="grid_U_3D" /> 142 146 <!-- uoce_eiv: available with key_traldf_eiv and key_diaeiv --> 143 147 <field id="uoce_eiv" long_name="EIV ocean current along i-axis" unit="m/s" grid_ref="grid_U_3D" /> … … 159 163 <field id="voce" long_name="ocean current along j-axis" unit="m/s" grid_ref="grid_V_3D" /> 160 164 <field id="vocetr_eff" long_name="Effective ocean transport along j-axis" unit="m3/s" grid_ref="grid_V_3D" /> 165 <!-- variables available with MLE --> 166 <field id="psiv_mle" long_name="MLE streamfunction along j-axis" unit="m3/s" grid_ref="grid_U_3D" /> 161 167 <!-- voce_eiv: available with key_traldf_eiv and key_diaeiv --> 162 168 <field id="voce_eiv" long_name="EIV ocean current along j-axis" unit="m/s" grid_ref="grid_V_3D" /> -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90
r3718 r3959 6 6 !! History : 2.0 ! 2005-11 (G. Madec) Original code 7 7 !! 3.3 ! 2010-09 (C. Ethe, G. Madec) merge TRC-TRA + switch from velocity to transport 8 !! 4.0 ! 2011-06 (G. Madec) Addition of Mixed Layer Eddy parameterisation 8 9 !!---------------------------------------------------------------------- 9 10 … … 21 22 USE traadv_qck ! QUICKEST scheme (tra_adv_qck routine) 22 23 USE traadv_eiv ! eddy induced velocity (tra_adv_eiv routine) 24 USE traadv_mle ! ML eddy induced velocity (tra_adv_eiv routine) 23 25 USE cla ! cross land advection (cla_traadv routine) 24 26 USE ldftra_oce ! lateral diffusion coefficient on tracers … … 99 101 & CALL tra_adv_eiv( kt, nit000, zun, zvn, zwn, 'TRA' ) ! add the eiv transport (if necessary) 100 102 ! 103 IF( ln_mle ) CALL tra_adv_mle( kt, nit000, zun, zvn, zwn, 'TRA' ) ! add the mle transport (if necessary) 101 104 CALL iom_put( "uocetr_eff", zun ) ! output effective transport 102 105 CALL iom_put( "vocetr_eff", zvn ) … … 206 209 ENDIF 207 210 ! 211 CALL tra_adv_mle_init ! initialisation of the Mixed Layer Eddy parametrisation (MLE) 212 ! 208 213 END SUBROUTINE tra_adv_init 209 214 -
branches/2013/dev_r3948_NOC_FK/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90
r3294 r3959 24 24 25 25 PUBLIC zdf_mxl ! called by step.F90 26 27 REAL(wp), PUBLIC :: rho_c = 0.01_wp ! density criterion for mixed layer depth 28 REAL(wp), PUBLIC :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth 26 29 27 30 INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: nmln !: number of level in the mixed layer (used by TOP) … … 63 66 !! ** Method : The mixed layer depth is the shallowest W depth with 64 67 !! the density of the corresponding T point (just bellow) bellow a 65 !! given value defined locally as rho(10m) + zrho_c68 !! given value defined locally as rho(10m) + rho_c 66 69 !! The turbocline depth is the depth at which the vertical 67 70 !! eddy diffusivity coefficient (resulting from the vertical physics … … 76 79 INTEGER :: iikn, iiki ! temporary integer within a do loop 77 80 INTEGER, POINTER, DIMENSION(:,:) :: imld ! temporary workspace 78 REAL(wp) :: zrho_c = 0.01_wp ! density criterion for mixed layer depth79 REAL(wp) :: zavt_c = 5.e-4_wp ! Kz criterion for the turbocline depth80 81 !!---------------------------------------------------------------------- 81 82 ! … … 98 99 DO jj = 1, jpj 99 100 DO ji = 1, jpi 100 IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + zrho_c ) nmln(ji,jj) = jk ! Mixed layer101 IF( avt (ji,jj,jk) < zavt_c ) imld(ji,jj) = jk ! Turbocline101 IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + rho_c ) nmln(ji,jj) = jk ! Mixed layer 102 IF( avt (ji,jj,jk) < avt_c ) imld(ji,jj) = jk ! Turbocline 102 103 END DO 103 104 END DO
Note: See TracChangeset
for help on using the changeset viewer.