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 | !!---------------------------------------------------------------------- |
---|
28 | #if defined key_zdftke || defined key_esopa |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | !! 'key_zdftke' TKE vertical physics |
---|
31 | !!---------------------------------------------------------------------- |
---|
32 | !! zdf_tke : update momentum and tracer Kz from a tke scheme |
---|
33 | !! tke_tke : tke time stepping: update tke at now time step (en) |
---|
34 | !! tke_avn : compute mixing length scale and deduce avm and avt |
---|
35 | !! tke_init : initialization, namelist read, and parameters control |
---|
36 | !! tke_rst : read/write tke restart in ocean restart file |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | USE oce ! ocean dynamics and active tracers |
---|
39 | USE dom_oce ! ocean space and time domain |
---|
40 | USE domvvl ! ocean space and time domain : variable volume layer |
---|
41 | USE zdf_oce ! ocean vertical physics |
---|
42 | USE sbc_oce ! surface boundary condition: ocean |
---|
43 | USE phycst ! physical constants |
---|
44 | USE zdfmxl ! mixed layer |
---|
45 | USE restart ! only for lrst_oce |
---|
46 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
47 | USE prtctl ! Print control |
---|
48 | USE in_out_manager ! I/O manager |
---|
49 | USE iom ! I/O manager library |
---|
50 | USE zdfbfr ! bottom friction |
---|
51 | |
---|
52 | IMPLICIT NONE |
---|
53 | PRIVATE |
---|
54 | |
---|
55 | PUBLIC zdf_tke ! routine called in step module |
---|
56 | PUBLIC tke_rst ! routine called in step module |
---|
57 | |
---|
58 | LOGICAL , PUBLIC, PARAMETER :: lk_zdftke = .TRUE. !: TKE vertical mixing flag |
---|
59 | |
---|
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 | #endif |
---|
65 | |
---|
66 | ! !!! ** Namelist namzdf_tke ** |
---|
67 | LOGICAL :: ln_mxl0 = .FALSE. ! mixing length scale surface value as function of wind stress or not |
---|
68 | INTEGER :: nn_mxl = 2 ! type of mixing length (=0/1/2/3) |
---|
69 | REAL(wp) :: rn_lmin0 = 0.4_wp ! surface min value of mixing length [m] |
---|
70 | REAL(wp) :: rn_lmin = 0.1_wp ! interior min value of mixing length [m] |
---|
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 ! background shear (>0) |
---|
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_efr = 1.0_wp ! fraction of TKE surface value which penetrates in the ocean |
---|
81 | LOGICAL :: ln_lc = .FALSE. ! Langmuir cells (LC) as a source term of TKE or not |
---|
82 | REAL(wp) :: rn_lc = 0.15_wp ! coef to compute vertical velocity of Langmuir cells |
---|
83 | |
---|
84 | REAL(wp) :: ri_cri ! critic Richardson number (deduced from rn_ediff and rn_ediss values) |
---|
85 | |
---|
86 | REAL(wp), DIMENSION(jpi,jpj) :: htau ! depth of tke penetration (nn_htau) |
---|
87 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: en ! now turbulent kinetic energy [m2/s2] |
---|
88 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: dissl ! now mixing lenght of dissipation |
---|
89 | |
---|
90 | !! * Substitutions |
---|
91 | # include "domzgr_substitute.h90" |
---|
92 | # include "vectopt_loop_substitute.h90" |
---|
93 | !!---------------------------------------------------------------------- |
---|
94 | !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009) |
---|
95 | !! $Id: zdftke2.F90 1201 2008-09-24 13:24:21Z rblod $ |
---|
96 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
97 | !!---------------------------------------------------------------------- |
---|
98 | |
---|
99 | CONTAINS |
---|
100 | |
---|
101 | SUBROUTINE zdf_tke( kt ) |
---|
102 | !!---------------------------------------------------------------------- |
---|
103 | !! *** ROUTINE zdf_tke *** |
---|
104 | !! |
---|
105 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
106 | !! coefficients using a turbulent closure scheme (TKE). |
---|
107 | !! |
---|
108 | !! ** Method : The time evolution of the turbulent kinetic energy (tke) |
---|
109 | !! is computed from a prognostic equation : |
---|
110 | !! d(en)/dt = avm (d(u)/dz)**2 ! shear production |
---|
111 | !! + d( avm d(en)/dz )/dz ! diffusion of tke |
---|
112 | !! + avt N^2 ! stratif. destruc. |
---|
113 | !! - rn_ediss / emxl en**(2/3) ! Kolmogoroff dissipation |
---|
114 | !! with the boundary conditions: |
---|
115 | !! surface: en = max( rn_emin0, rn_ebb * taum ) |
---|
116 | !! bottom : en = rn_emin |
---|
117 | !! The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff) |
---|
118 | !! |
---|
119 | !! The now Turbulent kinetic energy is computed using the following |
---|
120 | !! time stepping: implicit for vertical diffusion term, linearized semi |
---|
121 | !! implicit for kolmogoroff dissipation term, and explicit forward for |
---|
122 | !! both buoyancy and shear production terms. Therefore a tridiagonal |
---|
123 | !! linear system is solved. Note that buoyancy and shear terms are |
---|
124 | !! discretized in a energy conserving form (Bruchard 2002). |
---|
125 | !! |
---|
126 | !! The dissipative and mixing length scale are computed from en and |
---|
127 | !! the stratification (see tke_avn) |
---|
128 | !! |
---|
129 | !! The now vertical eddy vicosity and diffusivity coefficients are |
---|
130 | !! given by: |
---|
131 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
132 | !! avt = max( avmb, pdl * avm ) |
---|
133 | !! eav = max( avmb, avm ) |
---|
134 | !! where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and |
---|
135 | !! given by an empirical funtion of the localRichardson number if nn_pdl=1 |
---|
136 | !! |
---|
137 | !! ** Action : compute en (now turbulent kinetic energy) |
---|
138 | !! update avt, avmu, avmv (before vertical eddy coef.) |
---|
139 | !! |
---|
140 | !! References : Gaspar et al., JGR, 1990, |
---|
141 | !! Blanke and Delecluse, JPO, 1991 |
---|
142 | !! Mellor and Blumberg, JPO 2004 |
---|
143 | !! Axell, JGR, 2002 |
---|
144 | !! Bruchard OM 2002 |
---|
145 | !!---------------------------------------------------------------------- |
---|
146 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
147 | !!---------------------------------------------------------------------- |
---|
148 | ! |
---|
149 | IF( kt == nit000 ) CALL tke_init ! initialisation |
---|
150 | ! |
---|
151 | CALL tke_tke ! now tke (en) |
---|
152 | ! |
---|
153 | CALL tke_avn ! now avt, avm, avmu, avmv |
---|
154 | ! |
---|
155 | END SUBROUTINE zdf_tke |
---|
156 | |
---|
157 | |
---|
158 | SUBROUTINE tke_tke |
---|
159 | !!---------------------------------------------------------------------- |
---|
160 | !! *** ROUTINE tke_tke *** |
---|
161 | !! |
---|
162 | !! ** Purpose : Compute the now Turbulente Kinetic Energy (TKE) |
---|
163 | !! |
---|
164 | !! ** Method : - TKE surface boundary condition |
---|
165 | !! - source term due to Langmuir cells (ln_lc=T) |
---|
166 | !! - source term due to shear (saved in avmu, avmv arrays) |
---|
167 | !! - Now TKE : resolution of the TKE equation by inverting |
---|
168 | !! a tridiagonal linear system by a "methode de chasse" |
---|
169 | !! - increase TKE due to surface and internal wave breaking |
---|
170 | !! |
---|
171 | !! ** Action : - en : now turbulent kinetic energy) |
---|
172 | !! - avmu, avmv : production of TKE by shear at u and v-points |
---|
173 | !! (= Kz dz[Ub] * dz[Un] ) |
---|
174 | !! --------------------------------------------------------------------- |
---|
175 | USE oce, zdiag => ua ! use ua as workspace |
---|
176 | USE oce, zd_up => va ! use va as workspace |
---|
177 | USE oce, zd_lw => ta ! use ta as workspace |
---|
178 | !! |
---|
179 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
180 | REAL(wp) :: zbbrau, zesh2 ! temporary scalars |
---|
181 | REAL(wp) :: zfact1, zfact2, zfact3 ! - - |
---|
182 | REAL(wp) :: ztx2 , zty2 , zcof ! - - |
---|
183 | REAL(wp) :: zus , zwlc , zind ! - - |
---|
184 | REAL(wp) :: zzd_up, zzd_lw ! - - |
---|
185 | INTEGER :: ikbu, ikbv, ikbum1, ikbvm1 ! temporary scalar |
---|
186 | INTEGER :: ikbt, ikbumm1, ikbvmm1 ! temporary scalar |
---|
187 | REAL(wp) :: zebot ! temporary scalars |
---|
188 | INTEGER , DIMENSION(jpi,jpj) :: imlc ! 2D workspace |
---|
189 | REAL(wp), DIMENSION(jpi,jpj) :: zhlc ! - - |
---|
190 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: zpelc ! 3D workspace |
---|
191 | !!-------------------------------------------------------------------- |
---|
192 | ! |
---|
193 | zbbrau = rn_ebb / rau0 ! Local constant initialisation |
---|
194 | zfact1 = -.5 * rdt |
---|
195 | zfact2 = 1.5 * rdt * rn_ediss |
---|
196 | zfact3 = 0.5 * rn_ediss |
---|
197 | ! |
---|
198 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
199 | ! ! Surface boundary condition on tke |
---|
200 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
201 | DO jj = 2, jpjm1 ! en(1) = rn_ebb taum / rau0 (min value rn_emin0) |
---|
202 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
203 | en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1) |
---|
204 | END DO |
---|
205 | END DO |
---|
206 | ! |
---|
207 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
208 | ! ! Bottom boundary condition on tke |
---|
209 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
210 | ! en(bot) = 0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin) |
---|
211 | !CDIR NOVERRCHK |
---|
212 | DO jj = 2, jpjm1 |
---|
213 | !CDIR NOVERRCHK |
---|
214 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
215 | ikbu = MIN( mbathy(ji+1,jj), mbathy(ji,jj) ) |
---|
216 | ikbv = MIN( mbathy(ji,jj+1), mbathy(ji,jj) ) |
---|
217 | ikbum1 = MAX( ikbu-1, 1 ) |
---|
218 | ikbvm1 = MAX( ikbv-1, 1 ) |
---|
219 | ikbu = MIN( mbathy(ji,jj), mbathy(ji-1,jj) ) |
---|
220 | ikbv = MIN( mbathy(ji,jj), mbathy(ji,jj-1) ) |
---|
221 | ikbumm1 = MAX( ikbu-1, 1 ) |
---|
222 | ikbvmm1 = MAX( ikbv-1, 1 ) |
---|
223 | ikbt = MAX( mbathy(ji,jj), 1 ) |
---|
224 | ztx2 = bfrua(ji-1,jj) * fse3u(ji-1,jj ,ikbumm1) + & |
---|
225 | bfrua(ji ,jj) * fse3u(ji ,jj ,ikbum1 ) |
---|
226 | zty2 = bfrva(ji,jj ) * fse3v(ji ,jj ,ikbvm1) + & |
---|
227 | bfrva(ji,jj-1) * fse3v(ji ,jj-1,ikbvmm1 ) |
---|
228 | zebot = 0.25_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) |
---|
229 | en (ji,jj,ikbt) = MAX( zebot, rn_emin ) * tmask(ji,jj,1) |
---|
230 | END DO |
---|
231 | END DO |
---|
232 | ! |
---|
233 | ! |
---|
234 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
235 | IF( ln_lc ) THEN ! Langmuir circulation source term added to tke |
---|
236 | ! !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
237 | ! |
---|
238 | ! !* total energy produce by LC : cumulative sum over jk |
---|
239 | zpelc(:,:,1) = MAX( rn2b(:,:,1), 0. ) * fsdepw(:,:,1) * fse3w(:,:,1) |
---|
240 | DO jk = 2, jpk |
---|
241 | zpelc(:,:,jk) = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0. ) * fsdepw(:,:,jk) * fse3w(:,:,jk) |
---|
242 | END DO |
---|
243 | ! !* finite Langmuir Circulation depth |
---|
244 | imlc(:,:) = mbathy(:,:) ! Initialization to the number of w ocean point mbathy |
---|
245 | DO jk = jpkm1, 2, -1 |
---|
246 | DO jj = 1, jpj ! Last w-level at which zpelc>=0.5*us*us |
---|
247 | DO ji = 1, jpi ! with us=0.016*wind(starting from jpk-1) |
---|
248 | zus = 0.000128 * wndm(ji,jj) * wndm(ji,jj) |
---|
249 | IF( zpelc(ji,jj,jk) > zus ) imlc(ji,jj) = jk |
---|
250 | END DO |
---|
251 | END DO |
---|
252 | END DO |
---|
253 | ! ! finite LC depth |
---|
254 | # if defined key_vectopt_loop |
---|
255 | DO jj = 1, 1 |
---|
256 | DO ji = 1, jpij ! vector opt. (forced unrolling) |
---|
257 | # else |
---|
258 | DO jj = 1, jpj |
---|
259 | DO ji = 1, jpi |
---|
260 | # endif |
---|
261 | zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj)) |
---|
262 | END DO |
---|
263 | END DO |
---|
264 | # if defined key_c1d |
---|
265 | hlc(:,:) = zhlc(:,:) * tmask(:,:,1) ! c1d configuration: save finite Langmuir Circulation depth |
---|
266 | # endif |
---|
267 | !CDIR NOVERRCHK |
---|
268 | DO jk = 2, jpkm1 !* TKE Langmuir circulation source term added to en |
---|
269 | !CDIR NOVERRCHK |
---|
270 | DO jj = 2, jpjm1 |
---|
271 | !CDIR NOVERRCHK |
---|
272 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
273 | !!gm replace here wndn by a formulation with the stress module |
---|
274 | zus = 0.016 * wndm(ji,jj) ! Stokes drift |
---|
275 | ! ! vertical velocity due to LC |
---|
276 | zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) ) |
---|
277 | zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) ) |
---|
278 | ! ! TKE Langmuir circulation source term |
---|
279 | en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) * tmask(ji,jj,jk) |
---|
280 | END DO |
---|
281 | END DO |
---|
282 | END DO |
---|
283 | ! |
---|
284 | ENDIF |
---|
285 | ! |
---|
286 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
287 | ! ! Now Turbulent kinetic energy (output in en) |
---|
288 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
289 | ! ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
290 | ! ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ). |
---|
291 | ! ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal |
---|
292 | ! |
---|
293 | DO jk = 2, jpkm1 !* Shear production at uw- and vw-points (energy conserving form) |
---|
294 | DO jj = 1, jpj ! here avmu, avmv used as workspace |
---|
295 | DO ji = 1, jpi |
---|
296 | avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) ) & |
---|
297 | & * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) ) & |
---|
298 | & / ( fse3uw_n(ji,jj,jk) & |
---|
299 | & * fse3uw_b(ji,jj,jk) ) |
---|
300 | avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) ) & |
---|
301 | & * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) ) & |
---|
302 | & / ( fse3vw_n(ji,jj,jk) & |
---|
303 | & * fse3vw_b(ji,jj,jk) ) |
---|
304 | END DO |
---|
305 | END DO |
---|
306 | END DO |
---|
307 | ! |
---|
308 | DO jk = 2, jpkm1 !* Matrix and right hand side in en |
---|
309 | DO jj = 2, jpjm1 |
---|
310 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
311 | zcof = zfact1 * tmask(ji,jj,jk) |
---|
312 | zzd_up = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) & ! upper diagonal |
---|
313 | & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk ) ) |
---|
314 | zzd_lw = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) & ! lower diagonal |
---|
315 | & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) |
---|
316 | ! ! shear prod. at w-point weightened by mask |
---|
317 | zesh2 = ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1.e0 , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) & |
---|
318 | & + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1.e0 , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) ) |
---|
319 | ! |
---|
320 | zd_up(ji,jj,jk) = zzd_up ! Matrix (zdiag, zd_up, zd_lw) |
---|
321 | zd_lw(ji,jj,jk) = zzd_lw |
---|
322 | zdiag(ji,jj,jk) = 1.e0 - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk) |
---|
323 | ! |
---|
324 | ! ! right hand side in en |
---|
325 | en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( zesh2 - avt(ji,jj,jk) * rn2(ji,jj,jk) & |
---|
326 | & + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk) ) * tmask(ji,jj,jk) |
---|
327 | END DO |
---|
328 | END DO |
---|
329 | END DO |
---|
330 | ! !* Matrix inversion from level 2 (tke prescribed at level 1) |
---|
331 | DO jk = 3, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
332 | DO jj = 2, jpjm1 |
---|
333 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
334 | zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1) |
---|
335 | END DO |
---|
336 | END DO |
---|
337 | END DO |
---|
338 | DO jj = 2, jpjm1 ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
339 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
340 | zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1) ! Surface boudary conditions on tke |
---|
341 | END DO |
---|
342 | END DO |
---|
343 | DO jk = 3, jpkm1 |
---|
344 | DO jj = 2, jpjm1 |
---|
345 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
346 | 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) |
---|
347 | END DO |
---|
348 | END DO |
---|
349 | END DO |
---|
350 | DO jj = 2, jpjm1 ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
351 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
352 | en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1) |
---|
353 | END DO |
---|
354 | END DO |
---|
355 | DO jk = jpk-2, 2, -1 |
---|
356 | DO jj = 2, jpjm1 |
---|
357 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
358 | en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk) |
---|
359 | END DO |
---|
360 | END DO |
---|
361 | END DO |
---|
362 | DO jk = 2, jpkm1 ! set the minimum value of tke |
---|
363 | DO jj = 2, jpjm1 |
---|
364 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
365 | en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * tmask(ji,jj,jk) |
---|
366 | END DO |
---|
367 | END DO |
---|
368 | END DO |
---|
369 | |
---|
370 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
371 | ! ! TKE due to surface and internal wave breaking |
---|
372 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
373 | IF( nn_etau == 1 ) THEN !* penetration throughout the water column |
---|
374 | DO jk = 2, jpkm1 |
---|
375 | DO jj = 2, jpjm1 |
---|
376 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
377 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) ) & |
---|
378 | & * ( 1.e0 - fr_i(ji,jj) ) * tmask(ji,jj,jk) |
---|
379 | END DO |
---|
380 | END DO |
---|
381 | END DO |
---|
382 | ELSEIF( nn_etau == 2 ) THEN !* act only at the base of the mixed layer (jk=nmln) |
---|
383 | DO jj = 2, jpjm1 |
---|
384 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
385 | jk = nmln(ji,jj) |
---|
386 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) ) & |
---|
387 | & * ( 1.e0 - fr_i(ji,jj) ) * tmask(ji,jj,jk) |
---|
388 | END DO |
---|
389 | END DO |
---|
390 | ENDIF |
---|
391 | ! |
---|
392 | CALL lbc_lnk( en, 'W', 1. ) ! Lateral boundary conditions (sign unchanged) |
---|
393 | ! |
---|
394 | END SUBROUTINE tke_tke |
---|
395 | |
---|
396 | |
---|
397 | SUBROUTINE tke_avn |
---|
398 | !!---------------------------------------------------------------------- |
---|
399 | !! *** ROUTINE tke_avn *** |
---|
400 | !! |
---|
401 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
402 | !! |
---|
403 | !! ** Method : At this stage, en, the now TKE, is known (computed in |
---|
404 | !! the tke_tke routine). First, the now mixing lenth is |
---|
405 | !! computed from en and the strafification (N^2), then the mixings |
---|
406 | !! coefficients are computed. |
---|
407 | !! - Mixing length : a first evaluation of the mixing lengh |
---|
408 | !! scales is: |
---|
409 | !! mxl = sqrt(2*en) / N |
---|
410 | !! where N is the brunt-vaisala frequency, with a minimum value set |
---|
411 | !! to rn_lmin (rn_lmin0) in the interior (surface) ocean. |
---|
412 | !! The mixing and dissipative length scale are bound as follow : |
---|
413 | !! nn_mxl=0 : mxl bounded by the distance to surface and bottom. |
---|
414 | !! zmxld = zmxlm = mxl |
---|
415 | !! nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl |
---|
416 | !! nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is |
---|
417 | !! less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl |
---|
418 | !! nn_mxl=3 : mxl is bounded from the surface to the bottom usings |
---|
419 | !! |d/dz(xml)|<1 to obtain lup, and from the bottom to |
---|
420 | !! the surface to obtain ldown. the resulting length |
---|
421 | !! scales are: |
---|
422 | !! zmxld = sqrt( lup * ldown ) |
---|
423 | !! zmxlm = min ( lup , ldown ) |
---|
424 | !! - Vertical eddy viscosity and diffusivity: |
---|
425 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
426 | !! avt = max( avmb, pdlr * avm ) |
---|
427 | !! with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise. |
---|
428 | !! |
---|
429 | !! ** Action : - avt : now vertical eddy diffusivity (w-point) |
---|
430 | !! - avmu, avmv : now vertical eddy viscosity at uw- and vw-points |
---|
431 | !!---------------------------------------------------------------------- |
---|
432 | USE oce, zmpdl => ua ! use ua as workspace |
---|
433 | USE oce, zmxlm => va ! use va as workspace |
---|
434 | USE oce, zmxld => ta ! use ta as workspace |
---|
435 | !! |
---|
436 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
437 | REAL(wp) :: zrn2, zraug ! temporary scalars |
---|
438 | REAL(wp) :: zdku ! - - |
---|
439 | REAL(wp) :: zdkv ! - - |
---|
440 | REAL(wp) :: zcoef, zav ! - - |
---|
441 | REAL(wp) :: zpdlr, zri, zsqen ! - - |
---|
442 | REAL(wp) :: zemxl, zemlm, zemlp ! - - |
---|
443 | !!-------------------------------------------------------------------- |
---|
444 | |
---|
445 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
446 | ! ! Mixing length |
---|
447 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
448 | ! |
---|
449 | ! !* Buoyancy length scale: l=sqrt(2*e/n**2) |
---|
450 | ! |
---|
451 | IF( ln_mxl0 ) THEN ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g) |
---|
452 | !!gm this should be useless |
---|
453 | zmxlm(:,:,1) = 0.e0 |
---|
454 | !!gm end |
---|
455 | zraug = vkarmn * 2.e5 / ( rau0 * grav ) |
---|
456 | DO jj = 2, jpjm1 |
---|
457 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
458 | zmxlm(ji,jj,1) = MAX( rn_lmin0, zraug * taum(ji,jj) ) |
---|
459 | END DO |
---|
460 | END DO |
---|
461 | ELSE ! surface set to the minimum value |
---|
462 | zmxlm(:,:,1) = rn_lmin0 |
---|
463 | ENDIF |
---|
464 | zmxlm(:,:,jpk) = rn_lmin ! bottom set to the interior minium value |
---|
465 | ! |
---|
466 | !CDIR NOVERRCHK |
---|
467 | DO jk = 2, jpkm1 ! interior value : l=sqrt(2*e/n**2) |
---|
468 | !CDIR NOVERRCHK |
---|
469 | DO jj = 2, jpjm1 |
---|
470 | !CDIR NOVERRCHK |
---|
471 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
472 | zrn2 = MAX( rn2(ji,jj,jk), rsmall ) |
---|
473 | zmxlm(ji,jj,jk) = MAX( rn_lmin, SQRT( 2. * en(ji,jj,jk) / zrn2 ) ) |
---|
474 | END DO |
---|
475 | END DO |
---|
476 | END DO |
---|
477 | ! |
---|
478 | !!gm CAUTION: to be added here the bottom boundary condition on zmxlm |
---|
479 | ! |
---|
480 | ! !* Physical limits for the mixing length |
---|
481 | ! |
---|
482 | zmxld(:,:, 1 ) = zmxlm(:,:,1) ! surface set to the minimum value |
---|
483 | zmxld(:,:,jpk) = rn_lmin ! bottom set to the minimum value |
---|
484 | ! |
---|
485 | SELECT CASE ( nn_mxl ) |
---|
486 | ! |
---|
487 | CASE ( 0 ) ! bounded by the distance to surface and bottom |
---|
488 | DO jk = 2, jpkm1 |
---|
489 | DO jj = 2, jpjm1 |
---|
490 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
491 | zemxl = MIN( fsdepw(ji,jj,jk), zmxlm(ji,jj,jk), & |
---|
492 | & fsdepw(ji,jj,mbathy(ji,jj)) - fsdepw(ji,jj,jk) ) |
---|
493 | zmxlm(ji,jj,jk) = zemxl |
---|
494 | zmxld(ji,jj,jk) = zemxl |
---|
495 | END DO |
---|
496 | END DO |
---|
497 | END DO |
---|
498 | ! |
---|
499 | CASE ( 1 ) ! bounded by the vertical scale factor |
---|
500 | DO jk = 2, jpkm1 |
---|
501 | DO jj = 2, jpjm1 |
---|
502 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
503 | zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
504 | zmxlm(ji,jj,jk) = zemxl |
---|
505 | zmxld(ji,jj,jk) = zemxl |
---|
506 | END DO |
---|
507 | END DO |
---|
508 | END DO |
---|
509 | ! |
---|
510 | CASE ( 2 ) ! |dk[xml]| bounded by e3t : |
---|
511 | DO jk = 2, jpkm1 ! from the surface to the bottom : |
---|
512 | DO jj = 2, jpjm1 |
---|
513 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
514 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
515 | END DO |
---|
516 | END DO |
---|
517 | END DO |
---|
518 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : |
---|
519 | DO jj = 2, jpjm1 |
---|
520 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
521 | zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
522 | zmxlm(ji,jj,jk) = zemxl |
---|
523 | zmxld(ji,jj,jk) = zemxl |
---|
524 | END DO |
---|
525 | END DO |
---|
526 | END DO |
---|
527 | ! |
---|
528 | CASE ( 3 ) ! lup and ldown, |dk[xml]| bounded by e3t : |
---|
529 | DO jk = 2, jpkm1 ! from the surface to the bottom : lup |
---|
530 | DO jj = 2, jpjm1 |
---|
531 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
532 | zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
533 | END DO |
---|
534 | END DO |
---|
535 | END DO |
---|
536 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : ldown |
---|
537 | DO jj = 2, jpjm1 |
---|
538 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
539 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
540 | END DO |
---|
541 | END DO |
---|
542 | END DO |
---|
543 | !CDIR NOVERRCHK |
---|
544 | DO jk = 2, jpkm1 |
---|
545 | !CDIR NOVERRCHK |
---|
546 | DO jj = 2, jpjm1 |
---|
547 | !CDIR NOVERRCHK |
---|
548 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
549 | zemlm = MIN ( zmxld(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
550 | zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) ) |
---|
551 | zmxlm(ji,jj,jk) = zemlm |
---|
552 | zmxld(ji,jj,jk) = zemlp |
---|
553 | END DO |
---|
554 | END DO |
---|
555 | END DO |
---|
556 | ! |
---|
557 | END SELECT |
---|
558 | ! |
---|
559 | # if defined key_c1d |
---|
560 | e_dis(:,:,:) = zmxld(:,:,:) ! c1d configuration : save mixing and dissipation turbulent length scales |
---|
561 | e_mix(:,:,:) = zmxlm(:,:,:) |
---|
562 | # endif |
---|
563 | |
---|
564 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
565 | ! ! Vertical eddy viscosity and diffusivity (avmu, avmv, avt) |
---|
566 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
567 | !CDIR NOVERRCHK |
---|
568 | DO jk = 1, jpkm1 !* vertical eddy viscosity & diffivity at w-points |
---|
569 | !CDIR NOVERRCHK |
---|
570 | DO jj = 2, jpjm1 |
---|
571 | !CDIR NOVERRCHK |
---|
572 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
573 | zsqen = SQRT( en(ji,jj,jk) ) |
---|
574 | zav = rn_ediff * zmxlm(ji,jj,jk) * zsqen |
---|
575 | avm (ji,jj,jk) = MAX( zav, avmb(jk) ) * tmask(ji,jj,jk) |
---|
576 | avt (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk) |
---|
577 | dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk) |
---|
578 | END DO |
---|
579 | END DO |
---|
580 | END DO |
---|
581 | CALL lbc_lnk( avm, 'W', 1. ) ! Lateral boundary conditions (sign unchanged) |
---|
582 | ! |
---|
583 | DO jk = 2, jpkm1 !* vertical eddy viscosity at u- and v-points |
---|
584 | DO jj = 2, jpjm1 |
---|
585 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
586 | avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj ,jk) ) * umask(ji,jj,jk) |
---|
587 | avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji ,jj+1,jk) ) * vmask(ji,jj,jk) |
---|
588 | END DO |
---|
589 | END DO |
---|
590 | END DO |
---|
591 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions |
---|
592 | ! |
---|
593 | IF( nn_pdl == 1 ) THEN !* Prandtl number case: update avt |
---|
594 | DO jk = 2, jpkm1 |
---|
595 | DO jj = 2, jpjm1 |
---|
596 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
597 | zcoef = 0.5 / ( fse3w(ji,jj,jk) * fse3w(ji,jj,jk) ) |
---|
598 | ! ! shear |
---|
599 | zdku = avmu(ji-1,jj,jk) * ( un(ji-1,jj,jk-1) - un(ji-1,jj,jk) ) * ( ub(ji-1,jj,jk-1) - ub(ji-1,jj,jk) ) & |
---|
600 | & + avmu(ji ,jj,jk) * ( un(ji ,jj,jk-1) - un(ji ,jj,jk) ) * ( ub(ji ,jj,jk-1) - ub(ji ,jj,jk) ) |
---|
601 | zdkv = avmv(ji,jj-1,jk) * ( vn(ji,jj-1,jk-1) - vn(ji,jj-1,jk) ) * ( vb(ji,jj-1,jk-1) - vb(ji,jj-1,jk) ) & |
---|
602 | & + avmv(ji,jj ,jk) * ( vn(ji,jj ,jk-1) - vn(ji,jj ,jk) ) * ( vb(ji,jj ,jk-1) - vb(ji,jj ,jk) ) |
---|
603 | ! ! local Richardson number |
---|
604 | zri = MAX( rn2b(ji,jj,jk), 0. ) * avm(ji,jj,jk) / (zcoef * (zdku + zdkv + rn_bshear ) ) |
---|
605 | zpdlr = MAX( 0.1, 0.2 / MAX( 0.2 , zri ) ) |
---|
606 | !!gm and even better with the use of the "true" ri_crit=0.22222... (this change the results!) |
---|
607 | !!gm zpdlr = MAX( 0.1, ri_crit / MAX( ri_crit , zri ) ) |
---|
608 | avt(ji,jj,jk) = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk) |
---|
609 | # if defined key_c1d |
---|
610 | e_pdl(ji,jj,jk) = zpdlr * tmask(ji,jj,jk) ! c1d configuration : save masked Prandlt number |
---|
611 | e_ric(ji,jj,jk) = zri * tmask(ji,jj,jk) ! c1d config. : save Ri |
---|
612 | # endif |
---|
613 | END DO |
---|
614 | END DO |
---|
615 | END DO |
---|
616 | ENDIF |
---|
617 | CALL lbc_lnk( avt, 'W', 1. ) ! Lateral boundary conditions on avt (sign unchanged) |
---|
618 | |
---|
619 | IF(ln_ctl) THEN |
---|
620 | CALL prt_ctl( tab3d_1=en , clinfo1=' tke - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk) |
---|
621 | CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke - u: ', mask1=umask, & |
---|
622 | & tab3d_2=avmv, clinfo2= ' v: ', mask2=vmask, ovlap=1, kdim=jpk ) |
---|
623 | ENDIF |
---|
624 | ! |
---|
625 | END SUBROUTINE tke_avn |
---|
626 | |
---|
627 | |
---|
628 | SUBROUTINE tke_init |
---|
629 | !!---------------------------------------------------------------------- |
---|
630 | !! *** ROUTINE tke_init *** |
---|
631 | !! |
---|
632 | !! ** Purpose : Initialization of the vertical eddy diffivity and |
---|
633 | !! viscosity when using a tke turbulent closure scheme |
---|
634 | !! |
---|
635 | !! ** Method : Read the namzdf_tke namelist and check the parameters |
---|
636 | !! called at the first timestep (nit000) |
---|
637 | !! |
---|
638 | !! ** input : Namlist namzdf_tke |
---|
639 | !! |
---|
640 | !! ** Action : Increase by 1 the nstop flag is setting problem encounter |
---|
641 | !!---------------------------------------------------------------------- |
---|
642 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
643 | !! |
---|
644 | NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb, rn_emin, & |
---|
645 | & rn_emin0, rn_bshear, nn_mxl, ln_mxl0, & |
---|
646 | & rn_lmin , rn_lmin0 , nn_pdl, nn_etau, & |
---|
647 | & nn_htau , rn_efr , ln_lc , rn_lc |
---|
648 | !!---------------------------------------------------------------------- |
---|
649 | |
---|
650 | REWIND ( numnam ) !* Read Namelist namzdf_tke : Turbulente Kinetic Energy |
---|
651 | READ ( numnam, namzdf_tke ) |
---|
652 | |
---|
653 | ri_cri = 2. / ( 2. + rn_ediss / rn_ediff ) ! resulting critical Richardson number |
---|
654 | |
---|
655 | IF(lwp) THEN !* Control print |
---|
656 | WRITE(numout,*) |
---|
657 | WRITE(numout,*) 'zdf_tke : tke turbulent closure scheme - initialisation' |
---|
658 | WRITE(numout,*) '~~~~~~~~' |
---|
659 | WRITE(numout,*) ' Namelist namzdf_tke : set tke mixing parameters' |
---|
660 | WRITE(numout,*) ' coef. to compute avt rn_ediff = ', rn_ediff |
---|
661 | WRITE(numout,*) ' Kolmogoroff dissipation coef. rn_ediss = ', rn_ediss |
---|
662 | WRITE(numout,*) ' tke surface input coef. rn_ebb = ', rn_ebb |
---|
663 | WRITE(numout,*) ' minimum value of tke rn_emin = ', rn_emin |
---|
664 | WRITE(numout,*) ' surface minimum value of tke rn_emin0 = ', rn_emin0 |
---|
665 | WRITE(numout,*) ' background shear (>0) rn_bshear= ', rn_bshear |
---|
666 | WRITE(numout,*) ' mixing length type nn_mxl = ', nn_mxl |
---|
667 | WRITE(numout,*) ' prandl number flag nn_pdl = ', nn_pdl |
---|
668 | WRITE(numout,*) ' surface mixing length = F(stress) or not ln_mxl0 = ', ln_mxl0 |
---|
669 | WRITE(numout,*) ' surface mixing length minimum value rn_lmin0 = ', rn_lmin0 |
---|
670 | WRITE(numout,*) ' interior mixing length minimum value rn_lmin0 = ', rn_lmin0 |
---|
671 | WRITE(numout,*) ' test param. to add tke induced by wind nn_etau = ', nn_etau |
---|
672 | WRITE(numout,*) ' flag for computation of exp. tke profile nn_htau = ', nn_htau |
---|
673 | WRITE(numout,*) ' % of rn_emin0 which pene. the thermocline rn_efr = ', rn_efr |
---|
674 | WRITE(numout,*) ' flag to take into acc. Langmuir circ. ln_lc = ', ln_lc |
---|
675 | WRITE(numout,*) ' coef to compute verticla velocity of LC rn_lc = ', rn_lc |
---|
676 | WRITE(numout,*) |
---|
677 | WRITE(numout,*) ' critical Richardson nb with your parameters ri_cri = ', ri_cri |
---|
678 | ENDIF |
---|
679 | |
---|
680 | ! !* Check of some namelist values |
---|
681 | IF( nn_mxl < 0 .OR. nn_mxl > 3 ) CALL ctl_stop( 'bad flag: nn_mxl is 0, 1 or 2 ' ) |
---|
682 | IF( nn_pdl < 0 .OR. nn_pdl > 1 ) CALL ctl_stop( 'bad flag: nn_pdl is 0 or 1 ' ) |
---|
683 | IF( nn_htau < 0 .OR. nn_htau > 1 ) CALL ctl_stop( 'bad flag: nn_htau is 0 or 1 ' ) |
---|
684 | 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 ' ) |
---|
685 | |
---|
686 | IF( nn_etau == 2 ) CALL zdf_mxl( nit000 ) ! Initialization of nmln |
---|
687 | |
---|
688 | ! !* depth of penetration of surface tke |
---|
689 | IF( nn_etau /= 0 ) THEN |
---|
690 | SELECT CASE( nn_htau ) ! Choice of the depth of penetration |
---|
691 | CASE( 0 ) ! constant depth penetration (here 10 meters) |
---|
692 | htau(:,:) = 10.e0 |
---|
693 | CASE( 1 ) ! F(latitude) : 0.5m to 30m at high lat. |
---|
694 | DO jj = 1, jpj |
---|
695 | DO ji = 1, jpi |
---|
696 | htau(ji,jj) = MAX( 0.5, 3./4. * MIN( 40., 60.*ABS( SIN( rpi/180. * gphit(ji,jj) ) ) ) ) |
---|
697 | END DO |
---|
698 | END DO |
---|
699 | END SELECT |
---|
700 | ENDIF |
---|
701 | |
---|
702 | ! !* set vertical eddy coef. to the background value |
---|
703 | DO jk = 1, jpk |
---|
704 | avt (:,:,jk) = avtb(jk) * tmask(:,:,jk) |
---|
705 | avm (:,:,jk) = avmb(jk) * tmask(:,:,jk) |
---|
706 | avmu(:,:,jk) = avmb(jk) * umask(:,:,jk) |
---|
707 | avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk) |
---|
708 | END DO |
---|
709 | dissl(:,:,:) = 1.e-12 |
---|
710 | ! !* read or initialize all required files |
---|
711 | CALL tke_rst( nit000, 'READ' ) |
---|
712 | ! |
---|
713 | END SUBROUTINE tke_init |
---|
714 | |
---|
715 | |
---|
716 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
717 | !!--------------------------------------------------------------------- |
---|
718 | !! *** ROUTINE tke_rst *** |
---|
719 | !! |
---|
720 | !! ** Purpose : Read or write TKE file (en) in restart file |
---|
721 | !! |
---|
722 | !! ** Method : use of IOM library |
---|
723 | !! if the restart does not contain TKE, en is either |
---|
724 | !! set to rn_emin or recomputed |
---|
725 | !!---------------------------------------------------------------------- |
---|
726 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
727 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
728 | ! |
---|
729 | INTEGER :: jit, jk ! dummy loop indices |
---|
730 | INTEGER :: id1, id2, id3, id4, id5, id6 |
---|
731 | !!---------------------------------------------------------------------- |
---|
732 | ! |
---|
733 | IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise |
---|
734 | ! ! --------------- |
---|
735 | IF( ln_rstart ) THEN !* Read the restart file |
---|
736 | id1 = iom_varid( numror, 'en' , ldstop = .FALSE. ) |
---|
737 | id2 = iom_varid( numror, 'avt' , ldstop = .FALSE. ) |
---|
738 | id3 = iom_varid( numror, 'avm' , ldstop = .FALSE. ) |
---|
739 | id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. ) |
---|
740 | id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. ) |
---|
741 | id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. ) |
---|
742 | ! |
---|
743 | IF( id1 > 0 ) THEN ! 'en' exists |
---|
744 | CALL iom_get( numror, jpdom_autoglo, 'en', en ) |
---|
745 | IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN ! all required arrays exist |
---|
746 | CALL iom_get( numror, jpdom_autoglo, 'avt' , avt ) |
---|
747 | CALL iom_get( numror, jpdom_autoglo, 'avm' , avm ) |
---|
748 | CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu ) |
---|
749 | CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv ) |
---|
750 | CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl ) |
---|
751 | ELSE ! one at least array is missing |
---|
752 | CALL tke_avn ! compute avt, avm, avmu, avmv and dissl (approximation) |
---|
753 | ENDIF |
---|
754 | ELSE ! No TKE array found: initialisation |
---|
755 | IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop' |
---|
756 | en (:,:,:) = rn_emin * tmask(:,:,:) |
---|
757 | CALL tke_avn ! recompute avt, avm, avmu, avmv and dissl (approximation) |
---|
758 | DO jit = nit000 + 1, nit000 + 10 ; CALL zdf_tke( jit ) ; END DO |
---|
759 | ENDIF |
---|
760 | ELSE !* Start from rest |
---|
761 | en(:,:,:) = rn_emin * tmask(:,:,:) |
---|
762 | DO jk = 1, jpk ! set the Kz to the background value |
---|
763 | avt (:,:,jk) = avtb(jk) * tmask(:,:,jk) |
---|
764 | avm (:,:,jk) = avmb(jk) * tmask(:,:,jk) |
---|
765 | avmu(:,:,jk) = avmb(jk) * umask(:,:,jk) |
---|
766 | avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk) |
---|
767 | END DO |
---|
768 | ENDIF |
---|
769 | ! |
---|
770 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file |
---|
771 | ! ! ------------------- |
---|
772 | IF(lwp) WRITE(numout,*) '---- tke-rst ----' |
---|
773 | CALL iom_rstput( kt, nitrst, numrow, 'en' , en ) |
---|
774 | CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt ) |
---|
775 | CALL iom_rstput( kt, nitrst, numrow, 'avm' , avm ) |
---|
776 | CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu ) |
---|
777 | CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv ) |
---|
778 | CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl ) |
---|
779 | ! |
---|
780 | ENDIF |
---|
781 | ! |
---|
782 | END SUBROUTINE tke_rst |
---|
783 | |
---|
784 | #else |
---|
785 | !!---------------------------------------------------------------------- |
---|
786 | !! Dummy module : NO TKE scheme |
---|
787 | !!---------------------------------------------------------------------- |
---|
788 | LOGICAL, PUBLIC, PARAMETER :: lk_zdftke = .FALSE. !: TKE flag |
---|
789 | CONTAINS |
---|
790 | SUBROUTINE zdf_tke( kt ) ! Empty routine |
---|
791 | WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt |
---|
792 | END SUBROUTINE zdf_tke |
---|
793 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
794 | CHARACTER(len=*) :: cdrw |
---|
795 | WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr |
---|
796 | END SUBROUTINE tke_rst |
---|
797 | #endif |
---|
798 | |
---|
799 | !!====================================================================== |
---|
800 | END MODULE zdftke |
---|