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