1 |
! $Header: $ |
2 |
! $Name: $ |
3 |
|
4 |
MODULE MY25_TURB_MOD |
5 |
|
6 |
!======================================================================= |
7 |
! MELLOR-YAMADA LEVEL 2.5 TURBULENCE CLOSURE SCHEME - GFDL VERSION ! |
8 |
!======================================================================= |
9 |
|
10 |
! use Fms_Mod, ONLY: FILE_EXIST, OPEN_NAMELIST_FILE, ERROR_MESG, FATAL,& |
11 |
! read_data, write_data, CLOSE_FILE,& |
12 |
! check_nml_error,mpp_pe, mpp_root_pe, & |
13 |
! write_version_number, stdlog, open_file |
14 |
! use Tridiagonal_Mod, ONLY: TRI_INVERT, CLOSE_TRIDIAGONAL |
15 |
use constants_mod, only: grav, vonkarm |
16 |
use monin_obukhov_mod, only : mo_diff |
17 |
|
18 |
!--------------------------------------------------------------------- |
19 |
implicit none |
20 |
private |
21 |
!--------------------------------------------------------------------- |
22 |
|
23 |
public :: MY25_TURB, MY25_TURB_INIT, MY25_TURB_END, TKE_SURF |
24 |
|
25 |
!--------------------------------------------------------------------- |
26 |
! --- GLOBAL STORAGE |
27 |
!--------------------------------------------------------------------- |
28 |
|
29 |
real, public, allocatable, dimension(:,:,:) :: TKE |
30 |
|
31 |
!--------------------------------------------------------------------- |
32 |
|
33 |
character(len=128) :: version = '$Id: my25_turb_mod.F90,v 1.1 2012/09/11 03:53:05 jmc Exp $' |
34 |
character(len=128) :: tag = '$Name: $' |
35 |
|
36 |
logical :: do_init = .true. |
37 |
logical :: init_tke |
38 |
integer :: num_total_pts, pts_done |
39 |
|
40 |
!--------------------------------------------------------------------- |
41 |
! --- CONSTANTS |
42 |
!--------------------------------------------------------------------- |
43 |
|
44 |
real :: aa1, aa2, bb1, bb2, ccc |
45 |
real :: ckm1, ckm2, ckm3, ckm4, ckm5, ckm6, ckm7, ckm8 |
46 |
real :: ckh1, ckh2, ckh3, ckh4 |
47 |
real :: cvfq1, cvfq2, bcq |
48 |
|
49 |
real, parameter :: aa1_old = 0.78 |
50 |
real, parameter :: aa2_old = 0.79 |
51 |
real, parameter :: bb1_old = 15.0 |
52 |
real, parameter :: bb2_old = 8.0 |
53 |
real, parameter :: ccc_old = 0.056 |
54 |
real, parameter :: aa1_new = 0.92 |
55 |
real, parameter :: aa2_new = 0.74 |
56 |
real, parameter :: bb1_new = 16.0 |
57 |
real, parameter :: bb2_new = 10.0 |
58 |
real, parameter :: ccc_new = 0.08 |
59 |
real, parameter :: cc1 = 0.27 |
60 |
real, parameter :: t00 = 2.7248e2 |
61 |
real, parameter :: small = 1.0e-10 |
62 |
|
63 |
!--------------------------------------------------------------------- |
64 |
! --- NAMELIST |
65 |
!--------------------------------------------------------------------- |
66 |
|
67 |
real :: TKEmax = 5.0 |
68 |
real :: TKEmin = 0.0 |
69 |
real :: el0max = 1.0e6 |
70 |
real :: el0min = 0.0 |
71 |
real :: alpha_land = 0.10 |
72 |
real :: alpha_sea = 0.10 |
73 |
real :: akmax = 1.0e4 |
74 |
real :: akmin_land = 5.0 |
75 |
real :: akmin_sea = 0.0 |
76 |
integer :: nk_lim = 2 |
77 |
integer :: init_iters = 20 |
78 |
logical :: do_thv_stab = .true. |
79 |
logical :: use_old_cons = .false. |
80 |
real :: kcrit = 0.01 |
81 |
|
82 |
NAMELIST / my25_turb_nml / & |
83 |
TKEmax, TKEmin, init_iters, & |
84 |
akmax, akmin_land, akmin_sea, nk_lim, & |
85 |
el0max, el0min, alpha_land, alpha_sea, & |
86 |
do_thv_stab, use_old_cons, & |
87 |
kcrit |
88 |
|
89 |
!--------------------------------------------------------------------- |
90 |
|
91 |
contains |
92 |
|
93 |
!####################################################################### |
94 |
|
95 |
SUBROUTINE MY25_TURB( delt, fracland, phalf, pfull, theta, & |
96 |
um, vm, zhalf, zfull, z0, & |
97 |
TKE, el0, el, akm, akh, & |
98 |
mask, kbot, ustar, bstar, h ) |
99 |
|
100 |
!======================================================================= |
101 |
!--------------------------------------------------------------------- |
102 |
! Arguments (Intent in) |
103 |
! delt - Time step in seconds |
104 |
! fracland - Fractional amount of land beneath a grid box |
105 |
! phalf - Pressure at half levels |
106 |
! pfull - Pressure at full levels |
107 |
! theta - Potential temperature |
108 |
! um, vm - Wind components |
109 |
! zhalf - Height at half levels |
110 |
! zfull - Height at full levels |
111 |
! z0 - Roughness length |
112 |
! mask - OPTIONAL; floating point mask (0. or 1.) designating |
113 |
! where data is present |
114 |
! kbot - OPTIONAL;lowest model level index (integer); |
115 |
! at levels > kbot, mask = 0. |
116 |
! ustar - OPTIONAL:friction velocity (m/sec) |
117 |
! bstar - OPTIONAL:buoyancy scale (m/sec**2) |
118 |
!--------------------------------------------------------------------- |
119 |
real, intent(in) :: delt |
120 |
real, intent(in), dimension(:,:) :: fracland, z0 |
121 |
real, intent(in), dimension(:,:,:) :: phalf, pfull, zhalf, zfull |
122 |
real, intent(in), dimension(:,:,:) :: um, vm, theta |
123 |
|
124 |
integer, intent(in), OPTIONAL, dimension(:,:) :: kbot |
125 |
real, intent(in), OPTIONAL, dimension(:,:,:) :: mask |
126 |
real, intent(in), OPTIONAL, dimension(:,:) :: ustar, bstar |
127 |
|
128 |
!--------------------------------------------------------------------- |
129 |
! Arguments (Intent in/out) |
130 |
! TKE - turbulent kinetic energy |
131 |
!--------------------------------------------------------------------- |
132 |
real, intent(inout), dimension(:,:,:) :: TKE |
133 |
|
134 |
!--------------------------------------------------------------------- |
135 |
! Arguments (Intent out) |
136 |
! el0 - characteristic length scale |
137 |
! el - master length scale |
138 |
! akm - mixing coefficient for momentum |
139 |
! akh - mixing coefficient for heat and moisture |
140 |
! h - OPTIONAL, diagnosed depth of planetary boundary |
141 |
! layer (m) |
142 |
!--------------------------------------------------------------------- |
143 |
real, intent(out), dimension(:,:) :: el0 |
144 |
real, intent(out), dimension(:,:,:) :: akm, akh, el |
145 |
real, intent(out), OPTIONAL, dimension(:,:) :: h |
146 |
|
147 |
!--------------------------------------------------------------------- |
148 |
! (Intent local) |
149 |
!--------------------------------------------------------------------- |
150 |
integer :: ix, jx, kx, i, j, k |
151 |
integer :: kxp, kxm, klim, it, itermax |
152 |
real :: cvfqdt, dvfqdt |
153 |
|
154 |
real, dimension(SIZE(um,1),SIZE(um,2)) :: zsfc, x1, x2, akmin |
155 |
|
156 |
real, dimension(SIZE(um,1),SIZE(um,2),SIZE(um,3)-1) :: & |
157 |
dsdzh, shear, buoync, qm2, qm3, qm4, el2, & |
158 |
aaa, bbb, ccc, ddd, & |
159 |
xxm1, xxm2, xxm3, xxm4, xxm5 |
160 |
|
161 |
real, dimension(SIZE(um,1),SIZE(um,2),SIZE(um,3)) :: & |
162 |
dsdz, qm, xx1, xx2 |
163 |
|
164 |
!==================================================================== |
165 |
|
166 |
! --- Check to see if MY25_TURB has been initialized |
167 |
!del if( do_init ) CALL ERROR_MESG( ' MY25_TURB', & |
168 |
!del ' MY25_TURB_INIT has not been called',& |
169 |
!del FATAL ) |
170 |
|
171 |
! --- Set dimensions etc |
172 |
ix = SIZE( um, 1 ) |
173 |
jx = SIZE( um, 2 ) |
174 |
kx = SIZE( um, 3 ) |
175 |
kxp = kx + 1 |
176 |
kxm = kx - 1 |
177 |
|
178 |
!==================================================================== |
179 |
! --- SURFACE HEIGHT |
180 |
!==================================================================== |
181 |
|
182 |
if( PRESENT( kbot ) ) then |
183 |
do j = 1,jx |
184 |
do i = 1,ix |
185 |
k = kbot(i,j) + 1 |
186 |
zsfc(i,j) = zhalf(i,j,k) |
187 |
end do |
188 |
end do |
189 |
else |
190 |
zsfc(:,:) = zhalf(:,:,kxp) |
191 |
endif |
192 |
|
193 |
!==================================================================== |
194 |
! --- D( )/DZ OPERATORS: AT FULL LEVELS & AT HALF LEVELS |
195 |
!==================================================================== |
196 |
|
197 |
dsdz(:,:,1:kx) = 1.0 / ( zhalf(:,:,2:kxp) - zhalf(:,:,1:kx) ) |
198 |
dsdzh(:,:,1:kxm) = 1.0 / ( zfull(:,:,2:kx) - zfull(:,:,1:kxm) ) |
199 |
|
200 |
!==================================================================== |
201 |
! --- WIND SHEAR |
202 |
!==================================================================== |
203 |
|
204 |
xxm1(:,:,1:kxm) = dsdzh(:,:,1:kxm)*( um(:,:,2:kx) - um(:,:,1:kxm) ) |
205 |
xxm2(:,:,1:kxm) = dsdzh(:,:,1:kxm)*( vm(:,:,2:kx) - vm(:,:,1:kxm) ) |
206 |
|
207 |
shear = xxm1 * xxm1 + xxm2 * xxm2 |
208 |
|
209 |
!==================================================================== |
210 |
! --- BUOYANCY |
211 |
!==================================================================== |
212 |
|
213 |
xxm1(:,:,1:kxm) = theta(:,:,2:kx) - theta(:,:,1:kxm) |
214 |
|
215 |
if( do_thv_stab ) then |
216 |
xxm2(:,:,1:kxm) = 0.5*( theta(:,:,2:kx) + theta(:,:,1:kxm) ) |
217 |
else |
218 |
xxm2(:,:,1:kxm) = t00 |
219 |
end if |
220 |
|
221 |
buoync = grav * dsdzh * xxm1 / xxm2 |
222 |
|
223 |
!==================================================================== |
224 |
! --- MASK OUT UNDERGROUND VALUES FOR ETA COORDINATE |
225 |
!==================================================================== |
226 |
|
227 |
if( PRESENT( mask ) ) then |
228 |
where (mask(:,:,2:kx) < 0.1) ! assume mask(k=1) always eq 1 ? |
229 |
TKE(:,:,3:kxp) = 0.0 |
230 |
dsdz(:,:,2:kx ) = 0.0 |
231 |
dsdzh(:,:,1:kxm) = 0.0 |
232 |
shear(:,:,1:kxm) = 0.0 |
233 |
buoync(:,:,1:kxm) = 0.0 |
234 |
endwhere |
235 |
endif |
236 |
|
237 |
!==================================================================== |
238 |
! --- SET ITERATION LOOP IF INITALIZING TKE |
239 |
!==================================================================== |
240 |
|
241 |
itermax = 1 |
242 |
if (init_tke) then |
243 |
itermax = init_iters |
244 |
pts_done = pts_done + ix*jx |
245 |
if (pts_done >= num_total_pts) init_tke = .false. |
246 |
endif |
247 |
|
248 |
! $$$$$$$$$$$$$$$$$ |
249 |
do it = 1,itermax |
250 |
! $$$$$$$$$$$$$$$$$ |
251 |
|
252 |
!==================================================================== |
253 |
! --- SOME TKE STUFF |
254 |
!==================================================================== |
255 |
|
256 |
xx1(:,:,1:kx) = 2.0 * TKE(:,:,2:kxp) |
257 |
where (xx1(:,:,1:kx) > 0.0) |
258 |
qm(:,:,1:kx) = SQRT( xx1(:,:,1:kx) ) |
259 |
elsewhere |
260 |
qm(:,:,1:kx) = 0.0 |
261 |
endwhere |
262 |
|
263 |
qm2(:,:,1:kxm) = xx1(:,:,1:kxm) |
264 |
qm3(:,:,1:kxm) = qm(:,:,1:kxm) * qm2(:,:,1:kxm) |
265 |
qm4(:,:,1:kxm) = qm2(:,:,1:kxm) * qm2(:,:,1:kxm) |
266 |
|
267 |
!==================================================================== |
268 |
! --- CHARACTERISTIC LENGTH SCALE |
269 |
!==================================================================== |
270 |
|
271 |
xx1(:,:,1:kxm) = qm(:,:,1:kxm)*( pfull(:,:,2:kx) - pfull(:,:,1:kxm) ) |
272 |
do k = 1, kxm |
273 |
xx2(:,:,k) = xx1(:,:,k) * ( zhalf(:,:,k+1) - zsfc(:,:) ) |
274 |
end do |
275 |
|
276 |
if( PRESENT( kbot ) ) then |
277 |
xx1(:,:,kx) = 0.0 |
278 |
xx2(:,:,kx) = 0.0 |
279 |
do j = 1,jx |
280 |
do i = 1,ix |
281 |
k = kbot(i,j) |
282 |
xx1(i,j,k) = qm(i,j,k) * ( phalf(i,j,k+1) - pfull(i,j,k) ) |
283 |
xx2(i,j,k) = xx1(i,j,k) * z0(i,j) |
284 |
end do |
285 |
end do |
286 |
else |
287 |
xx1(:,:,kx) = qm(:,:,kx) * ( phalf(:,:,kxp) - pfull(:,:,kx) ) |
288 |
xx2(:,:,kx) = xx1(:,:,kx) * z0(:,:) |
289 |
endif |
290 |
|
291 |
if (PRESENT(mask)) then |
292 |
x1 = SUM( xx1, 3, mask=mask.gt.0.1 ) |
293 |
x2 = SUM( xx2, 3, mask=mask.gt.0.1 ) |
294 |
else |
295 |
x1 = SUM( xx1, 3 ) |
296 |
x2 = SUM( xx2, 3 ) |
297 |
endif |
298 |
|
299 |
!---- should never be equal to zero ---- |
300 |
!del if (count(x1 <= 0.0) > 0) CALL ERROR_MESG( ' MY25_TURB', & |
301 |
!del 'divid by zero, x1 <= 0.0', FATAL) |
302 |
el0 = x2 / x1 |
303 |
el0 = el0 * (alpha_land*fracland + alpha_sea*(1.-fracland)) |
304 |
|
305 |
el0 = MIN( el0, el0max ) |
306 |
el0 = MAX( el0, el0min ) |
307 |
|
308 |
!==================================================================== |
309 |
! --- MASTER LENGTH SCALE |
310 |
!==================================================================== |
311 |
|
312 |
do k = 1, kxm |
313 |
xx1(:,:,k) = vonkarm * ( zhalf(:,:,k+1) - zsfc(:,:) ) |
314 |
end do |
315 |
|
316 |
x1(:,:) = vonkarm * z0(:,:) |
317 |
|
318 |
if( PRESENT( kbot ) ) then |
319 |
do j = 1,jx |
320 |
do i = 1,ix |
321 |
do k = kbot(i,j), kx |
322 |
xx1(i,j,k) = x1(i,j) |
323 |
end do |
324 |
end do |
325 |
end do |
326 |
else |
327 |
xx1(:,:,kx) = x1(:,:) |
328 |
endif |
329 |
|
330 |
do k = 1,kx |
331 |
el(:,:,k+1) = xx1(:,:,k) / ( 1.0 + xx1(:,:,k) / el0(:,:) ) |
332 |
end do |
333 |
el(:,:,1) = el0(:,:) |
334 |
|
335 |
el2(:,:,1:kxm) = el(:,:,2:kx) * el(:,:,2:kx) |
336 |
|
337 |
!==================================================================== |
338 |
! --- MIXING COEFFICIENTS |
339 |
!==================================================================== |
340 |
|
341 |
xxm3(:,:,1:kxm) = el2(:,:,1:kxm)*buoync(:,:,1:kxm) |
342 |
xxm4(:,:,1:kxm) = el2(:,:,1:kxm)* shear(:,:,1:kxm) |
343 |
xxm5(:,:,1:kxm) = el(:,:,2:kx )* qm3(:,:,1:kxm) |
344 |
|
345 |
!------------------------------------------------------------------- |
346 |
! --- MOMENTUM |
347 |
!------------------------------------------------------------------- |
348 |
|
349 |
xxm1 = xxm5*( ckm1*qm2 + ckm2*xxm3 ) |
350 |
xxm2 = qm4 + ckm5*qm2*xxm4 + xxm3*( ckm6*xxm4 + ckm7*qm2 + ckm8*xxm3 ) |
351 |
|
352 |
xxm2 = MAX( xxm2, 0.2*qm4 ) |
353 |
xxm2 = MAX( xxm2, small ) |
354 |
|
355 |
akm(:,:,1) = 0.0 |
356 |
akm(:,:,2:kx) = xxm1(:,:,1:kxm) / xxm2(:,:,1:kxm) |
357 |
|
358 |
akm = MAX( akm, 0.0 ) |
359 |
|
360 |
!------------------------------------------------------------------- |
361 |
! --- HEAT AND MOISTURE |
362 |
!------------------------------------------------------------------- |
363 |
|
364 |
xxm1(:,:,1:kxm) = ckh1*xxm5(:,:,1:kxm) - ckh2*xxm4(:,:,1:kxm)*akm(:,:,2:kx) |
365 |
xxm2(:,:,1:kxm) = qm2(:,:,1:kxm) + ckh3*xxm3(:,:,1:kxm) |
366 |
|
367 |
xxm1 = MAX( xxm1, ckh4*xxm5 ) |
368 |
xxm2 = MAX( xxm2, 0.4*qm2 ) |
369 |
xxm2 = MAX( xxm2, small ) |
370 |
|
371 |
akh(:,:,1) = 0.0 |
372 |
akh(:,:,2:kx) = xxm1(:,:,1:kxm) / xxm2(:,:,1:kxm) |
373 |
|
374 |
!------------------------------------------------------------------- |
375 |
! --- BOUNDS |
376 |
!------------------------------------------------------------------- |
377 |
|
378 |
! --- UPPER BOUND |
379 |
akm = MIN( akm, akmax ) |
380 |
akh = MIN( akh, akmax ) |
381 |
|
382 |
! --- LOWER BOUND |
383 |
! where( akm(:,:,1:klim) < small ) akm(:,:,1:klim) = 0.0 |
384 |
! where( akh(:,:,1:klim) < small ) akh(:,:,1:klim) = 0.0 |
385 |
|
386 |
! --- LOWER BOUND NEAR SURFACE |
387 |
|
388 |
akmin = akmin_land*fracland + akmin_sea*(1.-fracland) |
389 |
|
390 |
if( PRESENT( kbot ) ) then |
391 |
do j = 1,jx |
392 |
do i = 1,ix |
393 |
klim = kbot(i,j) - nk_lim + 1 |
394 |
do k = klim,kbot(i,j) |
395 |
akm(i,j,k) = MAX( akm(i,j,k), akmin(i,j) ) |
396 |
akh(i,j,k) = MAX( akh(i,j,k), akmin(i,j) ) |
397 |
end do |
398 |
end do |
399 |
end do |
400 |
else |
401 |
klim = kx - nk_lim + 1 |
402 |
do k = klim,kx |
403 |
akm(:,:,k) = MAX( akm(:,:,k), akmin(:,:) ) |
404 |
akh(:,:,k) = MAX( akh(:,:,k), akmin(:,:) ) |
405 |
end do |
406 |
endif |
407 |
|
408 |
!------------------------------------------------------------------- |
409 |
! --- MASK OUT UNDERGROUND VALUES FOR ETA COORDINATE |
410 |
!------------------------------------------------------------------- |
411 |
|
412 |
if( PRESENT( mask ) ) then |
413 |
akm(:,:,1:kx) = akm(:,:,1:kx) * mask(:,:,1:kx) |
414 |
akh(:,:,1:kx) = akh(:,:,1:kx) * mask(:,:,1:kx) |
415 |
endif |
416 |
|
417 |
!==================================================================== |
418 |
! --- PROGNOSTICATE TURBULENT KE |
419 |
!==================================================================== |
420 |
|
421 |
cvfqdt = cvfq1 * delt |
422 |
dvfqdt = cvfq2 * delt * 2.0 |
423 |
|
424 |
!------------------------------------------------------------------- |
425 |
! --- PART OF LINEARIZED ENERGY DISIIPATION TERM |
426 |
!------------------------------------------------------------------- |
427 |
|
428 |
xxm1(:,:,1:kxm) = dvfqdt * qm(:,:,1:kxm) / el(:,:,2:kx) |
429 |
|
430 |
!------------------------------------------------------------------- |
431 |
! --- PART OF LINEARIZED VERTICAL DIFFUSION TERM |
432 |
!------------------------------------------------------------------- |
433 |
|
434 |
xx1(:,:,1:kx) = el(:,:,2:kxp) * qm(:,:,1:kx) |
435 |
|
436 |
xx2(:,:,1) = 0.5* xx1(:,:,1) |
437 |
xx2(:,:,2:kx) = 0.5*( xx1(:,:,2:kx) + xx1(:,:,1:kxm) ) |
438 |
|
439 |
xx1 = xx2 * dsdz |
440 |
|
441 |
!------------------------------------------------------------------- |
442 |
! --- IMPLICIT TIME DIFFERENCING FOR VERTICAL DIFFUSION |
443 |
! --- AND ENERGY DISSIPATION TERM |
444 |
!------------------------------------------------------------------- |
445 |
|
446 |
aaa(:,:,1:kxm) = -cvfqdt * xx1(:,:,2:kx ) * dsdzh(:,:,1:kxm) |
447 |
ccc(:,:,1:kxm) = -cvfqdt * xx1(:,:,1:kxm) * dsdzh(:,:,1:kxm) |
448 |
bbb(:,:,1:kxm) = 1.0 - aaa(:,:,1:kxm) - ccc(:,:,1:kxm) |
449 |
bbb(:,:,1:kxm) = bbb(:,:,1:kxm) + xxm1(:,:,1:kxm) |
450 |
ddd(:,:,1:kxm) = TKE(:,:,2:kx ) |
451 |
|
452 |
! correction for vertical diffusion of TKE surface boundary condition |
453 |
|
454 |
if (present(kbot)) then |
455 |
do j = 1,jx |
456 |
do i = 1,ix |
457 |
k = kbot(i,j) |
458 |
ddd(:,:,k-1) = ddd(:,:,k-1) - aaa(:,:,k-1) * TKE(:,:,k+1) |
459 |
enddo |
460 |
enddo |
461 |
else |
462 |
ddd(:,:,kxm) = ddd(:,:,kxm) - aaa(:,:,kxm) * TKE(:,:,kxp) |
463 |
endif |
464 |
|
465 |
! mask out terms below ground |
466 |
|
467 |
if (present(mask)) then |
468 |
where (mask(:,:,2:kx) < 0.1) ddd(:,:,1:kxm) = 0.0 |
469 |
endif |
470 |
|
471 |
|
472 |
!del CALL TRI_INVERT( xxm1, ddd, aaa, bbb, ccc ) |
473 |
!del CALL CLOSE_TRIDIAGONAL |
474 |
|
475 |
!------------------------------------------------------------------- |
476 |
! --- MASK OUT UNDERGROUND VALUES FOR ETA COORDINATE |
477 |
!------------------------------------------------------------------- |
478 |
|
479 |
if( PRESENT( mask ) ) then |
480 |
where (mask(:,:,2:kx) < 0.1) xxm1(:,:,1:kxm) = TKE(:,:,2:kx) |
481 |
endif |
482 |
|
483 |
!------------------------------------------------------------------- |
484 |
! --- SHEAR AND BUOYANCY TERMS |
485 |
!------------------------------------------------------------------- |
486 |
|
487 |
xxm2(:,:,1:kxm) = delt*( akm(:,:,2:kx)* shear(:,:,1:kxm) & |
488 |
- akh(:,:,2:kx)*buoync(:,:,1:kxm) ) |
489 |
|
490 |
!------------------------------------------------------------------- |
491 |
! --- UPDATE TURBULENT KINETIC ENERGY |
492 |
!------------------------------------------------------------------- |
493 |
|
494 |
TKE(:,:,1) = 0.0 |
495 |
TKE(:,:,2:kx) = xxm1(:,:,1:kxm) + xxm2(:,:,1:kxm) |
496 |
|
497 |
!==================================================================== |
498 |
! --- BOUND TURBULENT KINETIC ENERGY |
499 |
!==================================================================== |
500 |
|
501 |
TKE = MIN( TKE, TKEmax ) |
502 |
TKE = MAX( TKE, TKEmin ) |
503 |
|
504 |
if( PRESENT( mask ) ) then |
505 |
where (mask(:,:,1:kx) < 0.1) TKE(:,:,2:kxp) = 0.0 |
506 |
endif |
507 |
|
508 |
!==================================================================== |
509 |
! --- COMPUTE PBL DEPTH IF DESIRED |
510 |
!==================================================================== |
511 |
|
512 |
if (present(h)) then |
513 |
|
514 |
if (.not.present(ustar).or..not.present(bstar)) then |
515 |
!del CALL ERROR_MESG( ' MY25_TURB', & |
516 |
!del 'cannot request pbl depth diagnostic if ustar'// & |
517 |
!del ' and bstar are not also supplied', FATAL ) |
518 |
end if |
519 |
if (present(kbot)) then |
520 |
call k_pbl_depth(ustar,bstar,akm,akh,zsfc,zfull,zhalf,& |
521 |
h,kbot=kbot) |
522 |
else |
523 |
call k_pbl_depth(ustar,bstar,akm,akh,zsfc,zfull,zhalf,h) |
524 |
end if |
525 |
|
526 |
end if |
527 |
!==================================================================== |
528 |
|
529 |
! $$$$$$$$$$$$$$$$$ |
530 |
end do |
531 |
! $$$$$$$$$$$$$$$$$ |
532 |
|
533 |
!==================================================================== |
534 |
end SUBROUTINE MY25_TURB |
535 |
|
536 |
!####################################################################### |
537 |
|
538 |
SUBROUTINE MY25_TURB_INIT( ix, jx, kx ) |
539 |
|
540 |
!======================================================================= |
541 |
! ***** INITIALIZE MELLOR-YAMADA |
542 |
!======================================================================= |
543 |
!--------------------------------------------------------------------- |
544 |
! Arguments (Intent in) |
545 |
! ix, jx - Horizontal dimensions for global storage arrays |
546 |
! kx - Number of vertical levels in model |
547 |
!--------------------------------------------------------------------- |
548 |
integer, intent(in) :: ix, jx, kx |
549 |
!--------------------------------------------------------------------- |
550 |
! (Intent local) |
551 |
!--------------------------------------------------------------------- |
552 |
integer :: unit, io, ierr |
553 |
real :: actp, facm |
554 |
real, dimension(15) :: au, tem |
555 |
|
556 |
!===================================================================== |
557 |
|
558 |
!--------------------------------------------------------------------- |
559 |
! --- Read namelist |
560 |
!--------------------------------------------------------------------- |
561 |
|
562 |
!del if( FILE_EXIST( 'input.nml' ) ) then |
563 |
! ------------------------------------- |
564 |
!del unit = OPEN_NAMELIST_FILE () |
565 |
!del ierr = 1 |
566 |
!del do while( ierr .ne. 0 ) |
567 |
!del READ ( unit, nml = my25_turb_nml, iostat = io, end = 10 ) |
568 |
!del ierr = check_nml_error (io, 'my25_turb_nml') |
569 |
!del end do |
570 |
10 continue |
571 |
!del CALL CLOSE_FILE( unit ) |
572 |
! ------------------------------------- |
573 |
!del end if |
574 |
|
575 |
!--------------------------------------------------------------------- |
576 |
! --- Output version |
577 |
!--------------------------------------------------------------------- |
578 |
|
579 |
|
580 |
!del if ( mpp_pe() == mpp_root_pe() )write (stdlog(), nml=my25_turb_nml) |
581 |
!del call write_version_number(version, tag) |
582 |
|
583 |
!--------------------------------------------------------------------- |
584 |
! --- Initialize constants |
585 |
!--------------------------------------------------------------------- |
586 |
|
587 |
if( use_old_cons ) then |
588 |
aa1 = aa1_old |
589 |
aa2 = aa2_old |
590 |
bb1 = bb1_old |
591 |
bb2 = bb2_old |
592 |
ccc = ccc_old |
593 |
else |
594 |
aa1 = aa1_new |
595 |
aa2 = aa2_new |
596 |
bb1 = bb1_new |
597 |
bb2 = bb2_new |
598 |
ccc = ccc_new |
599 |
end if |
600 |
|
601 |
ckm1 = ( 1.0 - 3.0*ccc )*aa1 |
602 |
ckm3 = 3.0 * aa1*aa2* ( bb2 - 3.0*aa2 ) |
603 |
ckm4 = 9.0 * aa1*aa2*ccc*( bb2 + 4.0*aa1 ) |
604 |
ckm5 = 6.0 * aa1*aa1 |
605 |
ckm6 = 18.0 * aa1*aa1*aa2*( bb2 - 3.0*aa2 ) |
606 |
ckm7 = 3.0 * aa2* ( bb2 + 7.0*aa1 ) |
607 |
ckm8 = 27.0 * aa1*aa2*aa2*( bb2 + 4.0*aa1 ) |
608 |
ckm2 = ckm3 - ckm4 |
609 |
ckh1 = aa2 |
610 |
ckh2 = 6.0 * aa1*aa2 |
611 |
ckh3 = 3.0 * aa2*( bb2 + 4.0*aa1 ) |
612 |
ckh4 = 2.0e-6 * aa2 |
613 |
cvfq1 = 5.0 * cc1 / 3.0 |
614 |
cvfq2 = 1.0 / bb1 |
615 |
bcq = 0.5 * ( bb1**(2.0/3.0) ) |
616 |
|
617 |
!--------------------------------------------------------------------- |
618 |
! --- Allocate storage for TKE |
619 |
!--------------------------------------------------------------------- |
620 |
|
621 |
if( ALLOCATED( TKE ) ) DEALLOCATE( TKE ) |
622 |
ALLOCATE( TKE(ix,jx,kx+1) ) |
623 |
|
624 |
!------------------------------------------------------------------- |
625 |
do_init = .false. |
626 |
!--------------------------------------------------------------------- |
627 |
! --- Input TKE |
628 |
!--------------------------------------------------------------------- |
629 |
|
630 |
!del if( FILE_EXIST( 'INPUT/my25_turb.res' ) ) then |
631 |
|
632 |
!del unit = OPEN_FILE ( file = 'INPUT/my25_turb.res', & |
633 |
!del form = 'native', action = 'read' ) |
634 |
!del call read_data ( unit, TKE ) |
635 |
!del CALL CLOSE_FILE( unit ) |
636 |
|
637 |
init_tke = .false. |
638 |
|
639 |
!del else |
640 |
|
641 |
TKE = TKEmin |
642 |
|
643 |
init_tke = .true. |
644 |
num_total_pts = ix*jx |
645 |
pts_done = 0 |
646 |
|
647 |
!del endif |
648 |
|
649 |
|
650 |
!===================================================================== |
651 |
end SUBROUTINE MY25_TURB_INIT |
652 |
|
653 |
!####################################################################### |
654 |
|
655 |
SUBROUTINE MY25_TURB_END |
656 |
!======================================================================= |
657 |
integer :: unit |
658 |
!======================================================================= |
659 |
|
660 |
!del unit = OPEN_FILE ( file = 'RESTART/my25_turb.res', & |
661 |
!del form = 'native', action = 'write' ) |
662 |
!del call write_data ( unit, TKE ) |
663 |
!del CALL CLOSE_FILE ( unit ) |
664 |
|
665 |
!===================================================================== |
666 |
|
667 |
end SUBROUTINE MY25_TURB_END |
668 |
|
669 |
!####################################################################### |
670 |
|
671 |
SUBROUTINE K_PBL_DEPTH(ustar,bstar,akm,akh,zsfc,zfull,zhalf,h,kbot) |
672 |
!======================================================================= |
673 |
|
674 |
real, intent(in), dimension(:,:) :: ustar, bstar, zsfc |
675 |
real, intent(in), dimension(:,:,:) :: zhalf, zfull |
676 |
real, intent(in), dimension(:,:,:) :: akm, akh |
677 |
real, intent(out),dimension(:,:) :: h |
678 |
integer, intent(in), optional, dimension(:,:) :: kbot |
679 |
|
680 |
!======================================================================= |
681 |
real, dimension(size(zfull,1),size(zfull,2)) :: km_surf, kh_surf |
682 |
real, dimension(size(zfull,1),size(zfull,2)) :: zhalfhalf |
683 |
real, dimension(size(zfull,1),size(zfull,2),size(zfull,3)) :: zfull_ag |
684 |
real, dimension(size(zfull,1),size(zfull,2),size(zfull,3)+1):: zhalf_ag |
685 |
real, dimension(size(zfull,1),size(zfull,2),size(zfull,3)+1):: diff_tm |
686 |
integer, dimension(size(zfull,1),size(zfull,2)) :: ibot |
687 |
integer :: i,j,k,nlev,nlat,nlon |
688 |
|
689 |
|
690 |
nlev = size(zfull,3) |
691 |
nlat = size(zfull,2) |
692 |
nlon = size(zfull,1) |
693 |
|
694 |
!compute height of surface |
695 |
if (present(kbot)) then |
696 |
ibot=kbot |
697 |
else |
698 |
ibot(:,:) = nlev |
699 |
end if |
700 |
|
701 |
!compute density profile, and heights relative to surface |
702 |
do k = 1, nlev |
703 |
zfull_ag(:,:,k) = zfull(:,:,k) - zsfc(:,:) |
704 |
zhalf_ag(:,:,k) = zhalf(:,:,k) - zsfc(:,:) |
705 |
end do |
706 |
zhalf_ag(:,:,nlev+1) = zhalf(:,:,nlev+1) - zsfc(:,:) |
707 |
|
708 |
|
709 |
!compute height half way between surface and lowest model level |
710 |
zhalfhalf=0.5*MINVAL(MAX(zfull_ag,0.0)) |
711 |
|
712 |
!compute k's there by a call to mo_diff |
713 |
call mo_diff(zhalfhalf,ustar,bstar,km_surf,kh_surf,1) |
714 |
|
715 |
!create combined surface k's and diffusivity matrix |
716 |
diff_tm(:,:,nlev+1) = 0. |
717 |
diff_tm(:,:,1:nlev) = 0.5*(akm+akh) |
718 |
if (present(kbot)) then |
719 |
do j=1,nlat |
720 |
do i=1,nlon |
721 |
diff_tm(i,j,ibot(i,j)+1) = 0.5*(km_surf(i,j)+kh_surf(i,j)) |
722 |
zhalf_ag(i,j,ibot(i,j)+1) = zhalfhalf(i,j) |
723 |
enddo |
724 |
enddo |
725 |
else |
726 |
diff_tm(:,:,nlev+1) = 0.5*(km_surf(:,:)+kh_surf(:,:)) |
727 |
zhalf_ag(:,:,nlev+1) = zhalfhalf(:,:) |
728 |
end if |
729 |
|
730 |
|
731 |
!determine pbl depth as the height above ground where diff_tm |
732 |
!first falls beneath a critical value kcrit. If the value between |
733 |
!ground and level 1 does not exceed kcrit set pbl depth equal to zero. |
734 |
!kcrit is a namelist parameter. |
735 |
|
736 |
do j = 1,nlat |
737 |
do i = 1,nlon |
738 |
if (diff_tm(i,j,ibot(i,j)+1) .gt. kcrit) then |
739 |
k=ibot(i,j)+1 |
740 |
do while (k.gt. 2 .and. & |
741 |
diff_tm(i,j,k-1).gt.kcrit) |
742 |
k=k-1 |
743 |
enddo |
744 |
h(i,j) = zhalf_ag(i,j,k) + & |
745 |
(zhalf_ag(i,j,k-1)-zhalf_ag(i,j,k))* & |
746 |
(diff_tm(i,j,k)-kcrit) / & |
747 |
(diff_tm(i,j,k)-diff_tm(i,j,k-1)) |
748 |
else |
749 |
h(i,j) = 0. |
750 |
end if |
751 |
enddo |
752 |
enddo |
753 |
|
754 |
!----------------------------------------------------------------------- |
755 |
|
756 |
|
757 |
!===================================================================== |
758 |
end SUBROUTINE K_PBL_DEPTH |
759 |
|
760 |
!####################################################################### |
761 |
|
762 |
SUBROUTINE TKE_SURF ( u_star, TKE, kbot ) |
763 |
|
764 |
!======================================================================= |
765 |
!--------------------------------------------------------------------- |
766 |
! Arguments (Intent in) |
767 |
! u_star - surface friction velocity (m/s) |
768 |
! kbot - OPTIONAL;lowest model level index (integer); |
769 |
! at levels > Kbot, Mask = 0. |
770 |
!--------------------------------------------------------------------- |
771 |
real, intent(in), dimension(:,:) :: u_star |
772 |
|
773 |
integer, intent(in), OPTIONAL, dimension(:,:) :: kbot |
774 |
|
775 |
!--------------------------------------------------------------------- |
776 |
! Arguments (Intent inout) |
777 |
! TKE - turbulent kinetic energy |
778 |
!--------------------------------------------------------------------- |
779 |
real, intent(inout), dimension(:,:,:) :: TKE |
780 |
|
781 |
!--------------------------------------------------------------------- |
782 |
! (Intent local) |
783 |
!--------------------------------------------------------------------- |
784 |
real, dimension(SIZE(u_star,1),SIZE(u_star,2)) :: x1 |
785 |
integer :: ix, jx, kxp, i, j, k |
786 |
|
787 |
!======================================================================= |
788 |
|
789 |
ix = SIZE( u_star, 1 ) |
790 |
jx = SIZE( u_star, 2 ) |
791 |
kxp = SIZE( TKE, 3 ) |
792 |
|
793 |
!--------------------------------------------------------------------- |
794 |
|
795 |
x1 = bcq * u_star * u_star |
796 |
|
797 |
if( PRESENT( kbot ) ) then |
798 |
do j = 1,jx |
799 |
do i = 1,ix |
800 |
k = kbot(i,j) + 1 |
801 |
TKE(i,j,k) = x1(i,j) |
802 |
end do |
803 |
end do |
804 |
else |
805 |
TKE(:,:,kxp) = x1(:,:) |
806 |
endif |
807 |
|
808 |
!======================================================================= |
809 |
end SUBROUTINE TKE_SURF |
810 |
|
811 |
!####################################################################### |
812 |
end MODULE MY25_TURB_MOD |