1 | MODULE zdftke_old |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE zdftke_old *** |
---|
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 zdf_tke_init routine |
---|
16 | !! - ! 2002-08 (G. Madec) rn_cri and Free form, F90 |
---|
17 | !! - ! 2004-10 (C. Ethe ) 1D configuration |
---|
18 | !! 2.0 ! 2006-07 (S. Masson) distributed restart using iom |
---|
19 | !! 3.0 ! 2008-05 (C. Ethe, G.Madec) : update TKE physics: |
---|
20 | !! - tke penetration (wind steering) |
---|
21 | !! - suface condition for tke & mixing length |
---|
22 | !! - Langmuir cells |
---|
23 | !! - ! 2008-05 (J.-M. Molines, G. Madec) 2D form of avtb |
---|
24 | !! - ! 2008-06 (G. Madec) style + DOCTOR name for namelist parameters |
---|
25 | !!---------------------------------------------------------------------- |
---|
26 | #if defined key_zdftke_old || defined key_esopa |
---|
27 | !!---------------------------------------------------------------------- |
---|
28 | !! 'key_zdftke_old' TKE vertical physics |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | !!---------------------------------------------------------------------- |
---|
31 | !! zdf_tke_old : update momentum and tracer Kz from a tke scheme |
---|
32 | !! zdf_tke_init : initialization, namelist read, and parameters control |
---|
33 | !! tke_rst : read/write tke restart in ocean restart file |
---|
34 | !!---------------------------------------------------------------------- |
---|
35 | USE oce ! ocean dynamics and active tracers |
---|
36 | USE dom_oce ! ocean space and time domain |
---|
37 | USE zdf_oce ! ocean vertical physics |
---|
38 | USE sbc_oce ! surface boundary condition: ocean |
---|
39 | USE phycst ! physical constants |
---|
40 | USE zdfmxl ! mixed layer |
---|
41 | USE restart ! only for lrst_oce |
---|
42 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
43 | USE prtctl ! Print control |
---|
44 | USE in_out_manager ! I/O manager |
---|
45 | USE iom ! I/O manager library |
---|
46 | |
---|
47 | IMPLICIT NONE |
---|
48 | PRIVATE |
---|
49 | |
---|
50 | PUBLIC zdf_tke_old ! routine called in step module |
---|
51 | |
---|
52 | LOGICAL , PUBLIC, PARAMETER :: lk_zdftke_old = .TRUE. !: TKE vertical mixing flag |
---|
53 | REAL(wp), PUBLIC :: eboost !: multiplicative coeff of the shear product. |
---|
54 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: en !: now turbulent kinetic energy |
---|
55 | # if defined key_vectopt_memory |
---|
56 | ! !!! key_vectopt_memory |
---|
57 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: etmean !: coefficient used for horizontal smoothing |
---|
58 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: eumean, evmean !: at t-, u- and v-points |
---|
59 | # endif |
---|
60 | #if defined key_c1d |
---|
61 | ! !!! 1D cfg only |
---|
62 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: e_dis, e_mix !: dissipation and mixing turbulent lengh scales |
---|
63 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: e_pdl, e_ric !: prandl and local Richardson numbers |
---|
64 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: hlc !: save finite Langmuir Circulation depth |
---|
65 | #endif |
---|
66 | |
---|
67 | ! !!! ** Namelist namzdf_tke ** |
---|
68 | LOGICAL :: ln_mxl0 = .FALSE. ! mixing length scale surface value as function of wind stress or not |
---|
69 | LOGICAL :: ln_lc = .FALSE. ! Langmuir cells (LC) as a source term of TKE or not |
---|
70 | INTEGER :: nn_mxl = 2 ! type of mixing length (=0/1/2/3) |
---|
71 | INTEGER :: nn_pdl = 1 ! Prandtl number or not (ratio avt/avm) (=0/1) |
---|
72 | REAL(wp) :: rn_ediff = 0.1_wp ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e) |
---|
73 | REAL(wp) :: rn_ediss = 0.7_wp ! coefficient of the Kolmogoroff dissipation |
---|
74 | REAL(wp) :: rn_ebb = 3.75_wp ! coefficient of the surface input of tke |
---|
75 | REAL(wp) :: rn_emin = 0.7071e-6_wp ! minimum value of tke (m2/s2) |
---|
76 | REAL(wp) :: rn_emin0 = 1.e-4_wp ! surface minimum value of tke (m2/s2) |
---|
77 | REAL(wp) :: rn_bshear= 1.e-20_wp ! background shear (>0) (Not used in old TKE) |
---|
78 | INTEGER :: nn_etau = 0 ! type of depth penetration of surface tke (=0/1/2) |
---|
79 | INTEGER :: nn_htau = 0 ! type of tke profile of penetration (=0/1) |
---|
80 | REAL(wp) :: rn_lmin0 = 0.4_wp ! surface min value of mixing length |
---|
81 | REAL(wp) :: rn_lmin = 0.1_wp ! interior min value of mixing length |
---|
82 | REAL(wp) :: rn_efr = 1.0_wp ! fraction of TKE surface value which penetrates in the ocean |
---|
83 | REAL(wp) :: rn_addhft= 0.0_wp ! add offset applied to HF tau (Not used in old TKE) |
---|
84 | REAL(wp) :: rn_sclhft= 1.0_wp ! scale factor applied to HF tau (Not used in old TKE) |
---|
85 | REAL(wp) :: rn_lc = 0.15_wp ! coef to compute vertical velocity of Langmuir cells |
---|
86 | |
---|
87 | ! !! ** old namelist value: now hard coded ** |
---|
88 | INTEGER :: nn_ave = 1 ! horizontal average or not on avt, avmu, avmv (=0/1) |
---|
89 | REAL(wp) :: rn_efave = 1._wp ! coefficient for ave : ave=rn_efave*avm |
---|
90 | REAL(wp) :: rn_cri = 2._wp / 9._wp ! critic Richardson number |
---|
91 | |
---|
92 | !! * Substitutions |
---|
93 | # include "domzgr_substitute.h90" |
---|
94 | # include "vectopt_loop_substitute.h90" |
---|
95 | !!---------------------------------------------------------------------- |
---|
96 | !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008) |
---|
97 | !! $Id$ |
---|
98 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
99 | !!---------------------------------------------------------------------- |
---|
100 | |
---|
101 | CONTAINS |
---|
102 | |
---|
103 | SUBROUTINE zdf_tke_old( kt ) |
---|
104 | !!---------------------------------------------------------------------- |
---|
105 | !! *** ROUTINE zdf_tke_old *** |
---|
106 | !! |
---|
107 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
108 | !! coefficients using a 1.5 turbulent closure scheme. |
---|
109 | !! |
---|
110 | !! ** Method : The time evolution of the turbulent kinetic energy |
---|
111 | !! (tke) is computed from a prognostic equation : |
---|
112 | !! d(en)/dt = eboost eav (d(u)/dz)**2 ! shear production |
---|
113 | !! + d( rn_efave eav d(en)/dz )/dz ! diffusion of tke |
---|
114 | !! + grav/rau0 pdl eav d(rau)/dz ! stratif. destruc. |
---|
115 | !! - rn_ediss / emxl en**(2/3) ! dissipation |
---|
116 | !! with the boundary conditions: |
---|
117 | !! surface: en = max( rn_emin0, rn_ebb sqrt(utau^2 + vtau^2) ) |
---|
118 | !! bottom : en = rn_emin |
---|
119 | !! -1- The dissipation and mixing turbulent lengh scales are computed |
---|
120 | !! from the usual diagnostic buoyancy length scale: |
---|
121 | !! mxl= sqrt(2*en)/N where N is the brunt-vaisala frequency |
---|
122 | !! with mxl = rn_lmin at the bottom minimum value of 0.4 |
---|
123 | !! Four cases : |
---|
124 | !! nn_mxl=0 : mxl bounded by the distance to surface and bottom. |
---|
125 | !! zmxld = zmxlm = mxl |
---|
126 | !! nn_mxl=1 : mxl bounded by the vertical scale factor. |
---|
127 | !! zmxld = zmxlm = mxl |
---|
128 | !! nn_mxl=2 : mxl bounded such that the vertical derivative of mxl |
---|
129 | !! is less than 1 (|d/dz(xml)|<1). |
---|
130 | !! zmxld = zmxlm = mxl |
---|
131 | !! nn_mxl=3 : lup = mxl bounded using |d/dz(xml)|<1 from the surface |
---|
132 | !! to the bottom |
---|
133 | !! ldown = mxl bounded using |d/dz(xml)|<1 from the bottom |
---|
134 | !! to the surface |
---|
135 | !! zmxld = sqrt (lup*ldown) ; zmxlm = min(lup,ldown) |
---|
136 | !! -2- Compute the now Turbulent kinetic energy. The time differencing |
---|
137 | !! is implicit for vertical diffusion term, linearized for kolmo- |
---|
138 | !! goroff dissipation term, and explicit forward for both buoyancy |
---|
139 | !! and dynamic production terms. Thus a tridiagonal linear system is |
---|
140 | !! solved. |
---|
141 | !! Note that - the shear production is multiplied by eboost in order |
---|
142 | !! to set the critic richardson number to rn_cri (namelist parameter) |
---|
143 | !! - the destruction by stratification term is multiplied |
---|
144 | !! by the Prandtl number (defined by an empirical funtion of the local |
---|
145 | !! Richardson number) if nn_pdl=1 (namelist parameter) |
---|
146 | !! coefficient (zesh2): |
---|
147 | !! -3- Compute the now vertical eddy vicosity and diffusivity |
---|
148 | !! coefficients from en (before the time stepping) and zmxlm: |
---|
149 | !! avm = max( avtb, rn_ediff*zmxlm*en^1/2 ) |
---|
150 | !! avt = max( avmb, pdl*avm ) (pdl=1 if nn_pdl=0) |
---|
151 | !! eav = max( avmb, avm ) |
---|
152 | !! avt and avm are horizontally averaged to avoid numerical insta- |
---|
153 | !! bilities. |
---|
154 | !! N.B. The computation is done from jk=2 to jpkm1 except for |
---|
155 | !! en. Surface value of avt avmu avmv are set once a time to |
---|
156 | !! their background value in routine zdf_tke_init. |
---|
157 | !! |
---|
158 | !! ** Action : compute en (now turbulent kinetic energy) |
---|
159 | !! update avt, avmu, avmv (before vertical eddy coef.) |
---|
160 | !! |
---|
161 | !! References : Gaspar et al., JGR, 1990, |
---|
162 | !! Blanke and Delecluse, JPO, 1991 |
---|
163 | !! Mellor and Blumberg, JPO 2004 |
---|
164 | !! Axell, JGR, 2002 |
---|
165 | !!---------------------------------------------------------------------- |
---|
166 | USE oce, zwd => ua ! use ua as workspace |
---|
167 | USE oce, zmxlm => va ! use va as workspace |
---|
168 | USE oce, zmxld => ta ! use ta as workspace |
---|
169 | USE oce, ztkelc => sa ! use sa as workspace |
---|
170 | ! |
---|
171 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
172 | ! |
---|
173 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
174 | REAL(wp) :: zbbrau, zrn2, zesurf ! temporary scalars |
---|
175 | REAL(wp) :: zfact1, ztx2, zdku ! - - |
---|
176 | REAL(wp) :: zfact2, zty2, zdkv ! - - |
---|
177 | REAL(wp) :: zfact3, zcoef, zcof, zav ! - - |
---|
178 | REAL(wp) :: zsh2, zpdl, zri, zsqen, zesh2 ! - - |
---|
179 | REAL(wp) :: zemxl, zemlm, zemlp ! - - |
---|
180 | REAL(wp) :: zraug, zus, zwlc, zind ! - - |
---|
181 | INTEGER , DIMENSION(jpi,jpj) :: imlc ! 2D workspace |
---|
182 | REAL(wp), DIMENSION(jpi,jpj) :: zhtau ! - - |
---|
183 | REAL(wp), DIMENSION(jpi,jpj) :: zhlc ! - - |
---|
184 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: zpelc ! 3D workspace |
---|
185 | !!-------------------------------------------------------------------- |
---|
186 | |
---|
187 | IF( kt == nit000 ) CALL zdf_tke_init ! Initialization (first time-step only) |
---|
188 | |
---|
189 | ! ! Local constant initialization |
---|
190 | zbbrau = .5 * rn_ebb / rau0 |
---|
191 | zfact1 = -.5 * rdt * rn_efave |
---|
192 | zfact2 = 1.5 * rdt * rn_ediss |
---|
193 | zfact3 = 0.5 * rdt * rn_ediss |
---|
194 | |
---|
195 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
196 | ! I. Mixing length |
---|
197 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
198 | |
---|
199 | ! Buoyancy length scale: l=sqrt(2*e/n**2) |
---|
200 | ! --------------------- |
---|
201 | IF( ln_mxl0 ) THEN ! surface mixing length = F(stress) : l=vkarmn*2.e5*sqrt(utau^2 + vtau^2)/(rau0*g) |
---|
202 | !!gm this should be useless |
---|
203 | zmxlm(:,:,1) = 0.e0 |
---|
204 | !!gm end |
---|
205 | zraug = 0.5 * vkarmn * 2.e5 / ( rau0 * grav ) |
---|
206 | DO jj = 2, jpjm1 |
---|
207 | !CDIR NOVERRCHK |
---|
208 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
209 | ztx2 = utau(ji-1,jj ) + utau(ji,jj) |
---|
210 | zty2 = vtau(ji ,jj-1) + vtau(ji,jj) |
---|
211 | zmxlm(ji,jj,1) = MAX( rn_lmin0, zraug * SQRT( ztx2 * ztx2 + zty2 * zty2 ) ) |
---|
212 | END DO |
---|
213 | END DO |
---|
214 | ELSE ! surface set to the minimum value |
---|
215 | zmxlm(:,:,1) = rn_lmin0 |
---|
216 | ENDIF |
---|
217 | zmxlm(:,:,jpk) = rn_lmin ! bottom set to the interior minium value |
---|
218 | ! |
---|
219 | !CDIR NOVERRCHK |
---|
220 | DO jk = 2, jpkm1 ! interior value : l=sqrt(2*e/n**2) |
---|
221 | !CDIR NOVERRCHK |
---|
222 | DO jj = 2, jpjm1 |
---|
223 | !CDIR NOVERRCHK |
---|
224 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
225 | zrn2 = MAX( rn2b(ji,jj,jk), rsmall ) |
---|
226 | zmxlm(ji,jj,jk) = MAX( rn_lmin, SQRT( 2. * en(ji,jj,jk) / zrn2 ) ) |
---|
227 | END DO |
---|
228 | END DO |
---|
229 | END DO |
---|
230 | |
---|
231 | ! Physical limits for the mixing length |
---|
232 | ! ------------------------------------- |
---|
233 | zmxld(:,:, 1 ) = zmxlm(:,:,1) ! surface set to the minimum value |
---|
234 | zmxld(:,:,jpk) = rn_lmin ! bottom set to the minimum value |
---|
235 | |
---|
236 | SELECT CASE ( nn_mxl ) |
---|
237 | ! |
---|
238 | CASE ( 0 ) ! bounded by the distance to surface and bottom |
---|
239 | DO jk = 2, jpkm1 |
---|
240 | DO jj = 2, jpjm1 |
---|
241 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
242 | zemxl = MIN( fsdepw(ji,jj,jk), zmxlm(ji,jj,jk), & |
---|
243 | & fsdepw(ji,jj,mbathy(ji,jj)) - fsdepw(ji,jj,jk) ) |
---|
244 | zmxlm(ji,jj,jk) = zemxl |
---|
245 | zmxld(ji,jj,jk) = zemxl |
---|
246 | END DO |
---|
247 | END DO |
---|
248 | END DO |
---|
249 | ! |
---|
250 | CASE ( 1 ) ! bounded by the vertical scale factor |
---|
251 | DO jk = 2, jpkm1 |
---|
252 | DO jj = 2, jpjm1 |
---|
253 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
254 | zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
255 | zmxlm(ji,jj,jk) = zemxl |
---|
256 | zmxld(ji,jj,jk) = zemxl |
---|
257 | END DO |
---|
258 | END DO |
---|
259 | END DO |
---|
260 | ! |
---|
261 | CASE ( 2 ) ! |dk[xml]| bounded by e3t : |
---|
262 | DO jk = 2, jpkm1 ! from the surface to the bottom : |
---|
263 | DO jj = 2, jpjm1 |
---|
264 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
265 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
266 | END DO |
---|
267 | END DO |
---|
268 | END DO |
---|
269 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : |
---|
270 | DO jj = 2, jpjm1 |
---|
271 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
272 | zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
273 | zmxlm(ji,jj,jk) = zemxl |
---|
274 | zmxld(ji,jj,jk) = zemxl |
---|
275 | END DO |
---|
276 | END DO |
---|
277 | END DO |
---|
278 | ! |
---|
279 | CASE ( 3 ) ! lup and ldown, |dk[xml]| bounded by e3t : |
---|
280 | DO jk = 2, jpkm1 ! from the surface to the bottom : lup |
---|
281 | DO jj = 2, jpjm1 |
---|
282 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
283 | zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
284 | END DO |
---|
285 | END DO |
---|
286 | END DO |
---|
287 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : ldown |
---|
288 | DO jj = 2, jpjm1 |
---|
289 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
290 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
291 | END DO |
---|
292 | END DO |
---|
293 | END DO |
---|
294 | !CDIR NOVERRCHK |
---|
295 | DO jk = 2, jpkm1 |
---|
296 | !CDIR NOVERRCHK |
---|
297 | DO jj = 2, jpjm1 |
---|
298 | !CDIR NOVERRCHK |
---|
299 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
300 | zemlm = MIN ( zmxld(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
301 | zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) ) |
---|
302 | zmxlm(ji,jj,jk) = zemlm |
---|
303 | zmxld(ji,jj,jk) = zemlp |
---|
304 | END DO |
---|
305 | END DO |
---|
306 | END DO |
---|
307 | ! |
---|
308 | END SELECT |
---|
309 | |
---|
310 | # if defined key_c1d |
---|
311 | ! c1d configuration : save mixing and dissipation turbulent length scales |
---|
312 | e_dis(:,:,:) = zmxld(:,:,:) |
---|
313 | e_mix(:,:,:) = zmxlm(:,:,:) |
---|
314 | # endif |
---|
315 | |
---|
316 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
317 | ! II TKE Langmuir circulation source term |
---|
318 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
319 | IF( ln_lc ) THEN |
---|
320 | ! |
---|
321 | ! Computation of total energy produce by LC : cumulative sum over jk |
---|
322 | zpelc(:,:,1) = MAX( rn2b(:,:,1), 0. ) * fsdepw(:,:,1) * fse3w(:,:,1) |
---|
323 | DO jk = 2, jpk |
---|
324 | zpelc(:,:,jk) = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0. ) * fsdepw(:,:,jk) * fse3w(:,:,jk) |
---|
325 | END DO |
---|
326 | ! |
---|
327 | ! Computation of finite Langmuir Circulation depth |
---|
328 | ! Initialization to the number of w ocean point mbathy |
---|
329 | imlc(:,:) = mbathy(:,:) |
---|
330 | DO jk = jpkm1, 2, -1 |
---|
331 | DO jj = 1, jpj |
---|
332 | DO ji = 1, jpi |
---|
333 | ! Last w-level at which zpelc>=0.5*us*us with us=0.016*wind(starting from jpk-1) |
---|
334 | zus = 0.000128 * wndm(ji,jj) * wndm(ji,jj) |
---|
335 | IF( zpelc(ji,jj,jk) > zus ) imlc(ji,jj) = jk |
---|
336 | END DO |
---|
337 | END DO |
---|
338 | END DO |
---|
339 | ! |
---|
340 | ! finite LC depth |
---|
341 | DO jj = 1, jpj |
---|
342 | DO ji = 1, jpi |
---|
343 | zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj)) |
---|
344 | END DO |
---|
345 | END DO |
---|
346 | ! |
---|
347 | # if defined key_c1d |
---|
348 | hlc(:,:) = zhlc(:,:) * tmask(:,:,1) ! c1d configuration: save finite Langmuir Circulation depth |
---|
349 | # endif |
---|
350 | ! |
---|
351 | ! TKE Langmuir circulation source term |
---|
352 | !CDIR NOVERRCHK |
---|
353 | DO jk = 2, jpkm1 |
---|
354 | !CDIR NOVERRCHK |
---|
355 | DO jj = 2, jpjm1 |
---|
356 | !CDIR NOVERRCHK |
---|
357 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
358 | ! Stokes drift |
---|
359 | zus = 0.016 * wndm(ji,jj) |
---|
360 | ! computation of vertical velocity due to LC |
---|
361 | zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) ) |
---|
362 | zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) ) |
---|
363 | ! TKE Langmuir circulation source term |
---|
364 | ztkelc(ji,jj,jk) = ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) |
---|
365 | END DO |
---|
366 | END DO |
---|
367 | END DO |
---|
368 | ! |
---|
369 | ENDIF |
---|
370 | |
---|
371 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
372 | ! III Tubulent kinetic energy time stepping |
---|
373 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
374 | |
---|
375 | ! 1. Vertical eddy viscosity on tke (put in zmxlm) and first estimate of avt |
---|
376 | ! --------------------------------------------------------------------- |
---|
377 | !CDIR NOVERRCHK |
---|
378 | DO jk = 2, jpkm1 |
---|
379 | !CDIR NOVERRCHK |
---|
380 | DO jj = 2, jpjm1 |
---|
381 | !CDIR NOVERRCHK |
---|
382 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
383 | zsqen = SQRT( en(ji,jj,jk) ) |
---|
384 | zav = rn_ediff * zmxlm(ji,jj,jk) * zsqen |
---|
385 | avt (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk) |
---|
386 | zmxlm(ji,jj,jk) = MAX( zav, avmb(jk) ) * tmask(ji,jj,jk) |
---|
387 | zmxld(ji,jj,jk) = zsqen / zmxld(ji,jj,jk) |
---|
388 | END DO |
---|
389 | END DO |
---|
390 | END DO |
---|
391 | |
---|
392 | ! 2. Surface boundary condition on tke and its eddy viscosity (zmxlm) |
---|
393 | ! ------------------------------------------------- |
---|
394 | ! en(1) = rn_ebb sqrt(utau^2+vtau^2) / rau0 (min value rn_emin0) |
---|
395 | ! zmxlm(1) = avmb(1) and zmxlm(jpk) = 0. |
---|
396 | !CDIR NOVERRCHK |
---|
397 | DO jj = 2, jpjm1 |
---|
398 | !CDIR NOVERRCHK |
---|
399 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
400 | ztx2 = utau(ji-1,jj ) + utau(ji,jj) |
---|
401 | zty2 = vtau(ji ,jj-1) + vtau(ji,jj) |
---|
402 | zesurf = zbbrau * SQRT( ztx2 * ztx2 + zty2 * zty2 ) |
---|
403 | en (ji,jj,1) = MAX( zesurf, rn_emin0 ) * tmask(ji,jj,1) |
---|
404 | zav = rn_ediff * zmxlm(ji,jj,1) * SQRT( en(ji,jj,1) ) |
---|
405 | zmxlm(ji,jj,1 ) = MAX( zav, avmb(1) ) * tmask(ji,jj,1) |
---|
406 | avt (ji,jj,1 ) = MAX( zav, avtb(1) * avtb_2d(ji,jj) ) * tmask(ji,jj,1) |
---|
407 | zmxlm(ji,jj,jpk) = 0.e0 |
---|
408 | END DO |
---|
409 | END DO |
---|
410 | |
---|
411 | ! 3. Now Turbulent kinetic energy (output in en) |
---|
412 | ! ------------------------------- |
---|
413 | ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
414 | ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ). |
---|
415 | |
---|
416 | SELECT CASE ( nn_pdl ) |
---|
417 | ! |
---|
418 | CASE ( 0 ) ! No Prandtl number |
---|
419 | DO jk = 2, jpkm1 |
---|
420 | DO jj = 2, jpjm1 |
---|
421 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
422 | ! ! shear prod. - stratif. destruction |
---|
423 | zcoef = 0.5 / fse3w(ji,jj,jk) |
---|
424 | zdku = zcoef * ( ub(ji-1, jj ,jk-1) + ub(ji,jj,jk-1) & ! shear |
---|
425 | & - ub(ji-1, jj ,jk ) - ub(ji,jj,jk ) ) |
---|
426 | zdkv = zcoef * ( vb( ji ,jj-1,jk-1) + vb(ji,jj,jk-1) & |
---|
427 | & - vb( ji ,jj-1,jk ) - vb(ji,jj,jk ) ) |
---|
428 | zesh2 = eboost * ( zdku*zdku + zdkv*zdkv ) - rn2b(ji,jj,jk) ! coefficient (zesh2) |
---|
429 | ! |
---|
430 | ! ! Matrix |
---|
431 | zcof = zfact1 * tmask(ji,jj,jk) |
---|
432 | ! ! lower diagonal |
---|
433 | avmv(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk ) + zmxlm(ji,jj,jk-1) ) & |
---|
434 | & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) |
---|
435 | ! ! upper diagonal |
---|
436 | avmu(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk+1) + zmxlm(ji,jj,jk ) ) & |
---|
437 | & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) ) |
---|
438 | ! ! diagonal |
---|
439 | zwd(ji,jj,jk) = 1. - avmv(ji,jj,jk) - avmu(ji,jj,jk) + zfact2 * zmxld(ji,jj,jk) |
---|
440 | ! |
---|
441 | ! ! right hand side in en |
---|
442 | IF( .NOT. ln_lc ) THEN ! No Langmuir cells |
---|
443 | en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en(ji,jj,jk) & |
---|
444 | & + rdt * zmxlm (ji,jj,jk) * zesh2 |
---|
445 | ELSE ! Langmuir cell source term |
---|
446 | en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en(ji,jj,jk) & |
---|
447 | & + rdt * zmxlm (ji,jj,jk) * zesh2 & |
---|
448 | & + rdt * ztkelc(ji,jj,jk) |
---|
449 | ENDIF |
---|
450 | END DO |
---|
451 | END DO |
---|
452 | END DO |
---|
453 | ! |
---|
454 | CASE ( 1 ) ! Prandtl number |
---|
455 | DO jk = 2, jpkm1 |
---|
456 | DO jj = 2, jpjm1 |
---|
457 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
458 | ! ! shear prod. - stratif. destruction |
---|
459 | zcoef = 0.5 / fse3w(ji,jj,jk) |
---|
460 | zdku = zcoef * ( ub(ji-1,jj ,jk-1) + ub(ji,jj,jk-1) & ! shear |
---|
461 | & - ub(ji-1,jj ,jk ) - ub(ji,jj,jk ) ) |
---|
462 | zdkv = zcoef * ( vb(ji ,jj-1,jk-1) + vb(ji,jj,jk-1) & |
---|
463 | & - vb(ji ,jj-1,jk ) - vb(ji,jj,jk ) ) |
---|
464 | zsh2 = zdku * zdku + zdkv * zdkv ! square of shear |
---|
465 | zri = MAX( rn2b(ji,jj,jk), 0. ) / ( zsh2 + 1.e-20 ) ! local Richardson number |
---|
466 | # if defined key_c1d |
---|
467 | e_ric(ji,jj,jk) = zri * tmask(ji,jj,jk) ! c1d config. : save Ri |
---|
468 | # endif |
---|
469 | zpdl = 1.0 ! Prandtl number |
---|
470 | IF( zri >= 0.2 ) zpdl = 0.2 / zri |
---|
471 | zpdl = MAX( 0.1, zpdl ) |
---|
472 | zesh2 = eboost * zsh2 - zpdl * rn2b(ji,jj,jk) ! coefficient (esh2) |
---|
473 | ! |
---|
474 | ! ! Matrix |
---|
475 | zcof = zfact1 * tmask(ji,jj,jk) |
---|
476 | ! ! lower diagonal |
---|
477 | avmv(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk ) + zmxlm(ji,jj,jk-1) ) & |
---|
478 | & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) |
---|
479 | ! ! upper diagonal |
---|
480 | avmu(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk+1) + zmxlm(ji,jj,jk ) ) & |
---|
481 | & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) ) |
---|
482 | ! ! diagonal |
---|
483 | zwd(ji,jj,jk) = 1. - avmv(ji,jj,jk) - avmu(ji,jj,jk) + zfact2 * zmxld(ji,jj,jk) |
---|
484 | ! |
---|
485 | ! ! right hand side in en |
---|
486 | IF( .NOT. ln_lc ) THEN ! No Langmuir cells |
---|
487 | en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en (ji,jj,jk) & |
---|
488 | & + rdt * zmxlm (ji,jj,jk) * zesh2 |
---|
489 | ELSE ! Langmuir cell source term |
---|
490 | en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en (ji,jj,jk) & |
---|
491 | & + rdt * zmxlm (ji,jj,jk) * zesh2 & |
---|
492 | & + rdt * ztkelc(ji,jj,jk) |
---|
493 | ENDIF |
---|
494 | zmxld(ji,jj,jk) = zpdl * tmask(ji,jj,jk) ! store masked Prandlt number in zmxld array |
---|
495 | END DO |
---|
496 | END DO |
---|
497 | END DO |
---|
498 | ! |
---|
499 | END SELECT |
---|
500 | |
---|
501 | # if defined key_c1d |
---|
502 | e_pdl(:,:,2:jpkm1) = zmxld(:,:,2:jpkm1) ! c1d configuration : save masked Prandlt number |
---|
503 | e_pdl(:,:, 1) = e_pdl(:,:, 2) |
---|
504 | e_pdl(:,:, jpk) = e_pdl(:,:, jpkm1) |
---|
505 | # endif |
---|
506 | |
---|
507 | ! 4. Matrix inversion from level 2 (tke prescribed at level 1) |
---|
508 | !!-------------------------------- |
---|
509 | DO jk = 3, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
510 | DO jj = 2, jpjm1 |
---|
511 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
512 | zwd(ji,jj,jk) = zwd(ji,jj,jk) - avmv(ji,jj,jk) * avmu(ji,jj,jk-1) / zwd(ji,jj,jk-1) |
---|
513 | END DO |
---|
514 | END DO |
---|
515 | END DO |
---|
516 | DO jj = 2, jpjm1 ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
517 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
518 | avmv(ji,jj,2) = en(ji,jj,2) - avmv(ji,jj,2) * en(ji,jj,1) ! Surface boudary conditions on tke |
---|
519 | END DO |
---|
520 | END DO |
---|
521 | DO jk = 3, jpkm1 |
---|
522 | DO jj = 2, jpjm1 |
---|
523 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
524 | avmv(ji,jj,jk) = en(ji,jj,jk) - avmv(ji,jj,jk) / zwd(ji,jj,jk-1) *avmv(ji,jj,jk-1) |
---|
525 | END DO |
---|
526 | END DO |
---|
527 | END DO |
---|
528 | DO jj = 2, jpjm1 ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
529 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
530 | en(ji,jj,jpkm1) = avmv(ji,jj,jpkm1) / zwd(ji,jj,jpkm1) |
---|
531 | END DO |
---|
532 | END DO |
---|
533 | DO jk = jpk-2, 2, -1 |
---|
534 | DO jj = 2, jpjm1 |
---|
535 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
536 | en(ji,jj,jk) = ( avmv(ji,jj,jk) - avmu(ji,jj,jk) * en(ji,jj,jk+1) ) / zwd(ji,jj,jk) |
---|
537 | END DO |
---|
538 | END DO |
---|
539 | END DO |
---|
540 | DO jk = 2, jpkm1 ! set the minimum value of tke |
---|
541 | DO jj = 2, jpjm1 |
---|
542 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
543 | en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * tmask(ji,jj,jk) |
---|
544 | END DO |
---|
545 | END DO |
---|
546 | END DO |
---|
547 | |
---|
548 | ! 5. Add extra TKE due to surface and internal wave breaking (nn_etau /= 0) |
---|
549 | !!---------------------------------------------------------- |
---|
550 | IF( nn_etau /= 0 ) THEN ! extra tke : en = en + rn_efr * en(1) * exp( -z/zhtau ) |
---|
551 | ! |
---|
552 | SELECT CASE( nn_htau ) ! Choice of the depth of penetration |
---|
553 | CASE( 0 ) ! constant depth penetration (here 10 meters) |
---|
554 | DO jj = 2, jpjm1 |
---|
555 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
556 | zhtau(ji,jj) = 10. |
---|
557 | END DO |
---|
558 | END DO |
---|
559 | CASE( 1 ) ! meridional profile 1 |
---|
560 | DO jj = 2, jpjm1 ! ( 0.5m in the tropics to a maximum of 30 m at high lat.) |
---|
561 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
562 | zhtau(ji,jj) = MAX( 0.5, 3./4. * MIN( 40., 60.*ABS( SIN( rpi/180. * gphit(ji,jj) ) ) ) ) |
---|
563 | END DO |
---|
564 | END DO |
---|
565 | END SELECT |
---|
566 | ! |
---|
567 | IF( nn_etau == 1 ) THEN ! extra term throughout the water column |
---|
568 | DO jk = 2, jpkm1 |
---|
569 | DO jj = 2, jpjm1 |
---|
570 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
571 | en(ji,jj,jk) = en(ji,jj,jk) & |
---|
572 | & + rn_efr * en(ji,jj,1)*EXP( -fsdepw(ji,jj,jk) / zhtau(ji,jj) ) & |
---|
573 | & * ( 1.e0 - fr_i(ji,jj) ) * tmask(ji,jj,jk) |
---|
574 | END DO |
---|
575 | END DO |
---|
576 | END DO |
---|
577 | ELSEIF( nn_etau == 2 ) THEN ! extra term only at the base of the mixed layer |
---|
578 | DO jj = 2, jpjm1 |
---|
579 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
580 | jk = nmln(ji,jj) |
---|
581 | en(ji,jj,jk) = en(ji,jj,jk) & |
---|
582 | & + rn_efr * en(ji,jj,1)*EXP( -fsdepw(ji,jj,jk) / zhtau(ji,jj) ) & |
---|
583 | & * ( 1.e0 - fr_i(ji,jj) ) * tmask(ji,jj,jk) |
---|
584 | END DO |
---|
585 | END DO |
---|
586 | ENDIF |
---|
587 | ! |
---|
588 | ENDIF |
---|
589 | |
---|
590 | |
---|
591 | ! Lateral boundary conditions (sign unchanged) |
---|
592 | CALL lbc_lnk( en, 'W', 1. ) ; CALL lbc_lnk( avt, 'W', 1. ) |
---|
593 | |
---|
594 | |
---|
595 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
596 | ! IV. Before vertical eddy vicosity and diffusivity coefficients |
---|
597 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
598 | ! |
---|
599 | SELECT CASE ( nn_ave ) |
---|
600 | CASE ( 0 ) ! no horizontal average |
---|
601 | DO jk = 2, jpkm1 ! only vertical eddy viscosity |
---|
602 | DO jj = 2, jpjm1 |
---|
603 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
604 | avmu(ji,jj,jk) = ( avt (ji,jj,jk) + avt (ji+1,jj ,jk) ) * umask(ji,jj,jk) & |
---|
605 | & / MAX( 1., tmask(ji,jj,jk) + tmask(ji+1,jj ,jk) ) |
---|
606 | avmv(ji,jj,jk) = ( avt (ji,jj,jk) + avt (ji ,jj+1,jk) ) * vmask(ji,jj,jk) & |
---|
607 | & / MAX( 1., tmask(ji,jj,jk) + tmask(ji ,jj+1,jk) ) |
---|
608 | END DO |
---|
609 | END DO |
---|
610 | END DO |
---|
611 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions |
---|
612 | ! |
---|
613 | CASE ( 1 ) ! horizontal average ( 1/2 1/2 ) |
---|
614 | ! ! Vertical eddy viscosity avmu = 1/4 ( 1 1 ) |
---|
615 | ! ! ( 1/2 1/2 ) |
---|
616 | ! ! |
---|
617 | ! ! ( 1/2 1 1/2 ) |
---|
618 | ! ! avmv = 1/4 ( 1/2 1 1/2 ) |
---|
619 | DO jk = 2, jpkm1 |
---|
620 | DO jj = 2, jpjm1 |
---|
621 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
622 | # if defined key_vectopt_memory |
---|
623 | ! ! caution : vectopt_memory change the solution |
---|
624 | ! ! (last digit of the solver stat) |
---|
625 | avmu(ji,jj,jk) = ( avt(ji,jj ,jk) + avt(ji+1,jj ,jk) & |
---|
626 | & +.5*( avt(ji,jj-1,jk) + avt(ji+1,jj-1,jk) & |
---|
627 | & +avt(ji,jj+1,jk) + avt(ji+1,jj+1,jk) ) ) * eumean(ji,jj,jk) |
---|
628 | |
---|
629 | avmv(ji,jj,jk) = ( avt(ji ,jj,jk) + avt(ji ,jj+1,jk) & |
---|
630 | & +.5*( avt(ji-1,jj,jk) + avt(ji-1,jj+1,jk) & |
---|
631 | & +avt(ji+1,jj,jk) + avt(ji+1,jj+1,jk) ) ) * evmean(ji,jj,jk) |
---|
632 | # else |
---|
633 | avmu(ji,jj,jk) = ( avt (ji,jj ,jk) + avt (ji+1,jj ,jk) & |
---|
634 | & +.5*( avt (ji,jj-1,jk) + avt (ji+1,jj-1,jk) & |
---|
635 | & +avt (ji,jj+1,jk) + avt (ji+1,jj+1,jk) ) ) * umask(ji,jj,jk) & |
---|
636 | & / MAX( 1., tmask(ji,jj ,jk) + tmask(ji+1,jj ,jk) & |
---|
637 | & +.5*( tmask(ji,jj-1,jk) + tmask(ji+1,jj-1,jk) & |
---|
638 | & +tmask(ji,jj+1,jk) + tmask(ji+1,jj+1,jk) ) ) |
---|
639 | |
---|
640 | avmv(ji,jj,jk) = ( avt (ji ,jj,jk) + avt (ji ,jj+1,jk) & |
---|
641 | & +.5*( avt (ji-1,jj,jk) + avt (ji-1,jj+1,jk) & |
---|
642 | & +avt (ji+1,jj,jk) + avt (ji+1,jj+1,jk) ) ) * vmask(ji,jj,jk) & |
---|
643 | & / MAX( 1., tmask(ji ,jj,jk) + tmask(ji ,jj+1,jk) & |
---|
644 | & +.5*( tmask(ji-1,jj,jk) + tmask(ji-1,jj+1,jk) & |
---|
645 | & +tmask(ji+1,jj,jk) + tmask(ji+1,jj+1,jk) ) ) |
---|
646 | # endif |
---|
647 | END DO |
---|
648 | END DO |
---|
649 | END DO |
---|
650 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions |
---|
651 | ! |
---|
652 | ! ! Vertical eddy diffusivity (1 2 1) |
---|
653 | ! ! avt = 1/16 (2 4 2) |
---|
654 | ! ! (1 2 1) |
---|
655 | DO jk = 2, jpkm1 |
---|
656 | DO jj = 2, jpjm1 |
---|
657 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
658 | # if defined key_vectopt_memory |
---|
659 | avt(ji,jj,jk) = ( avmu(ji,jj,jk) + avmu(ji-1,jj ,jk) & |
---|
660 | & + avmv(ji,jj,jk) + avmv(ji ,jj-1,jk) ) * etmean(ji,jj,jk) |
---|
661 | # else |
---|
662 | avt(ji,jj,jk) = ( avmu (ji,jj,jk) + avmu (ji-1,jj ,jk) & |
---|
663 | & + avmv (ji,jj,jk) + avmv (ji ,jj-1,jk) ) * tmask(ji,jj,jk) & |
---|
664 | & / MAX( 1., umask(ji,jj,jk) + umask(ji-1,jj ,jk) & |
---|
665 | & + vmask(ji,jj,jk) + vmask(ji ,jj-1,jk) ) |
---|
666 | # endif |
---|
667 | END DO |
---|
668 | END DO |
---|
669 | END DO |
---|
670 | ! |
---|
671 | END SELECT |
---|
672 | ! |
---|
673 | IF( nn_pdl == 1 ) THEN ! Ponderation by the Prandtl number (nn_pdl=1) |
---|
674 | DO jk = 2, jpkm1 |
---|
675 | DO jj = 2, jpjm1 |
---|
676 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
677 | zpdl = zmxld(ji,jj,jk) |
---|
678 | avt(ji,jj,jk) = MAX( zpdl * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk) |
---|
679 | END DO |
---|
680 | END DO |
---|
681 | END DO |
---|
682 | ENDIF |
---|
683 | ! |
---|
684 | DO jk = 2, jpkm1 ! Minimum value on the eddy viscosity |
---|
685 | avmu(:,:,jk) = MAX( avmu(:,:,jk), avmb(jk) ) * umask(:,:,jk) |
---|
686 | avmv(:,:,jk) = MAX( avmv(:,:,jk), avmb(jk) ) * vmask(:,:,jk) |
---|
687 | END DO |
---|
688 | |
---|
689 | CALL lbc_lnk( avt, 'W', 1. ) ! Lateral boundary conditions on avt (sign unchanged) |
---|
690 | |
---|
691 | IF( lrst_oce ) CALL tke_rst( kt, 'WRITE' ) ! write en in restart file |
---|
692 | |
---|
693 | IF(ln_ctl) THEN |
---|
694 | CALL prt_ctl( tab3d_1=en , clinfo1=' tke - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk) |
---|
695 | CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke - u: ', mask1=umask, & |
---|
696 | & tab3d_2=avmv, clinfo2= ' v: ', mask2=vmask, ovlap=1, kdim=jpk ) |
---|
697 | ENDIF |
---|
698 | ! |
---|
699 | END SUBROUTINE zdf_tke_old |
---|
700 | |
---|
701 | |
---|
702 | SUBROUTINE zdf_tke_init |
---|
703 | !!---------------------------------------------------------------------- |
---|
704 | !! *** ROUTINE zdf_tke_init *** |
---|
705 | !! |
---|
706 | !! ** Purpose : Initialization of the vertical eddy diffivity and |
---|
707 | !! viscosity when using a tke turbulent closure scheme |
---|
708 | !! |
---|
709 | !! ** Method : Read the namzdf_tke namelist and check the parameters |
---|
710 | !! called at the first timestep (nit000) |
---|
711 | !! |
---|
712 | !! ** input : Namlist namzdf_tke |
---|
713 | !! |
---|
714 | !! ** Action : Increase by 1 the nstop flag is setting problem encounter |
---|
715 | !! |
---|
716 | !!---------------------------------------------------------------------- |
---|
717 | USE dynzdf_exp |
---|
718 | USE trazdf_exp |
---|
719 | ! |
---|
720 | # if defined key_vectopt_memory |
---|
721 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
722 | # else |
---|
723 | INTEGER :: jk ! dummy loop indices |
---|
724 | # endif |
---|
725 | !! |
---|
726 | NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin , & |
---|
727 | & rn_emin0, rn_bshear, nn_mxl , ln_mxl0 , & |
---|
728 | & rn_lmin , rn_lmin0 , nn_pdl , nn_etau , & |
---|
729 | & nn_htau , rn_efr , rn_addhft, rn_sclhft, & |
---|
730 | & ln_lc , rn_lc |
---|
731 | !!---------------------------------------------------------------------- |
---|
732 | |
---|
733 | ! Read Namelist namzdf_tke : Turbulente Kinetic Energy |
---|
734 | ! -------------------- |
---|
735 | REWIND ( numnam ) |
---|
736 | READ ( numnam, namzdf_tke ) |
---|
737 | |
---|
738 | ! Compute boost associated with the Richardson critic |
---|
739 | ! (control values: rn_cri = 0.3 ==> eboost=1.25 for nn_pdl=1) |
---|
740 | ! ( rn_cri = 0.222 ==> eboost=1. ) |
---|
741 | eboost = rn_cri * ( 2. + rn_ediss / rn_ediff ) / 2. |
---|
742 | |
---|
743 | |
---|
744 | |
---|
745 | ! Parameter control and print |
---|
746 | ! --------------------------- |
---|
747 | ! Control print |
---|
748 | IF(lwp) THEN |
---|
749 | WRITE(numout,*) |
---|
750 | WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme (old scheme)' |
---|
751 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
752 | WRITE(numout,*) ' Namelist namzdf_tke : set tke mixing parameters' |
---|
753 | WRITE(numout,*) ' coef. to compute avt rn_ediff = ', rn_ediff |
---|
754 | WRITE(numout,*) ' Kolmogoroff dissipation coef. rn_ediss = ', rn_ediss |
---|
755 | WRITE(numout,*) ' tke surface input coef. rn_ebb = ', rn_ebb |
---|
756 | WRITE(numout,*) ' minimum value of tke rn_emin = ', rn_emin |
---|
757 | WRITE(numout,*) ' surface minimum value of tke rn_emin0 = ', rn_emin0 |
---|
758 | WRITE(numout,*) ' mixing length type nn_mxl = ', nn_mxl |
---|
759 | WRITE(numout,*) ' prandl number flag nn_pdl = ', nn_pdl |
---|
760 | WRITE(numout,*) ' horizontal average flag nn_ave = ', nn_ave |
---|
761 | WRITE(numout,*) ' critic Richardson nb rn_cri = ', rn_cri |
---|
762 | WRITE(numout,*) ' and its associated coeff. eboost = ', eboost |
---|
763 | WRITE(numout,*) ' surface mixing length = F(stress) or not ln_mxl0 = ', ln_mxl0 |
---|
764 | WRITE(numout,*) ' surface mixing length minimum value rn_lmin0 = ', rn_lmin0 |
---|
765 | WRITE(numout,*) ' interior mixing length minimum value rn_lmin0 = ', rn_lmin0 |
---|
766 | WRITE(numout,*) ' test param. to add tke induced by wind nn_etau = ', nn_etau |
---|
767 | WRITE(numout,*) ' flag for computation of exp. tke profile nn_htau = ', nn_htau |
---|
768 | WRITE(numout,*) ' % of rn_emin0 which pene. the thermocline rn_efr = ', rn_efr |
---|
769 | WRITE(numout,*) ' flag to take into acc. Langmuir circ. ln_lc = ', ln_lc |
---|
770 | WRITE(numout,*) ' coef to compute verticla velocity of LC rn_lc = ', rn_lc |
---|
771 | WRITE(numout,*) |
---|
772 | ENDIF |
---|
773 | |
---|
774 | ! Check of some namelist values |
---|
775 | IF( nn_mxl < 0 .OR. nn_mxl > 3 ) CALL ctl_stop( 'bad flag: nn_mxl is 0, 1 or 2 ' ) |
---|
776 | IF( nn_pdl < 0 .OR. nn_pdl > 1 ) CALL ctl_stop( 'bad flag: nn_pdl is 0 or 1 ' ) |
---|
777 | IF( nn_ave < 0 .OR. nn_ave > 1 ) CALL ctl_stop( 'bad flag: nn_ave is 0 or 1 ' ) |
---|
778 | IF( nn_htau < 0 .OR. nn_htau > 1 ) CALL ctl_stop( 'bad flag: nn_htau is 0 or 1 ' ) |
---|
779 | IF( rn_lc < 0.15 .OR. rn_lc > 0.2 ) CALL ctl_stop( 'bad value: rn_lc must be between 0.15 and 0.2 ' ) |
---|
780 | |
---|
781 | IF( nn_etau == 2 ) CALL zdf_mxl( nit000 ) ! Initialization of nmln |
---|
782 | |
---|
783 | ! Horizontal average : initialization of weighting arrays |
---|
784 | ! ------------------- |
---|
785 | ! |
---|
786 | SELECT CASE ( nn_ave ) |
---|
787 | ! Notice: mean arrays have not to by defined at local domain boundaries due to the vertical nature of TKE |
---|
788 | ! |
---|
789 | CASE ( 0 ) ! no horizontal average |
---|
790 | IF(lwp) WRITE(numout,*) ' no horizontal average on avt, avmu, avmv' |
---|
791 | IF(lwp) WRITE(numout,*) ' only in very high horizontal resolution !' |
---|
792 | # if defined key_vectopt_memory |
---|
793 | ! caution vectopt_memory change the solution (last digit of the solver stat) |
---|
794 | ! weighting mean arrays etmean, eumean and evmean |
---|
795 | ! ( 1 1 ) ( 1 ) |
---|
796 | ! avt = 1/4 ( 1 1 ) avmu = 1/2 ( 1 1 ) avmv= 1/2 ( 1 ) |
---|
797 | ! |
---|
798 | etmean(:,:,:) = 0.e0 |
---|
799 | eumean(:,:,:) = 0.e0 |
---|
800 | evmean(:,:,:) = 0.e0 |
---|
801 | ! |
---|
802 | DO jk = 1, jpkm1 |
---|
803 | DO jj = 2, jpjm1 |
---|
804 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
805 | etmean(ji,jj,jk) = tmask(ji,jj,jk) / MAX( 1., umask(ji,jj,jk) + umask(ji-1,jj ,jk) & |
---|
806 | & + vmask(ji,jj,jk) + vmask(ji ,jj-1,jk) ) |
---|
807 | eumean(ji,jj,jk) = umask(ji,jj,jk) / MAX( 1., tmask(ji,jj,jk) + tmask(ji+1,jj ,jk) ) |
---|
808 | evmean(ji,jj,jk) = vmask(ji,jj,jk) / MAX( 1., tmask(ji,jj,jk) + tmask(ji ,jj+1,jk) ) |
---|
809 | END DO |
---|
810 | END DO |
---|
811 | END DO |
---|
812 | # endif |
---|
813 | ! |
---|
814 | CASE ( 1 ) ! horizontal average |
---|
815 | IF(lwp) WRITE(numout,*) ' horizontal average on avt, avmu, avmv' |
---|
816 | # if defined key_vectopt_memory |
---|
817 | ! caution vectopt_memory change the solution (last digit of the solver stat) |
---|
818 | ! weighting mean arrays etmean, eumean and evmean |
---|
819 | ! ( 1 1 ) ( 1/2 1/2 ) ( 1/2 1 1/2 ) |
---|
820 | ! avt = 1/4 ( 1 1 ) avmu = 1/4 ( 1 1 ) avmv= 1/4 ( 1/2 1 1/2 ) |
---|
821 | ! ( 1/2 1/2 ) |
---|
822 | etmean(:,:,:) = 0.e0 |
---|
823 | eumean(:,:,:) = 0.e0 |
---|
824 | evmean(:,:,:) = 0.e0 |
---|
825 | ! |
---|
826 | DO jk = 1, jpkm1 |
---|
827 | DO jj = 2, jpjm1 |
---|
828 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
829 | etmean(ji,jj,jk) = tmask(ji,jj,jk) / MAX( 1., umask(ji-1,jj ,jk) + umask(ji ,jj ,jk) & |
---|
830 | & + vmask(ji ,jj-1,jk) + vmask(ji ,jj ,jk) ) |
---|
831 | eumean(ji,jj,jk) = umask(ji,jj,jk) / MAX( 1., tmask(ji ,jj ,jk) + tmask(ji+1,jj ,jk) & |
---|
832 | & +.5 * ( tmask(ji ,jj-1,jk) + tmask(ji+1,jj-1,jk) & |
---|
833 | & + tmask(ji ,jj+1,jk) + tmask(ji+1,jj+1,jk) ) ) |
---|
834 | evmean(ji,jj,jk) = vmask(ji,jj,jk) / MAX( 1., tmask(ji ,jj ,jk) + tmask(ji ,jj+1,jk) & |
---|
835 | & +.5 * ( tmask(ji-1,jj ,jk) + tmask(ji-1,jj+1,jk) & |
---|
836 | & + tmask(ji+1,jj ,jk) + tmask(ji+1,jj+1,jk) ) ) |
---|
837 | END DO |
---|
838 | END DO |
---|
839 | END DO |
---|
840 | # endif |
---|
841 | ! |
---|
842 | END SELECT |
---|
843 | |
---|
844 | ! Initialization of vertical eddy coef. to the background value |
---|
845 | ! ------------------------------------------------------------- |
---|
846 | DO jk = 1, jpk |
---|
847 | avt (:,:,jk) = avtb(jk) * tmask(:,:,jk) |
---|
848 | avmu(:,:,jk) = avmb(jk) * umask(:,:,jk) |
---|
849 | avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk) |
---|
850 | END DO |
---|
851 | |
---|
852 | |
---|
853 | ! read or initialize turbulent kinetic energy ( en ) |
---|
854 | ! ------------------------------------------------- |
---|
855 | CALL tke_rst( nit000, 'READ' ) |
---|
856 | ! |
---|
857 | END SUBROUTINE zdf_tke_init |
---|
858 | |
---|
859 | |
---|
860 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
861 | !!--------------------------------------------------------------------- |
---|
862 | !! *** ROUTINE ts_rst *** |
---|
863 | !! |
---|
864 | !! ** Purpose : Read or write TKE file (en) in restart file |
---|
865 | !! |
---|
866 | !! ** Method : use of IOM library |
---|
867 | !! if the restart does not contain TKE, en is either |
---|
868 | !! set to rn_emin or recomputed (nn_itke/=0) |
---|
869 | !!---------------------------------------------------------------------- |
---|
870 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
871 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
872 | ! |
---|
873 | INTEGER :: jit ! dummy loop indices |
---|
874 | !!---------------------------------------------------------------------- |
---|
875 | ! |
---|
876 | IF( TRIM(cdrw) == 'READ' ) THEN |
---|
877 | IF( ln_rstart ) THEN |
---|
878 | IF( iom_varid( numror, 'en', ldstop = .FALSE. ) > 0 ) THEN |
---|
879 | CALL iom_get( numror, jpdom_autoglo, 'en', en ) |
---|
880 | ELSE |
---|
881 | IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme' |
---|
882 | IF(lwp) WRITE(numout,*) ' ===>>>> : en is set by iterative loop' |
---|
883 | en (:,:,:) = rn_emin * tmask(:,:,:) |
---|
884 | DO jit = 2, 51 |
---|
885 | CALL zdf_tke_old( jit ) |
---|
886 | END DO |
---|
887 | ENDIF |
---|
888 | ELSE |
---|
889 | en(:,:,:) = rn_emin * tmask(:,:,:) ! no restart: en set to its minimum value |
---|
890 | ENDIF |
---|
891 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN |
---|
892 | CALL iom_rstput( kt, nitrst, numrow, 'en', en ) |
---|
893 | ENDIF |
---|
894 | ! |
---|
895 | END SUBROUTINE tke_rst |
---|
896 | |
---|
897 | #else |
---|
898 | !!---------------------------------------------------------------------- |
---|
899 | !! Dummy module : NO TKE scheme |
---|
900 | !!---------------------------------------------------------------------- |
---|
901 | LOGICAL, PUBLIC, PARAMETER :: lk_zdftke_old = .FALSE. !: TKE flag |
---|
902 | CONTAINS |
---|
903 | SUBROUTINE zdf_tke_old( kt ) ! Empty routine |
---|
904 | WRITE(*,*) 'zdf_tke_old: You should not have seen this print! error?', kt |
---|
905 | END SUBROUTINE zdf_tke_old |
---|
906 | #endif |
---|
907 | |
---|
908 | !!====================================================================== |
---|
909 | END MODULE zdftke_old |
---|