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