[1531] | 1 | MODULE zdftke |
---|
[1239] | 2 | !!====================================================================== |
---|
[1531] | 3 | !! *** MODULE zdftke *** |
---|
[1239] | 4 | !! Ocean physics: vertical mixing coefficient computed from the tke |
---|
| 5 | !! turbulent closure parameterization |
---|
| 6 | !!===================================================================== |
---|
[1492] | 7 | !! History : OPA ! 1991-03 (b. blanke) Original code |
---|
| 8 | !! 7.0 ! 1991-11 (G. Madec) bug fix |
---|
| 9 | !! 7.1 ! 1992-10 (G. Madec) new mixing length and eav |
---|
| 10 | !! 7.2 ! 1993-03 (M. Guyon) symetrical conditions |
---|
| 11 | !! 7.3 ! 1994-08 (G. Madec, M. Imbard) nn_pdl flag |
---|
| 12 | !! 7.5 ! 1996-01 (G. Madec) s-coordinates |
---|
| 13 | !! 8.0 ! 1997-07 (G. Madec) lbc |
---|
| 14 | !! 8.1 ! 1999-01 (E. Stretta) new option for the mixing length |
---|
| 15 | !! NEMO 1.0 ! 2002-06 (G. Madec) add tke_init routine |
---|
| 16 | !! - ! 2004-10 (C. Ethe ) 1D configuration |
---|
| 17 | !! 2.0 ! 2006-07 (S. Masson) distributed restart using iom |
---|
| 18 | !! 3.0 ! 2008-05 (C. Ethe, G.Madec) : update TKE physics: |
---|
| 19 | !! ! - tke penetration (wind steering) |
---|
| 20 | !! ! - suface condition for tke & mixing length |
---|
| 21 | !! ! - Langmuir cells |
---|
| 22 | !! - ! 2008-05 (J.-M. Molines, G. Madec) 2D form of avtb |
---|
| 23 | !! - ! 2008-06 (G. Madec) style + DOCTOR name for namelist parameters |
---|
| 24 | !! - ! 2008-12 (G. Reffray) stable discretization of the production term |
---|
| 25 | !! 3.2 ! 2009-06 (G. Madec, S. Masson) TKE restart compatible with key_cpl |
---|
| 26 | !! ! + cleaning of the parameters + bugs correction |
---|
[2528] | 27 | !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase |
---|
[5120] | 28 | !! 3.6 ! 2014-11 (P. Mathiot) add ice shelf capability |
---|
[9019] | 29 | !! 4.0 ! 2017-04 (G. Madec) remove CPP ddm key & avm at t-point only |
---|
| 30 | !! - ! 2017-05 (G. Madec) add top/bottom friction as boundary condition (ln_drg) |
---|
[1239] | 31 | !!---------------------------------------------------------------------- |
---|
[9019] | 32 | |
---|
[1239] | 33 | !!---------------------------------------------------------------------- |
---|
[3625] | 34 | !! zdf_tke : update momentum and tracer Kz from a tke scheme |
---|
| 35 | !! tke_tke : tke time stepping: update tke at now time step (en) |
---|
| 36 | !! tke_avn : compute mixing length scale and deduce avm and avt |
---|
| 37 | !! zdf_tke_init : initialization, namelist read, and parameters control |
---|
| 38 | !! tke_rst : read/write tke restart in ocean restart file |
---|
[1239] | 39 | !!---------------------------------------------------------------------- |
---|
[2528] | 40 | USE oce ! ocean: dynamics and active tracers variables |
---|
| 41 | USE phycst ! physical constants |
---|
| 42 | USE dom_oce ! domain: ocean |
---|
| 43 | USE domvvl ! domain: variable volume layer |
---|
[1492] | 44 | USE sbc_oce ! surface boundary condition: ocean |
---|
[9019] | 45 | USE zdfdrg ! vertical physics: top/bottom drag coef. |
---|
[2528] | 46 | USE zdfmxl ! vertical physics: mixed layer |
---|
[9019] | 47 | ! |
---|
[1492] | 48 | USE in_out_manager ! I/O manager |
---|
| 49 | USE iom ! I/O manager library |
---|
[2715] | 50 | USE lib_mpp ! MPP library |
---|
[9019] | 51 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 52 | USE prtctl ! Print control |
---|
[3625] | 53 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
[12991] | 54 | USE sbcwave ! Surface boundary waves |
---|
[1239] | 55 | |
---|
| 56 | IMPLICIT NONE |
---|
| 57 | PRIVATE |
---|
| 58 | |
---|
[2528] | 59 | PUBLIC zdf_tke ! routine called in step module |
---|
| 60 | PUBLIC zdf_tke_init ! routine called in opa module |
---|
| 61 | PUBLIC tke_rst ! routine called in step module |
---|
[1239] | 62 | |
---|
[4147] | 63 | ! !!** Namelist namzdf_tke ** |
---|
| 64 | LOGICAL :: ln_mxl0 ! mixing length scale surface value as function of wind stress or not |
---|
[12991] | 65 | LOGICAL :: ln_mxhsw ! mixing length scale surface value as a fonction of wave height |
---|
[4147] | 66 | INTEGER :: nn_mxl ! type of mixing length (=0/1/2/3) |
---|
| 67 | REAL(wp) :: rn_mxl0 ! surface min value of mixing length (kappa*z_o=0.4*0.1 m) [m] |
---|
| 68 | INTEGER :: nn_pdl ! Prandtl number or not (ratio avt/avm) (=0/1) |
---|
| 69 | REAL(wp) :: rn_ediff ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e) |
---|
| 70 | REAL(wp) :: rn_ediss ! coefficient of the Kolmogoroff dissipation |
---|
| 71 | REAL(wp) :: rn_ebb ! coefficient of the surface input of tke |
---|
| 72 | REAL(wp) :: rn_emin ! minimum value of tke [m2/s2] |
---|
| 73 | REAL(wp) :: rn_emin0 ! surface minimum value of tke [m2/s2] |
---|
| 74 | REAL(wp) :: rn_bshear ! background shear (>0) currently a numerical threshold (do not change it) |
---|
[9019] | 75 | LOGICAL :: ln_drg ! top/bottom friction forcing flag |
---|
[4147] | 76 | INTEGER :: nn_etau ! type of depth penetration of surface tke (=0/1/2/3) |
---|
[9019] | 77 | INTEGER :: nn_htau ! type of tke profile of penetration (=0/1) |
---|
[12991] | 78 | INTEGER :: nn_bc_surf! surface condition (0/1=Dir/Neum) ! Only applicable for wave coupling |
---|
| 79 | INTEGER :: nn_bc_bot ! surface condition (0/1=Dir/Neum) ! Only applicable for wave coupling |
---|
[9019] | 80 | REAL(wp) :: rn_efr ! fraction of TKE surface value which penetrates in the ocean |
---|
[9546] | 81 | REAL(wp) :: rn_eice ! =0 ON below sea-ice, =4 OFF when ice fraction > 1/4 |
---|
[4147] | 82 | LOGICAL :: ln_lc ! Langmuir cells (LC) as a source term of TKE or not |
---|
[9019] | 83 | REAL(wp) :: rn_lc ! coef to compute vertical velocity of Langmuir cells |
---|
[1239] | 84 | |
---|
[4147] | 85 | REAL(wp) :: ri_cri ! critic Richardson number (deduced from rn_ediff and rn_ediss values) |
---|
| 86 | REAL(wp) :: rmxl_min ! minimum mixing length value (deduced from rn_ediff and rn_emin values) [m] |
---|
[2528] | 87 | REAL(wp) :: rhftau_add = 1.e-3_wp ! add offset applied to HF part of taum (nn_etau=3) |
---|
| 88 | REAL(wp) :: rhftau_scl = 1.0_wp ! scale factor applied to HF part of taum (nn_etau=3) |
---|
[1239] | 89 | |
---|
[9019] | 90 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: htau ! depth of tke penetration (nn_htau) |
---|
| 91 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: dissl ! now mixing lenght of dissipation |
---|
| 92 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: apdlr ! now mixing lenght of dissipation |
---|
[1492] | 93 | |
---|
[1239] | 94 | !! * Substitutions |
---|
[12377] | 95 | # include "do_loop_substitute.h90" |
---|
[1239] | 96 | !!---------------------------------------------------------------------- |
---|
[9598] | 97 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
[2528] | 98 | !! $Id$ |
---|
[10068] | 99 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[1239] | 100 | !!---------------------------------------------------------------------- |
---|
| 101 | CONTAINS |
---|
| 102 | |
---|
[2715] | 103 | INTEGER FUNCTION zdf_tke_alloc() |
---|
| 104 | !!---------------------------------------------------------------------- |
---|
| 105 | !! *** FUNCTION zdf_tke_alloc *** |
---|
| 106 | !!---------------------------------------------------------------------- |
---|
[9019] | 107 | ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) , STAT= zdf_tke_alloc ) |
---|
| 108 | ! |
---|
[10425] | 109 | CALL mpp_sum ( 'zdftke', zdf_tke_alloc ) |
---|
| 110 | IF( zdf_tke_alloc /= 0 ) CALL ctl_stop( 'STOP', 'zdf_tke_alloc: failed to allocate arrays' ) |
---|
[2715] | 111 | ! |
---|
| 112 | END FUNCTION zdf_tke_alloc |
---|
| 113 | |
---|
| 114 | |
---|
[12377] | 115 | SUBROUTINE zdf_tke( kt, Kbb, Kmm, p_sh2, p_avm, p_avt ) |
---|
[1239] | 116 | !!---------------------------------------------------------------------- |
---|
[1531] | 117 | !! *** ROUTINE zdf_tke *** |
---|
[1239] | 118 | !! |
---|
| 119 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
[1492] | 120 | !! coefficients using a turbulent closure scheme (TKE). |
---|
[1239] | 121 | !! |
---|
[1492] | 122 | !! ** Method : The time evolution of the turbulent kinetic energy (tke) |
---|
| 123 | !! is computed from a prognostic equation : |
---|
| 124 | !! d(en)/dt = avm (d(u)/dz)**2 ! shear production |
---|
| 125 | !! + d( avm d(en)/dz )/dz ! diffusion of tke |
---|
| 126 | !! + avt N^2 ! stratif. destruc. |
---|
| 127 | !! - rn_ediss / emxl en**(2/3) ! Kolmogoroff dissipation |
---|
[1239] | 128 | !! with the boundary conditions: |
---|
[1695] | 129 | !! surface: en = max( rn_emin0, rn_ebb * taum ) |
---|
[1239] | 130 | !! bottom : en = rn_emin |
---|
[1492] | 131 | !! The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff) |
---|
| 132 | !! |
---|
| 133 | !! The now Turbulent kinetic energy is computed using the following |
---|
| 134 | !! time stepping: implicit for vertical diffusion term, linearized semi |
---|
| 135 | !! implicit for kolmogoroff dissipation term, and explicit forward for |
---|
| 136 | !! both buoyancy and shear production terms. Therefore a tridiagonal |
---|
| 137 | !! linear system is solved. Note that buoyancy and shear terms are |
---|
| 138 | !! discretized in a energy conserving form (Bruchard 2002). |
---|
| 139 | !! |
---|
| 140 | !! The dissipative and mixing length scale are computed from en and |
---|
| 141 | !! the stratification (see tke_avn) |
---|
| 142 | !! |
---|
| 143 | !! The now vertical eddy vicosity and diffusivity coefficients are |
---|
| 144 | !! given by: |
---|
| 145 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
| 146 | !! avt = max( avmb, pdl * avm ) |
---|
[1239] | 147 | !! eav = max( avmb, avm ) |
---|
[1492] | 148 | !! where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and |
---|
| 149 | !! given by an empirical funtion of the localRichardson number if nn_pdl=1 |
---|
[1239] | 150 | !! |
---|
| 151 | !! ** Action : compute en (now turbulent kinetic energy) |
---|
[9019] | 152 | !! update avt, avm (before vertical eddy coef.) |
---|
[1239] | 153 | !! |
---|
| 154 | !! References : Gaspar et al., JGR, 1990, |
---|
| 155 | !! Blanke and Delecluse, JPO, 1991 |
---|
| 156 | !! Mellor and Blumberg, JPO 2004 |
---|
| 157 | !! Axell, JGR, 2002 |
---|
[1492] | 158 | !! Bruchard OM 2002 |
---|
[1239] | 159 | !!---------------------------------------------------------------------- |
---|
[9019] | 160 | INTEGER , INTENT(in ) :: kt ! ocean time step |
---|
[12377] | 161 | INTEGER , INTENT(in ) :: Kbb, Kmm ! ocean time level indices |
---|
[9019] | 162 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: p_sh2 ! shear production term |
---|
| 163 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: p_avm, p_avt ! momentum and tracer Kz (w-points) |
---|
[1492] | 164 | !!---------------------------------------------------------------------- |
---|
[1481] | 165 | ! |
---|
[12377] | 166 | CALL tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt ) ! now tke (en) |
---|
[5656] | 167 | ! |
---|
[12377] | 168 | CALL tke_avn( Kbb, Kmm, p_avm, p_avt ) ! now avt, avm, dissl |
---|
[3632] | 169 | ! |
---|
[5656] | 170 | END SUBROUTINE zdf_tke |
---|
[1239] | 171 | |
---|
[1492] | 172 | |
---|
[12377] | 173 | SUBROUTINE tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt ) |
---|
[1239] | 174 | !!---------------------------------------------------------------------- |
---|
[1492] | 175 | !! *** ROUTINE tke_tke *** |
---|
| 176 | !! |
---|
| 177 | !! ** Purpose : Compute the now Turbulente Kinetic Energy (TKE) |
---|
| 178 | !! |
---|
| 179 | !! ** Method : - TKE surface boundary condition |
---|
[2528] | 180 | !! - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T) |
---|
[9019] | 181 | !! - source term due to shear (= Kz dz[Ub] * dz[Un] ) |
---|
[1492] | 182 | !! - Now TKE : resolution of the TKE equation by inverting |
---|
| 183 | !! a tridiagonal linear system by a "methode de chasse" |
---|
| 184 | !! - increase TKE due to surface and internal wave breaking |
---|
[9019] | 185 | !! NB: when sea-ice is present, both LC parameterization |
---|
| 186 | !! and TKE penetration are turned off when the ice fraction |
---|
| 187 | !! is smaller than 0.25 |
---|
[1492] | 188 | !! |
---|
| 189 | !! ** Action : - en : now turbulent kinetic energy) |
---|
[1239] | 190 | !! --------------------------------------------------------------------- |
---|
[9019] | 191 | USE zdf_oce , ONLY : en ! ocean vertical physics |
---|
| 192 | !! |
---|
[12377] | 193 | INTEGER , INTENT(in ) :: Kbb, Kmm ! ocean time level indices |
---|
[9019] | 194 | REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: p_sh2 ! shear production term |
---|
| 195 | REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: p_avm, p_avt ! vertical eddy viscosity & diffusivity (w-points) |
---|
| 196 | ! |
---|
| 197 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
| 198 | REAL(wp) :: zetop, zebot, zmsku, zmskv ! local scalars |
---|
| 199 | REAL(wp) :: zrhoa = 1.22 ! Air density kg/m3 |
---|
| 200 | REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient |
---|
| 201 | REAL(wp) :: zbbrau, zri ! local scalars |
---|
| 202 | REAL(wp) :: zfact1, zfact2, zfact3 ! - - |
---|
| 203 | REAL(wp) :: ztx2 , zty2 , zcof ! - - |
---|
| 204 | REAL(wp) :: ztau , zdif ! - - |
---|
| 205 | REAL(wp) :: zus , zwlc , zind ! - - |
---|
| 206 | REAL(wp) :: zzd_up, zzd_lw ! - - |
---|
[12991] | 207 | REAL(wp) :: ztaui, ztauj, z1_norm |
---|
[9019] | 208 | INTEGER , DIMENSION(jpi,jpj) :: imlc |
---|
[12991] | 209 | REAL(wp), DIMENSION(jpi,jpj) :: zhlc, zfr_i, zWlc2 |
---|
[9019] | 210 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: zpelc, zdiag, zd_up, zd_lw |
---|
[1239] | 211 | !!-------------------------------------------------------------------- |
---|
[1492] | 212 | ! |
---|
[12489] | 213 | zbbrau = rn_ebb / rho0 ! Local constant initialisation |
---|
| 214 | zfact1 = -.5_wp * rn_Dt |
---|
| 215 | zfact2 = 1.5_wp * rn_Dt * rn_ediss |
---|
[2528] | 216 | zfact3 = 0.5_wp * rn_ediss |
---|
[1492] | 217 | ! |
---|
[12991] | 218 | zpelc(:,:,:) = 0._wp ! need to be initialised in case ln_lc is not used |
---|
| 219 | ! |
---|
[1492] | 220 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[9019] | 221 | ! ! Surface/top/bottom boundary condition on tke |
---|
[1492] | 222 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[12698] | 223 | ! |
---|
[12377] | 224 | DO_2D_00_00 |
---|
| 225 | en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1) |
---|
| 226 | END_2D |
---|
[9019] | 227 | ! |
---|
[1492] | 228 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 229 | ! ! Bottom boundary condition on tke |
---|
| 230 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[1719] | 231 | ! |
---|
[12489] | 232 | ! en(bot) = (ebb0/rho0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin) |
---|
[9019] | 233 | ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75) |
---|
| 234 | ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2 |
---|
[1492] | 235 | ! |
---|
[9019] | 236 | IF( ln_drg ) THEN !== friction used as top/bottom boundary condition on TKE |
---|
| 237 | ! |
---|
[12377] | 238 | DO_2D_00_00 |
---|
| 239 | zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) ) |
---|
| 240 | zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) ) |
---|
[12489] | 241 | ! ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0) |
---|
[12377] | 242 | zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT( ( zmsku*( uu(ji,jj,mbkt(ji,jj),Kbb)+uu(ji-1,jj,mbkt(ji,jj),Kbb) ) )**2 & |
---|
| 243 | & + ( zmskv*( vv(ji,jj,mbkt(ji,jj),Kbb)+vv(ji,jj-1,mbkt(ji,jj),Kbb) ) )**2 ) |
---|
| 244 | en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj) |
---|
| 245 | END_2D |
---|
[9019] | 246 | IF( ln_isfcav ) THEN ! top friction |
---|
[12377] | 247 | DO_2D_00_00 |
---|
| 248 | zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) ) |
---|
| 249 | zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) ) |
---|
[12489] | 250 | ! ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0) |
---|
[12377] | 251 | zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT( ( zmsku*( uu(ji,jj,mikt(ji,jj),Kbb)+uu(ji-1,jj,mikt(ji,jj),Kbb) ) )**2 & |
---|
| 252 | & + ( zmskv*( vv(ji,jj,mikt(ji,jj),Kbb)+vv(ji,jj-1,mikt(ji,jj),Kbb) ) )**2 ) |
---|
[12702] | 253 | ! (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj) = 1 where ice shelves are present |
---|
| 254 | en(ji,jj,mikt(ji,jj)) = en(ji,jj,1) * tmask(ji,jj,1) & |
---|
| 255 | & + MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj) |
---|
[12377] | 256 | END_2D |
---|
[9019] | 257 | ENDIF |
---|
| 258 | ! |
---|
| 259 | ENDIF |
---|
| 260 | ! |
---|
[1492] | 261 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[12991] | 262 | IF( ln_lc ) THEN ! Langmuir circulation source term added to tke (Axell JGR 2002) |
---|
[1492] | 263 | ! !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[1239] | 264 | ! |
---|
[12991] | 265 | ! !* Langmuir velocity scale |
---|
| 266 | ! |
---|
| 267 | IF ( cpl_sdrftx ) THEN ! Surface Stokes Drift available |
---|
| 268 | ! ! Craik-Leibovich velocity scale Wlc = ( u* u_s )^1/2 with u* = (taum/rho0)^1/2 |
---|
| 269 | ! ! associated kinetic energy : 1/2 (Wlc)^2 = u* u_s |
---|
| 270 | ! ! more precisely, it is the dot product that must be used : |
---|
| 271 | ! ! 1/2 (W_lc)^2 = MAX( u* u_s + v* v_s , 0 ) only the positive part |
---|
| 272 | !!gm ! PS: currently we don't have neither the 2 stress components at t-point !nor the angle between u* and u_s |
---|
| 273 | !!gm ! so we will overestimate the LC velocity.... !!gm I will do the work if !LC have an effect ! |
---|
| 274 | DO_2D_00_00 |
---|
| 275 | !!XC zWlc2(ji,jj) = 0.5_wp * SQRT( taum(ji,jj) * r1_rho0 * ( ut0sd(ji,jj)**2 +vt0sd(ji,jj)**2 ) ) |
---|
| 276 | zWlc2(ji,jj) = 0.5_wp * ( ut0sd(ji,jj)**2 +vt0sd(ji,jj)**2 ) |
---|
| 277 | END_2D |
---|
| 278 | ! |
---|
| 279 | ! Projection of Stokes drift in the wind stress direction |
---|
| 280 | ! |
---|
| 281 | DO_2D_00_00 |
---|
| 282 | ztaui = 0.5_wp * ( utau(ji,jj) + utau(ji-1,jj) ) |
---|
| 283 | ztauj = 0.5_wp * ( vtau(ji,jj) + vtau(ji,jj-1) ) |
---|
| 284 | z1_norm = 1._wp / MAX( SQRT(ztaui*ztaui+ztauj*ztauj), 1.e-12 ) * tmask(ji,jj,1) |
---|
| 285 | zWlc2(ji,jj) = 0.5_wp * z1_norm * ( MAX( ut0sd(ji,jj)*ztaui + vt0sd(ji,jj)*ztauj, 0._wp ) )**2 |
---|
| 286 | END_2D |
---|
| 287 | CALL lbc_lnk ( 'zdftke', zWlc2, 'T', 1. ) |
---|
| 288 | ! |
---|
| 289 | ELSE ! Surface Stokes drift deduced from surface stress |
---|
| 290 | ! ! Wlc = u_s with u_s = 0.016*U_10m, the surface stokes drift (Axell 2002, Eq.44) |
---|
| 291 | ! ! using |tau| = rho_air Cd |U_10m|^2 , it comes: |
---|
| 292 | ! ! Wlc = 0.016 * [|tau|/(rho_air Cdrag) ]^1/2 and thus: |
---|
| 293 | ! ! 1/2 Wlc^2 = 0.5 * 0.016 * 0.016 |tau| /( rho_air Cdrag ) |
---|
| 294 | zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag ) ! to convert stress in 10m wind using a constant drag |
---|
| 295 | DO_2D_11_11 |
---|
| 296 | zWlc2(ji,jj) = zcof * taum(ji,jj) |
---|
| 297 | END_2D |
---|
| 298 | ! |
---|
| 299 | ENDIF |
---|
| 300 | ! |
---|
| 301 | ! !* Depth of the LC circulation (Axell 2002, Eq.47) |
---|
| 302 | ! !- LHS of Eq.47 |
---|
[12377] | 303 | zpelc(:,:,1) = MAX( rn2b(:,:,1), 0._wp ) * gdepw(:,:,1,Kmm) * e3w(:,:,1,Kmm) |
---|
[1239] | 304 | DO jk = 2, jpk |
---|
[12377] | 305 | zpelc(:,:,jk) = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * gdepw(:,:,jk,Kmm) * e3w(:,:,jk,Kmm) |
---|
[1239] | 306 | END DO |
---|
[12991] | 307 | ! |
---|
| 308 | ! !- compare LHS to RHS of Eq.47 |
---|
[7753] | 309 | imlc(:,:) = mbkt(:,:) + 1 ! Initialization to the number of w ocean point (=2 over land) |
---|
[12377] | 310 | DO_3DS_11_11( jpkm1, 2, -1 ) |
---|
[12991] | 311 | IF( zpelc(ji,jj,jk) > zWlc2(ji,jj) ) imlc(ji,jj) = jk |
---|
[12377] | 312 | END_3D |
---|
[1492] | 313 | ! ! finite LC depth |
---|
[12377] | 314 | DO_2D_11_11 |
---|
| 315 | zhlc(ji,jj) = gdepw(ji,jj,imlc(ji,jj),Kmm) |
---|
| 316 | END_2D |
---|
[12991] | 317 | ! |
---|
[1705] | 318 | zcof = 0.016 / SQRT( zrhoa * zcdrag ) |
---|
[12377] | 319 | DO_2D_00_00 |
---|
[12991] | 320 | zus = SQRT( 2. * zWlc2(ji,jj) ) ! Stokes drift |
---|
[12377] | 321 | zfr_i(ji,jj) = ( 1._wp - 4._wp * fr_i(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok |
---|
| 322 | IF (zfr_i(ji,jj) < 0. ) zfr_i(ji,jj) = 0. |
---|
| 323 | END_2D |
---|
| 324 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 325 | IF ( zfr_i(ji,jj) /= 0. ) THEN |
---|
| 326 | ! vertical velocity due to LC |
---|
| 327 | IF ( gdepw(ji,jj,jk,Kmm) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN |
---|
| 328 | ! ! vertical velocity due to LC |
---|
| 329 | zwlc = rn_lc * SIN( rpi * gdepw(ji,jj,jk,Kmm) / zhlc(ji,jj) ) ! warning: optimization: zus^3 is in zfr_i |
---|
| 330 | ! ! TKE Langmuir circulation source term |
---|
[12489] | 331 | en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * zfr_i(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) |
---|
[12377] | 332 | ENDIF |
---|
| 333 | ENDIF |
---|
| 334 | END_3D |
---|
[1239] | 335 | ! |
---|
| 336 | ENDIF |
---|
[1492] | 337 | ! |
---|
| 338 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 339 | ! ! Now Turbulent kinetic energy (output in en) |
---|
| 340 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 341 | ! ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
| 342 | ! ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ). |
---|
| 343 | ! ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal |
---|
| 344 | ! |
---|
[9019] | 345 | IF( nn_pdl == 1 ) THEN !* Prandtl number = F( Ri ) |
---|
[12377] | 346 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 347 | ! ! local Richardson number |
---|
| 348 | zri = MAX( rn2b(ji,jj,jk), 0._wp ) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear ) |
---|
| 349 | ! ! inverse of Prandtl number |
---|
| 350 | apdlr(ji,jj,jk) = MAX( 0.1_wp, ri_cri / MAX( ri_cri , zri ) ) |
---|
| 351 | END_3D |
---|
[5656] | 352 | ENDIF |
---|
[5836] | 353 | ! |
---|
[12377] | 354 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 355 | zcof = zfact1 * tmask(ji,jj,jk) |
---|
| 356 | ! ! A minimum of 2.e-5 m2/s is imposed on TKE vertical |
---|
| 357 | ! ! eddy coefficient (ensure numerical stability) |
---|
| 358 | zzd_up = zcof * MAX( p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk ) , 2.e-5_wp ) & ! upper diagonal |
---|
| 359 | & / ( e3t(ji,jj,jk ,Kmm) * e3w(ji,jj,jk ,Kmm) ) |
---|
| 360 | zzd_lw = zcof * MAX( p_avm(ji,jj,jk ) + p_avm(ji,jj,jk-1) , 2.e-5_wp ) & ! lower diagonal |
---|
| 361 | & / ( e3t(ji,jj,jk-1,Kmm) * e3w(ji,jj,jk ,Kmm) ) |
---|
| 362 | ! |
---|
| 363 | zd_up(ji,jj,jk) = zzd_up ! Matrix (zdiag, zd_up, zd_lw) |
---|
| 364 | zd_lw(ji,jj,jk) = zzd_lw |
---|
| 365 | zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk) |
---|
| 366 | ! |
---|
| 367 | ! ! right hand side in en |
---|
[12489] | 368 | en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * ( p_sh2(ji,jj,jk) & ! shear |
---|
[12377] | 369 | & - p_avt(ji,jj,jk) * rn2(ji,jj,jk) & ! stratification |
---|
| 370 | & + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk) & ! dissipation |
---|
| 371 | & ) * wmask(ji,jj,jk) |
---|
| 372 | END_3D |
---|
[12991] | 373 | ! |
---|
| 374 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 375 | ! ! choose to keep coherence with previous estimation of |
---|
| 376 | ! ! Surface boundary condition on tke |
---|
| 377 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 378 | ! [EC] Should we keep this?? |
---|
| 379 | IF ( ln_isfcav ) THEN |
---|
| 380 | DO_2D_11_11 ! en(mikt(ji,jj)) = rn_emin |
---|
| 381 | en(ji,jj,mikt(ji,jj))=rn_emin * tmask(ji,jj,1) |
---|
| 382 | END_2D |
---|
| 383 | END IF |
---|
| 384 | |
---|
| 385 | IF ( cpl_phioc .and. ln_phioc ) THEN |
---|
| 386 | SELECT CASE (nn_bc_surf) !! Dirichlet Boundary Condition using surface TKE flux from waves |
---|
| 387 | |
---|
| 388 | CASE ( 0 ) |
---|
| 389 | |
---|
| 390 | DO_2D_00_00 ! en(1) = rn_ebb taum / rho0 (min value rn_emin0) |
---|
| 391 | IF ( phioc(ji,jj) < 0 ) phioc(ji,jj) = 0._wp |
---|
| 392 | en(ji,jj,1) = MAX( rn_emin0, .5 * ( 15.8 * phioc(ji,jj) / rho0 )**(2./3.) ) * tmask(ji,jj,1) |
---|
| 393 | zdiag(ji,jj,1) = 1._wp/en(ji,jj,1) ! choose to keep coherence with former estimation of |
---|
| 394 | zd_lw(ji,jj,1) = 1._wp ! zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1) |
---|
| 395 | zd_up(ji,jj,1) = 0._wp |
---|
| 396 | END_2D |
---|
| 397 | |
---|
| 398 | CASE ( 1 ) |
---|
| 399 | DO_2D_00_00 |
---|
| 400 | IF ( phioc(ji,jj) < 0 ) phioc(ji,jj) = 0._wp |
---|
| 401 | en(ji,jj,2) = en(ji,jj,2) + ( rn_Dt * phioc(ji,jj) / rho0 ) /e3w(ji,jj,2,Kmm) |
---|
| 402 | en(ji,jj,1) = en(ji,jj,2) + (2 * e3t(ji,jj,1,Kmm) * phioc(ji,jj)) / ( p_avm(ji,jj,1) + p_avm(ji,jj,2) ) |
---|
| 403 | zdiag(ji,jj,2) = zdiag(ji,jj,2) + zd_lw(ji,jj,2) |
---|
| 404 | zd_lw(ji,jj,2) = 0._wp |
---|
| 405 | zd_up(ji,jj,1) = 0._wp |
---|
| 406 | END_2D |
---|
| 407 | ! |
---|
| 408 | END SELECT |
---|
| 409 | |
---|
| 410 | ELSE ! TKE Dirichlet boundary condition (without wave coupling) |
---|
| 411 | |
---|
| 412 | DO_2D_00_00 ! en(1) = rn_ebb taum / rho0 (min value rn_emin0) |
---|
| 413 | en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1) |
---|
| 414 | zdiag(ji,jj,1) = 1._wp/en(ji,jj,1) ! choose to keep coherence with former estimation of |
---|
| 415 | zd_lw(ji,jj,1) = 1._wp ! zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1) |
---|
| 416 | zd_up(ji,jj,1) = 0._wp |
---|
| 417 | END_2D |
---|
| 418 | |
---|
| 419 | ENDIF |
---|
| 420 | ! |
---|
[5120] | 421 | ! !* Matrix inversion from level 2 (tke prescribed at level 1) |
---|
[12991] | 422 | ! DO_3D_00_00( 3, jpkm1 ) |
---|
| 423 | DO_3D_00_00( 2, jpkm1 ) ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
[12377] | 424 | zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1) |
---|
| 425 | END_3D |
---|
[12991] | 426 | !XC : commented to allow for neumann boundary condition |
---|
| 427 | ! DO_2D_00_00 |
---|
| 428 | ! zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1) ! Surface boudary conditions on tke |
---|
| 429 | ! END_2D |
---|
| 430 | ! DO_3D_00_00( 3, jpkm1 ) |
---|
| 431 | DO_3D_00_00( 2, jpkm1 ) ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
[12377] | 432 | zd_lw(ji,jj,jk) = en(ji,jj,jk) - zd_lw(ji,jj,jk) / zdiag(ji,jj,jk-1) *zd_lw(ji,jj,jk-1) |
---|
| 433 | END_3D |
---|
| 434 | DO_2D_00_00 |
---|
| 435 | en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1) |
---|
| 436 | END_2D |
---|
[12991] | 437 | DO_3DS_00_00( jpk-2, 2, -1 ) ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
[12377] | 438 | en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk) |
---|
| 439 | END_3D |
---|
[12991] | 440 | DO_3D_00_00( 2, jpkm1 ) ! set the minimum value of tke |
---|
[12377] | 441 | en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk) |
---|
| 442 | END_3D |
---|
[9019] | 443 | ! |
---|
[1492] | 444 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 445 | ! ! TKE due to surface and internal wave breaking |
---|
| 446 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[6140] | 447 | !!gm BUG : in the exp remove the depth of ssh !!! |
---|
[12377] | 448 | !!gm i.e. use gde3w in argument (gdepw(:,:,:,Kmm)) |
---|
[6140] | 449 | |
---|
| 450 | |
---|
[2528] | 451 | IF( nn_etau == 1 ) THEN !* penetration below the mixed layer (rn_efr fraction) |
---|
[12377] | 452 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 453 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) ) & |
---|
| 454 | & * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
| 455 | END_3D |
---|
[2528] | 456 | ELSEIF( nn_etau == 2 ) THEN !* act only at the base of the mixed layer (jk=nmln) (rn_efr fraction) |
---|
[12377] | 457 | DO_2D_00_00 |
---|
| 458 | jk = nmln(ji,jj) |
---|
| 459 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) ) & |
---|
| 460 | & * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
| 461 | END_2D |
---|
[2528] | 462 | ELSEIF( nn_etau == 3 ) THEN !* penetration belox the mixed layer (HF variability) |
---|
[12377] | 463 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 464 | ztx2 = utau(ji-1,jj ) + utau(ji,jj) |
---|
| 465 | zty2 = vtau(ji ,jj-1) + vtau(ji,jj) |
---|
| 466 | ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1) ! module of the mean stress |
---|
| 467 | zdif = taum(ji,jj) - ztau ! mean of modulus - modulus of the mean |
---|
| 468 | zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add ) ! apply some modifications... |
---|
| 469 | en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) ) & |
---|
| 470 | & * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
| 471 | END_3D |
---|
[1239] | 472 | ENDIF |
---|
[1492] | 473 | ! |
---|
[1239] | 474 | END SUBROUTINE tke_tke |
---|
| 475 | |
---|
[1492] | 476 | |
---|
[12377] | 477 | SUBROUTINE tke_avn( Kbb, Kmm, p_avm, p_avt ) |
---|
[1239] | 478 | !!---------------------------------------------------------------------- |
---|
[1492] | 479 | !! *** ROUTINE tke_avn *** |
---|
[1239] | 480 | !! |
---|
[1492] | 481 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
| 482 | !! |
---|
| 483 | !! ** Method : At this stage, en, the now TKE, is known (computed in |
---|
| 484 | !! the tke_tke routine). First, the now mixing lenth is |
---|
| 485 | !! computed from en and the strafification (N^2), then the mixings |
---|
| 486 | !! coefficients are computed. |
---|
| 487 | !! - Mixing length : a first evaluation of the mixing lengh |
---|
| 488 | !! scales is: |
---|
| 489 | !! mxl = sqrt(2*en) / N |
---|
| 490 | !! where N is the brunt-vaisala frequency, with a minimum value set |
---|
[2528] | 491 | !! to rmxl_min (rn_mxl0) in the interior (surface) ocean. |
---|
[1492] | 492 | !! The mixing and dissipative length scale are bound as follow : |
---|
| 493 | !! nn_mxl=0 : mxl bounded by the distance to surface and bottom. |
---|
| 494 | !! zmxld = zmxlm = mxl |
---|
| 495 | !! nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl |
---|
| 496 | !! nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is |
---|
| 497 | !! less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl |
---|
| 498 | !! nn_mxl=3 : mxl is bounded from the surface to the bottom usings |
---|
| 499 | !! |d/dz(xml)|<1 to obtain lup, and from the bottom to |
---|
| 500 | !! the surface to obtain ldown. the resulting length |
---|
| 501 | !! scales are: |
---|
| 502 | !! zmxld = sqrt( lup * ldown ) |
---|
| 503 | !! zmxlm = min ( lup , ldown ) |
---|
| 504 | !! - Vertical eddy viscosity and diffusivity: |
---|
| 505 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
| 506 | !! avt = max( avmb, pdlr * avm ) |
---|
| 507 | !! with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise. |
---|
| 508 | !! |
---|
[9019] | 509 | !! ** Action : - avt, avm : now vertical eddy diffusivity and viscosity (w-point) |
---|
[1239] | 510 | !!---------------------------------------------------------------------- |
---|
[9019] | 511 | USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d ! ocean vertical physics |
---|
| 512 | !! |
---|
[12377] | 513 | INTEGER , INTENT(in ) :: Kbb, Kmm ! ocean time level indices |
---|
[9019] | 514 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_avm, p_avt ! vertical eddy viscosity & diffusivity (w-points) |
---|
| 515 | ! |
---|
[2715] | 516 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[9019] | 517 | REAL(wp) :: zrn2, zraug, zcoef, zav ! local scalars |
---|
| 518 | REAL(wp) :: zdku, zdkv, zsqen ! - - |
---|
| 519 | REAL(wp) :: zemxl, zemlm, zemlp ! - - |
---|
| 520 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: zmxlm, zmxld ! 3D workspace |
---|
[1239] | 521 | !!-------------------------------------------------------------------- |
---|
[3294] | 522 | ! |
---|
[1492] | 523 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 524 | ! ! Mixing length |
---|
| 525 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 526 | ! |
---|
| 527 | ! !* Buoyancy length scale: l=sqrt(2*e/n**2) |
---|
| 528 | ! |
---|
[5120] | 529 | ! initialisation of interior minimum value (avoid a 2d loop with mikt) |
---|
[7753] | 530 | zmxlm(:,:,:) = rmxl_min |
---|
| 531 | zmxld(:,:,:) = rmxl_min |
---|
[5120] | 532 | ! |
---|
[12991] | 533 | IF(ln_sdw .AND. ln_mxhsw) THEN |
---|
| 534 | zmxlm(:,:,1)= vkarmn * MAX ( 1.6 * hsw(:,:) , 0.02 ) ! surface mixing length = F(wave height) |
---|
| 535 | ! from terray et al 1999 and mellor and blumberg 2004 it should be 0.85 and not 1.6 |
---|
| 536 | zcoef = vkarmn * ( (rn_ediff*rn_ediss)**0.25 ) / rn_ediff |
---|
| 537 | zmxlm(:,:,1)= zcoef * MAX ( 1.6 * hsw(:,:) , 0.02 ) ! surface mixing length = F(wave height) |
---|
| 538 | ELSE |
---|
| 539 | IF( ln_mxl0 ) THEN ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rho0*g) |
---|
| 540 | zraug = vkarmn * 2.e5_wp / ( rho0 * grav ) |
---|
| 541 | DO_2D_00_00 |
---|
| 542 | zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) ) |
---|
| 543 | END_2D |
---|
| 544 | ELSE |
---|
| 545 | zmxlm(:,:,1) = rn_mxl0 |
---|
| 546 | ENDIF |
---|
[1239] | 547 | ENDIF |
---|
| 548 | ! |
---|
[12377] | 549 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 550 | zrn2 = MAX( rn2(ji,jj,jk), rsmall ) |
---|
| 551 | zmxlm(ji,jj,jk) = MAX( rmxl_min, SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) ) |
---|
| 552 | END_3D |
---|
[1492] | 553 | ! |
---|
| 554 | ! !* Physical limits for the mixing length |
---|
| 555 | ! |
---|
[7753] | 556 | zmxld(:,:, 1 ) = zmxlm(:,:,1) ! surface set to the minimum value |
---|
| 557 | zmxld(:,:,jpk) = rmxl_min ! last level set to the minimum value |
---|
[1492] | 558 | ! |
---|
[1239] | 559 | SELECT CASE ( nn_mxl ) |
---|
| 560 | ! |
---|
[5836] | 561 | !!gm Not sure of that coding for ISF.... |
---|
[12377] | 562 | ! where wmask = 0 set zmxlm == e3w(:,:,:,Kmm) |
---|
[1239] | 563 | CASE ( 0 ) ! bounded by the distance to surface and bottom |
---|
[12377] | 564 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 565 | zemxl = MIN( gdepw(ji,jj,jk,Kmm) - gdepw(ji,jj,mikt(ji,jj),Kmm), zmxlm(ji,jj,jk), & |
---|
| 566 | & gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) - gdepw(ji,jj,jk,Kmm) ) |
---|
| 567 | ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj) |
---|
| 568 | zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk)) |
---|
| 569 | zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk)) |
---|
| 570 | END_3D |
---|
[1239] | 571 | ! |
---|
| 572 | CASE ( 1 ) ! bounded by the vertical scale factor |
---|
[12377] | 573 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 574 | zemxl = MIN( e3w(ji,jj,jk,Kmm), zmxlm(ji,jj,jk) ) |
---|
| 575 | zmxlm(ji,jj,jk) = zemxl |
---|
| 576 | zmxld(ji,jj,jk) = zemxl |
---|
| 577 | END_3D |
---|
[1239] | 578 | ! |
---|
| 579 | CASE ( 2 ) ! |dk[xml]| bounded by e3t : |
---|
[12377] | 580 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 581 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) ) |
---|
| 582 | END_3D |
---|
| 583 | DO_3DS_00_00( jpkm1, 2, -1 ) |
---|
| 584 | zemxl = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) ) |
---|
| 585 | zmxlm(ji,jj,jk) = zemxl |
---|
| 586 | zmxld(ji,jj,jk) = zemxl |
---|
| 587 | END_3D |
---|
[1239] | 588 | ! |
---|
| 589 | CASE ( 3 ) ! lup and ldown, |dk[xml]| bounded by e3t : |
---|
[12377] | 590 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 591 | zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) ) |
---|
| 592 | END_3D |
---|
| 593 | DO_3DS_00_00( jpkm1, 2, -1 ) |
---|
| 594 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) ) |
---|
| 595 | END_3D |
---|
| 596 | DO_3D_00_00( 2, jpkm1 ) |
---|
| 597 | zemlm = MIN ( zmxld(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
| 598 | zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) ) |
---|
| 599 | zmxlm(ji,jj,jk) = zemlm |
---|
| 600 | zmxld(ji,jj,jk) = zemlp |
---|
| 601 | END_3D |
---|
[1239] | 602 | ! |
---|
| 603 | END SELECT |
---|
[1492] | 604 | ! |
---|
| 605 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[9019] | 606 | ! ! Vertical eddy viscosity and diffusivity (avm and avt) |
---|
[1492] | 607 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[12377] | 608 | DO_3D_00_00( 1, jpkm1 ) |
---|
| 609 | zsqen = SQRT( en(ji,jj,jk) ) |
---|
| 610 | zav = rn_ediff * zmxlm(ji,jj,jk) * zsqen |
---|
| 611 | p_avm(ji,jj,jk) = MAX( zav, avmb(jk) ) * wmask(ji,jj,jk) |
---|
| 612 | p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk) |
---|
| 613 | dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk) |
---|
| 614 | END_3D |
---|
[1492] | 615 | ! |
---|
| 616 | ! |
---|
| 617 | IF( nn_pdl == 1 ) THEN !* Prandtl number case: update avt |
---|
[12377] | 618 | DO_3D_00_00( 2, jpkm1 ) |
---|
[12698] | 619 | p_avt(ji,jj,jk) = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk) |
---|
[12377] | 620 | END_3D |
---|
[1239] | 621 | ENDIF |
---|
[9019] | 622 | ! |
---|
[12377] | 623 | IF(sn_cfctl%l_prtctl) THEN |
---|
[9440] | 624 | CALL prt_ctl( tab3d_1=en , clinfo1=' tke - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk) |
---|
| 625 | CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke - m: ', kdim=jpk ) |
---|
[1239] | 626 | ENDIF |
---|
| 627 | ! |
---|
[1492] | 628 | END SUBROUTINE tke_avn |
---|
[1239] | 629 | |
---|
[1492] | 630 | |
---|
[12377] | 631 | SUBROUTINE zdf_tke_init( Kmm ) |
---|
[1239] | 632 | !!---------------------------------------------------------------------- |
---|
[2528] | 633 | !! *** ROUTINE zdf_tke_init *** |
---|
[1239] | 634 | !! |
---|
| 635 | !! ** Purpose : Initialization of the vertical eddy diffivity and |
---|
[1492] | 636 | !! viscosity when using a tke turbulent closure scheme |
---|
[1239] | 637 | !! |
---|
[1601] | 638 | !! ** Method : Read the namzdf_tke namelist and check the parameters |
---|
[1492] | 639 | !! called at the first timestep (nit000) |
---|
[1239] | 640 | !! |
---|
[1601] | 641 | !! ** input : Namlist namzdf_tke |
---|
[1239] | 642 | !! |
---|
| 643 | !! ** Action : Increase by 1 the nstop flag is setting problem encounter |
---|
| 644 | !!---------------------------------------------------------------------- |
---|
[9019] | 645 | USE zdf_oce , ONLY : ln_zdfiwm ! Internal Wave Mixing flag |
---|
| 646 | !! |
---|
[12377] | 647 | INTEGER, INTENT(in) :: Kmm ! time level index |
---|
| 648 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 649 | INTEGER :: ios |
---|
[1239] | 650 | !! |
---|
[9019] | 651 | NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin , & |
---|
| 652 | & rn_emin0, rn_bshear, nn_mxl , ln_mxl0 , & |
---|
| 653 | & rn_mxl0 , nn_pdl , ln_drg , ln_lc , rn_lc, & |
---|
[12991] | 654 | & nn_etau , nn_htau , rn_efr , rn_eice , & |
---|
| 655 | & nn_bc_surf, nn_bc_bot, ln_mxhsw |
---|
[1239] | 656 | !!---------------------------------------------------------------------- |
---|
[2715] | 657 | ! |
---|
[4147] | 658 | READ ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901) |
---|
[11536] | 659 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' ) |
---|
[4147] | 660 | |
---|
| 661 | READ ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 ) |
---|
[11536] | 662 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' ) |
---|
[4624] | 663 | IF(lwm) WRITE ( numond, namzdf_tke ) |
---|
[2715] | 664 | ! |
---|
[2528] | 665 | ri_cri = 2._wp / ( 2._wp + rn_ediss / rn_ediff ) ! resulting critical Richardson number |
---|
[2715] | 666 | ! |
---|
[1492] | 667 | IF(lwp) THEN !* Control print |
---|
[1239] | 668 | WRITE(numout,*) |
---|
[2528] | 669 | WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation' |
---|
| 670 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
[1601] | 671 | WRITE(numout,*) ' Namelist namzdf_tke : set tke mixing parameters' |
---|
[1705] | 672 | WRITE(numout,*) ' coef. to compute avt rn_ediff = ', rn_ediff |
---|
| 673 | WRITE(numout,*) ' Kolmogoroff dissipation coef. rn_ediss = ', rn_ediss |
---|
| 674 | WRITE(numout,*) ' tke surface input coef. rn_ebb = ', rn_ebb |
---|
| 675 | WRITE(numout,*) ' minimum value of tke rn_emin = ', rn_emin |
---|
| 676 | WRITE(numout,*) ' surface minimum value of tke rn_emin0 = ', rn_emin0 |
---|
[9019] | 677 | WRITE(numout,*) ' prandl number flag nn_pdl = ', nn_pdl |
---|
[1705] | 678 | WRITE(numout,*) ' background shear (>0) rn_bshear = ', rn_bshear |
---|
| 679 | WRITE(numout,*) ' mixing length type nn_mxl = ', nn_mxl |
---|
[9019] | 680 | WRITE(numout,*) ' surface mixing length = F(stress) or not ln_mxl0 = ', ln_mxl0 |
---|
| 681 | WRITE(numout,*) ' surface mixing length minimum value rn_mxl0 = ', rn_mxl0 |
---|
| 682 | WRITE(numout,*) ' top/bottom friction forcing flag ln_drg = ', ln_drg |
---|
| 683 | WRITE(numout,*) ' Langmuir cells parametrization ln_lc = ', ln_lc |
---|
| 684 | WRITE(numout,*) ' coef to compute vertical velocity of LC rn_lc = ', rn_lc |
---|
[1705] | 685 | WRITE(numout,*) ' test param. to add tke induced by wind nn_etau = ', nn_etau |
---|
[9019] | 686 | WRITE(numout,*) ' type of tke penetration profile nn_htau = ', nn_htau |
---|
| 687 | WRITE(numout,*) ' fraction of TKE that penetrates rn_efr = ', rn_efr |
---|
[9546] | 688 | WRITE(numout,*) ' below sea-ice: =0 ON rn_eice = ', rn_eice |
---|
| 689 | WRITE(numout,*) ' =4 OFF when ice fraction > 1/4 ' |
---|
[9019] | 690 | IF( ln_drg ) THEN |
---|
[9169] | 691 | WRITE(numout,*) |
---|
[9019] | 692 | WRITE(numout,*) ' Namelist namdrg_top/_bot: used values:' |
---|
| 693 | WRITE(numout,*) ' top ocean cavity roughness (m) rn_z0(_top)= ', r_z0_top |
---|
| 694 | WRITE(numout,*) ' Bottom seafloor roughness (m) rn_z0(_bot)= ', r_z0_bot |
---|
| 695 | ENDIF |
---|
| 696 | WRITE(numout,*) |
---|
[9190] | 697 | WRITE(numout,*) ' ==>>> critical Richardson nb with your parameters ri_cri = ', ri_cri |
---|
[9019] | 698 | WRITE(numout,*) |
---|
[1239] | 699 | ENDIF |
---|
[2715] | 700 | ! |
---|
[9019] | 701 | IF( ln_zdfiwm ) THEN ! Internal wave-driven mixing |
---|
| 702 | rn_emin = 1.e-10_wp ! specific values of rn_emin & rmxl_min are used |
---|
| 703 | rmxl_min = 1.e-03_wp ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s) |
---|
[9190] | 704 | IF(lwp) WRITE(numout,*) ' ==>>> Internal wave-driven mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3' |
---|
[9019] | 705 | ELSE ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s) |
---|
| 706 | rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) ) ! resulting minimum length to recover molecular viscosity |
---|
[9190] | 707 | IF(lwp) WRITE(numout,*) ' ==>>> minimum mixing length with your parameters rmxl_min = ', rmxl_min |
---|
[9019] | 708 | ENDIF |
---|
| 709 | ! |
---|
[2715] | 710 | ! ! allocate tke arrays |
---|
| 711 | IF( zdf_tke_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' ) |
---|
| 712 | ! |
---|
[1492] | 713 | ! !* Check of some namelist values |
---|
[4990] | 714 | IF( nn_mxl < 0 .OR. nn_mxl > 3 ) CALL ctl_stop( 'bad flag: nn_mxl is 0, 1 or 2 ' ) |
---|
| 715 | IF( nn_pdl < 0 .OR. nn_pdl > 1 ) CALL ctl_stop( 'bad flag: nn_pdl is 0 or 1 ' ) |
---|
| 716 | IF( nn_htau < 0 .OR. nn_htau > 1 ) CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' ) |
---|
[5407] | 717 | IF( nn_etau == 3 .AND. .NOT. ln_cpl ) CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) |
---|
[9019] | 718 | ! |
---|
[2528] | 719 | IF( ln_mxl0 ) THEN |
---|
[9169] | 720 | IF(lwp) WRITE(numout,*) |
---|
[9190] | 721 | IF(lwp) WRITE(numout,*) ' ==>>> use a surface mixing length = F(stress) : set rn_mxl0 = rmxl_min' |
---|
[2528] | 722 | rn_mxl0 = rmxl_min |
---|
| 723 | ENDIF |
---|
| 724 | |
---|
[12377] | 725 | IF( nn_etau == 2 ) CALL zdf_mxl( nit000, Kmm ) ! Initialization of nmln |
---|
[1239] | 726 | |
---|
[1492] | 727 | ! !* depth of penetration of surface tke |
---|
| 728 | IF( nn_etau /= 0 ) THEN |
---|
[1601] | 729 | SELECT CASE( nn_htau ) ! Choice of the depth of penetration |
---|
[2528] | 730 | CASE( 0 ) ! constant depth penetration (here 10 meters) |
---|
[7753] | 731 | htau(:,:) = 10._wp |
---|
[2528] | 732 | CASE( 1 ) ! F(latitude) : 0.5m to 30m poleward of 40 degrees |
---|
[7753] | 733 | htau(:,:) = MAX( 0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) ) ) |
---|
[1492] | 734 | END SELECT |
---|
| 735 | ENDIF |
---|
[9019] | 736 | ! !* read or initialize all required files |
---|
| 737 | CALL tke_rst( nit000, 'READ' ) ! (en, avt_k, avm_k, dissl) |
---|
[1239] | 738 | ! |
---|
[9367] | 739 | IF( lwxios ) THEN |
---|
| 740 | CALL iom_set_rstw_var_active('en') |
---|
| 741 | CALL iom_set_rstw_var_active('avt_k') |
---|
| 742 | CALL iom_set_rstw_var_active('avm_k') |
---|
| 743 | CALL iom_set_rstw_var_active('dissl') |
---|
| 744 | ENDIF |
---|
[2528] | 745 | END SUBROUTINE zdf_tke_init |
---|
[1239] | 746 | |
---|
| 747 | |
---|
[1531] | 748 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
[9019] | 749 | !!--------------------------------------------------------------------- |
---|
| 750 | !! *** ROUTINE tke_rst *** |
---|
| 751 | !! |
---|
| 752 | !! ** Purpose : Read or write TKE file (en) in restart file |
---|
| 753 | !! |
---|
| 754 | !! ** Method : use of IOM library |
---|
| 755 | !! if the restart does not contain TKE, en is either |
---|
| 756 | !! set to rn_emin or recomputed |
---|
| 757 | !!---------------------------------------------------------------------- |
---|
| 758 | USE zdf_oce , ONLY : en, avt_k, avm_k ! ocean vertical physics |
---|
| 759 | !! |
---|
| 760 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 761 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
| 762 | ! |
---|
| 763 | INTEGER :: jit, jk ! dummy loop indices |
---|
| 764 | INTEGER :: id1, id2, id3, id4 ! local integers |
---|
| 765 | !!---------------------------------------------------------------------- |
---|
| 766 | ! |
---|
| 767 | IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise |
---|
| 768 | ! ! --------------- |
---|
| 769 | IF( ln_rstart ) THEN !* Read the restart file |
---|
| 770 | id1 = iom_varid( numror, 'en' , ldstop = .FALSE. ) |
---|
| 771 | id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. ) |
---|
| 772 | id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. ) |
---|
| 773 | id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. ) |
---|
| 774 | ! |
---|
| 775 | IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN ! fields exist |
---|
[9367] | 776 | CALL iom_get( numror, jpdom_autoglo, 'en' , en , ldxios = lrxios ) |
---|
| 777 | CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k, ldxios = lrxios ) |
---|
| 778 | CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k, ldxios = lrxios ) |
---|
| 779 | CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, ldxios = lrxios ) |
---|
[9019] | 780 | ELSE ! start TKE from rest |
---|
[9169] | 781 | IF(lwp) WRITE(numout,*) |
---|
[9190] | 782 | IF(lwp) WRITE(numout,*) ' ==>>> previous run without TKE scheme, set en to background values' |
---|
[9019] | 783 | en (:,:,:) = rn_emin * wmask(:,:,:) |
---|
| 784 | dissl(:,:,:) = 1.e-12_wp |
---|
| 785 | ! avt_k, avm_k already set to the background value in zdf_phy_init |
---|
| 786 | ENDIF |
---|
| 787 | ELSE !* Start from rest |
---|
[9169] | 788 | IF(lwp) WRITE(numout,*) |
---|
[9190] | 789 | IF(lwp) WRITE(numout,*) ' ==>>> start from rest: set en to the background value' |
---|
[9019] | 790 | en (:,:,:) = rn_emin * wmask(:,:,:) |
---|
| 791 | dissl(:,:,:) = 1.e-12_wp |
---|
| 792 | ! avt_k, avm_k already set to the background value in zdf_phy_init |
---|
| 793 | ENDIF |
---|
| 794 | ! |
---|
| 795 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file |
---|
| 796 | ! ! ------------------- |
---|
[9169] | 797 | IF(lwp) WRITE(numout,*) '---- tke_rst ----' |
---|
[9367] | 798 | IF( lwxios ) CALL iom_swap( cwxios_context ) |
---|
| 799 | CALL iom_rstput( kt, nitrst, numrow, 'en' , en , ldxios = lwxios ) |
---|
| 800 | CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k, ldxios = lwxios ) |
---|
| 801 | CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k, ldxios = lwxios ) |
---|
| 802 | CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl, ldxios = lwxios ) |
---|
| 803 | IF( lwxios ) CALL iom_swap( cxios_context ) |
---|
[9019] | 804 | ! |
---|
| 805 | ENDIF |
---|
| 806 | ! |
---|
[1531] | 807 | END SUBROUTINE tke_rst |
---|
[1239] | 808 | |
---|
| 809 | !!====================================================================== |
---|
[1531] | 810 | END MODULE zdftke |
---|