Line data Source code
1 : module zm_conv_momtran
2 :
3 : use ccpp_kinds, only: kind_phys
4 :
5 : implicit none
6 :
7 : save
8 : private ! Make default type private to the module
9 : public zm_conv_momtran_run ! convective momentum transport
10 : integer, parameter, private :: num_winds=2 ! Number of wind directions (for historical purposes)
11 :
12 :
13 : contains
14 :
15 : !===============================================================================
16 : !> \section arg_table_zm_conv_momtran_run Argument Table
17 : !! \htmlinclude zm_conv_momtran_run.html
18 : !!
19 99072 : subroutine zm_conv_momtran_run(ncol, pver, pverp, &
20 99072 : domomtran,windu, windv, mu, md, &
21 : momcu, momcd, &
22 99072 : du, eu, ed, dp, dsubcld , &
23 99072 : jt, mx, ideep , il1g, il2g, &
24 99072 : nstep, windu_tend, windv_tend, pguallu, pguallv, pgdallu, pgdallv, &
25 297216 : icwuu, icwuv, icwdu, icwdv, dt, seten, scheme_name, errmsg, errflg)
26 : !-----------------------------------------------------------------------
27 : !
28 : ! Purpose:
29 : ! Convective transport of momentum
30 : !
31 : ! Mixing ratios may be with respect to either dry or moist air
32 : !
33 : ! Method:
34 : ! Based on the convtran subroutine by P. Rasch
35 : ! <Also include any applicable external references.>
36 : !
37 : ! Author: J. Richter and P. Rasch
38 : !
39 : !-----------------------------------------------------------------------
40 :
41 : implicit none
42 : !-----------------------------------------------------------------------
43 : !
44 : ! Input arguments
45 : !
46 : integer, intent(in) :: ncol ! number of atmospheric columns
47 : integer, intent(in) :: pver, pverp
48 : logical, intent(in) :: domomtran ! flag for doing convective transport
49 : real(kind_phys), intent(in) :: windu(:,:) ! U Wind array (ncol,pver)
50 : real(kind_phys), intent(in) :: windv(:,:) ! V Wind array (ncol,pver)
51 : real(kind_phys), intent(in) :: mu(:,:) ! Mass flux up (ncol,pver)
52 : real(kind_phys), intent(in) :: md(:,:) ! Mass flux down (ncol,pver)
53 : real(kind_phys), intent(in) :: momcu
54 : real(kind_phys), intent(in) :: momcd
55 : real(kind_phys), intent(in) :: du(:,:) ! Mass detraining from updraft (ncol,pver)
56 : real(kind_phys), intent(in) :: eu(:,:) ! Mass entraining from updraft (ncol,pver)
57 : real(kind_phys), intent(in) :: ed(:,:) ! Mass entraining from downdraft (ncol,pver)
58 : real(kind_phys), intent(in) :: dp(:,:) ! Delta pressure between interfaces (ncol,pver)
59 : real(kind_phys), intent(in) :: dsubcld(:) ! Delta pressure from cloud base to sfc (ncol)
60 : real(kind_phys), intent(in) :: dt ! time step in seconds
61 :
62 : integer, intent(in) :: jt(:) ! Index of cloud top for each column (ncol)
63 : integer, intent(in) :: mx(:) ! Index of cloud top for each column (ncol)
64 : integer, intent(in) :: ideep(:) ! Gathering array (ncol)
65 : integer, intent(in) :: il1g ! Gathered min lon indices over which to operate
66 : integer, intent(in) :: il2g ! Gathered max lon indices over which to operate
67 : integer, intent(in) :: nstep ! Time step index
68 :
69 :
70 :
71 : ! input/output
72 :
73 : real(kind_phys), intent(out) :: windu_tend(:,:) ! U wind tendency
74 : real(kind_phys), intent(out) :: windv_tend(:,:) ! V wind tendency
75 :
76 : character(len=512), intent(out) :: errmsg
77 : integer, intent(out) :: errflg
78 : character(len=40), intent(out) :: scheme_name
79 :
80 : !--------------------------Local Variables------------------------------
81 :
82 : integer i ! Work index
83 : integer k ! Work index
84 : integer kbm ! Highest altitude index of cloud base
85 : integer kk ! Work index
86 : integer kkp1 ! Work index
87 : integer kkm1 ! Work index
88 : integer km1 ! Work index
89 : integer kp1 ! Work index
90 : integer ktm ! Highest altitude index of cloud top
91 : integer m ! Work index
92 : integer ii ! Work index
93 :
94 : real(kind_phys) cabv ! Mix ratio of constituent above
95 : real(kind_phys) cbel ! Mix ratio of constituent below
96 : real(kind_phys) cdifr ! Normalized diff between cabv and cbel
97 198144 : real(kind_phys) chat(ncol,pver) ! Mix ratio in env at interfaces
98 198144 : real(kind_phys) cond(ncol,pver) ! Mix ratio in downdraft at interfaces
99 198144 : real(kind_phys) const(ncol,pver) ! Gathered wind array
100 198144 : real(kind_phys) conu(ncol,pver) ! Mix ratio in updraft at interfaces
101 198144 : real(kind_phys) dcondt(ncol,pver) ! Gathered tend array
102 : real(kind_phys) mbsth ! Threshold for mass fluxes
103 : real(kind_phys) mupdudp ! A work variable
104 : real(kind_phys) minc ! A work variable
105 : real(kind_phys) maxc ! A work variable
106 : real(kind_phys) fluxin ! A work variable
107 : real(kind_phys) fluxout ! A work variable
108 : real(kind_phys) netflux ! A work variable
109 :
110 :
111 : real(kind_phys) sum ! sum
112 : real(kind_phys) sum2 ! sum2
113 :
114 198144 : real(kind_phys) mududp(ncol,pver) ! working variable
115 198144 : real(kind_phys) mddudp(ncol,pver) ! working variable
116 :
117 198144 : real(kind_phys) pgu(ncol,pver) ! Pressure gradient term for updraft
118 198144 : real(kind_phys) pgd(ncol,pver) ! Pressure gradient term for downdraft
119 :
120 : real(kind_phys),intent(out) :: pguallu(:,:) ! Apparent force from updraft PG on U winds ! (ncol,pver)
121 : real(kind_phys),intent(out) :: pguallv(:,:) ! Apparent force from updraft PG on V winds ! (ncol,pver)
122 : real(kind_phys),intent(out) :: pgdallu(:,:) ! Apparent force from downdraft PG on U winds! (ncol,pver)
123 : real(kind_phys),intent(out) :: pgdallv(:,:) ! Apparent force from downdraft PG on V winds! (ncol,pver)
124 :
125 : real(kind_phys),intent(out) :: icwuu(:,:) ! In-cloud U winds in updraft ! (ncol,pver)
126 : real(kind_phys),intent(out) :: icwuv(:,:) ! In-cloud V winds in updraft ! (ncol,pver)
127 : real(kind_phys),intent(out) :: icwdu(:,:) ! In-cloud U winds in downdraft ! (ncol,pver)
128 : real(kind_phys),intent(out) :: icwdv(:,:) ! In-cloud V winds in downdraft ! (ncol,pver)
129 :
130 : real(kind_phys),intent(out) :: seten(:,:) ! Dry static energy tendency ! (ncol,pver)
131 198144 : real(kind_phys) gseten(ncol,pver) ! Gathered dry static energy tendency
132 :
133 198144 : real(kind_phys) :: winds(ncol,pver,num_winds) ! combined winds array
134 198144 : real(kind_phys) :: wind_tends(ncol,pver,num_winds) ! combined tendency array
135 198144 : real(kind_phys) :: pguall(ncol,pver,num_winds) ! Combined apparent force from updraft PG on U winds
136 198144 : real(kind_phys) :: pgdall(ncol,pver,num_winds) ! Combined apparent force from downdraft PG on U winds
137 198144 : real(kind_phys) :: icwu(ncol,pver,num_winds) ! Combined In-cloud winds in updraft
138 198144 : real(kind_phys) :: icwd(ncol,pver,num_winds) ! Combined In-cloud winds in downdraft
139 :
140 198144 : real(kind_phys) mflux(ncol,pverp,num_winds) ! Gathered momentum flux
141 :
142 198144 : real(kind_phys) wind0(ncol,pver,num_winds) ! gathered wind before time step
143 99072 : real(kind_phys) windf(ncol,pver,num_winds) ! gathered wind after time step
144 : real(kind_phys) fkeb, fket, ketend_cons, ketend, utop, ubot, vtop, vbot, gset2
145 :
146 :
147 : !-----------------------------------------------------------------------
148 99072 : scheme_name = "zm_conv_momtran_run"
149 99072 : errmsg = ''
150 99072 : errflg = 0
151 :
152 : ! Combine winds in single array
153 154045440 : winds(:,:,1) = windu(:,:)
154 153946368 : winds(:,:,2) = windv(:,:)
155 :
156 : ! Initialize outgoing fields
157 307991808 : pguall(:,:,:) = 0.0_kind_phys
158 307991808 : pgdall(:,:,:) = 0.0_kind_phys
159 : ! Initialize in-cloud winds to environmental wind
160 307991808 : icwu(:ncol,:,:) = winds(:ncol,:,:)
161 307991808 : icwd(:ncol,:,:) = winds(:ncol,:,:)
162 :
163 : ! Initialize momentum flux and final winds
164 311300352 : mflux(:,:,:) = 0.0_kind_phys
165 307991808 : wind0(:,:,:) = 0.0_kind_phys
166 307991808 : windf(:,:,:) = 0.0_kind_phys
167 :
168 : ! Initialize dry static energy
169 :
170 153946368 : seten(:,:) = 0.0_kind_phys
171 153946368 : gseten(:,:) = 0.0_kind_phys
172 :
173 : ! mbsth is the threshold below which we treat the mass fluxes as zero (in mb/s)
174 99072 : mbsth = 1.e-15_kind_phys
175 :
176 : ! Find the highest level top and bottom levels of convection
177 99072 : ktm = pver
178 99072 : kbm = pver
179 360519 : do i = il1g, il2g
180 261447 : ktm = min(ktm,jt(i))
181 360519 : kbm = min(kbm,mx(i))
182 : end do
183 :
184 : ! Loop ever each wind component
185 297216 : do m = 1, num_winds !start at m = 1 to transport momentum
186 297216 : if (domomtran) then
187 :
188 : ! Gather up the winds and set tend to zero
189 18625536 : do k = 1,pver
190 67254678 : do i =il1g,il2g
191 48629142 : const(i,k) = winds(ideep(i),k,m)
192 67056534 : wind0(i,k,m) = const(i,k)
193 : end do
194 : end do
195 :
196 :
197 : ! From now on work only with gathered data
198 :
199 : ! Interpolate winds to interfaces
200 :
201 18625536 : do k = 1,pver
202 18427392 : km1 = max(1,k-1)
203 67254678 : do i = il1g, il2g
204 :
205 : ! use arithmetic mean
206 48629142 : chat(i,k) = 0.5_kind_phys* (const(i,k)+const(i,km1))
207 :
208 : ! Provisional up and down draft values
209 48629142 : conu(i,k) = chat(i,k)
210 48629142 : cond(i,k) = chat(i,k)
211 :
212 : ! provisional tends
213 67056534 : dcondt(i,k) = 0._kind_phys
214 :
215 : end do
216 : end do
217 :
218 :
219 : !
220 : ! Pressure Perturbation Term
221 : !
222 :
223 : !Top boundary: assume mu is zero
224 :
225 198144 : k=1
226 919182 : pgu(:il2g,k) = 0.0_kind_phys
227 721038 : pgd(:il2g,k) = 0.0_kind_phys
228 :
229 18229248 : do k=2,pver-1
230 18031104 : km1 = max(1,k-1)
231 18031104 : kp1 = min(pver,k+1)
232 65812602 : do i = il1g,il2g
233 :
234 : !interior points
235 :
236 142750062 : mududp(i,k) = ( mu(i,k) * (const(i,k)- const(i,km1))/dp(i,km1) &
237 190333416 : + mu(i,kp1) * (const(i,kp1) - const(i,k))/dp(i,k))
238 :
239 47583354 : pgu(i,k) = - momcu * 0.5_kind_phys * mududp(i,k)
240 :
241 :
242 47583354 : mddudp(i,k) = ( md(i,k) * (const(i,k)- const(i,km1))/dp(i,km1) &
243 47583354 : + md(i,kp1) * (const(i,kp1) - const(i,k))/dp(i,k))
244 :
245 65614458 : pgd(i,k) = - momcd * 0.5_kind_phys * mddudp(i,k)
246 :
247 :
248 : end do
249 : end do
250 :
251 : ! bottom boundary
252 198144 : k = pver
253 198144 : km1 = max(1,k-1)
254 721038 : do i=il1g,il2g
255 :
256 522894 : mududp(i,k) = mu(i,k) * (const(i,k)- const(i,km1))/dp(i,km1)
257 522894 : pgu(i,k) = - momcu * mududp(i,k)
258 :
259 522894 : mddudp(i,k) = md(i,k) * (const(i,k)- const(i,km1))/dp(i,km1)
260 :
261 721038 : pgd(i,k) = - momcd * mddudp(i,k)
262 :
263 : end do
264 :
265 :
266 : !
267 : ! In-cloud velocity calculations
268 : !
269 :
270 : ! Do levels adjacent to top and bottom
271 721038 : k = 2
272 721038 : km1 = 1
273 721038 : kk = pver
274 721038 : kkm1 = max(1,kk-1)
275 721038 : do i = il1g,il2g
276 522894 : mupdudp = mu(i,kk) + du(i,kk)*dp(i,kk)
277 522894 : if (mupdudp > mbsth) then
278 :
279 0 : conu(i,kk) = (+eu(i,kk)*const(i,kk)*dp(i,kk)+pgu(i,kk)*dp(i,kk))/mupdudp
280 : endif
281 721038 : if (md(i,k) < -mbsth) then
282 0 : cond(i,k) = (-ed(i,km1)*const(i,km1)*dp(i,km1))-pgd(i,km1)*dp(i,km1)/md(i,k)
283 : endif
284 :
285 :
286 : end do
287 :
288 :
289 :
290 : ! Updraft from bottom to top
291 18427392 : do kk = pver-1,1,-1
292 18229248 : kkm1 = max(1,kk-1)
293 18229248 : kkp1 = min(pver,kk+1)
294 66533640 : do i = il1g,il2g
295 48106248 : mupdudp = mu(i,kk) + du(i,kk)*dp(i,kk)
296 66335496 : if (mupdudp > mbsth) then
297 :
298 17628312 : conu(i,kk) = ( mu(i,kkp1)*conu(i,kkp1)+eu(i,kk)* &
299 17628312 : const(i,kk)*dp(i,kk)+pgu(i,kk)*dp(i,kk))/mupdudp
300 : endif
301 : end do
302 :
303 : end do
304 :
305 :
306 : ! Downdraft from top to bottom
307 18229248 : do k = 3,pver
308 18031104 : km1 = max(1,k-1)
309 65812602 : do i = il1g,il2g
310 65614458 : if (md(i,k) < -mbsth) then
311 :
312 13810460 : cond(i,k) = ( md(i,km1)*cond(i,km1)-ed(i,km1)*const(i,km1) &
313 13810460 : *dp(i,km1)-pgd(i,km1)*dp(i,km1) )/md(i,k)
314 :
315 : endif
316 : end do
317 : end do
318 :
319 :
320 : sum = 0._kind_phys
321 : sum2 = 0._kind_phys
322 :
323 :
324 1822032 : do k = ktm,pver
325 1623888 : km1 = max(1,k-1)
326 1623888 : kp1 = min(pver,k+1)
327 14304946 : do i = il1g,il2g
328 12482914 : ii = ideep(i)
329 :
330 : ! version 1 hard to check for roundoff errors
331 12482914 : dcondt(i,k) = &
332 12482914 : +(mu(i,kp1)* (conu(i,kp1)-chat(i,kp1)) &
333 12482914 : -mu(i,k)* (conu(i,k)-chat(i,k)) &
334 12482914 : +md(i,kp1)* (cond(i,kp1)-chat(i,kp1)) &
335 12482914 : -md(i,k)* (cond(i,k)-chat(i,k)) &
336 39072630 : )/dp(i,k)
337 :
338 : end do
339 : end do
340 :
341 : ! dcont for bottom layer
342 : !
343 669230 : do k = kbm,pver
344 3134218 : km1 = max(1,k-1)
345 3332362 : do i = il1g,il2g
346 3134218 : if (k == mx(i)) then
347 :
348 : ! version 1
349 522894 : dcondt(i,k) = (1._kind_phys/dp(i,k))* &
350 522894 : (-mu(i,k)*(conu(i,k)-chat(i,k)) &
351 522894 : -md(i,k)*(cond(i,k)-chat(i,k)) &
352 1045788 : )
353 : end if
354 : end do
355 : end do
356 :
357 : ! Initialize to zero everywhere, then scatter tendency back to full array
358 307892736 : wind_tends(:,:,m) = 0._kind_phys
359 :
360 18625536 : do k = 1,pver
361 67254678 : do i = il1g,il2g
362 48629142 : ii = ideep(i)
363 48629142 : wind_tends(ii,k,m) = dcondt(i,k)
364 : ! Output apparent force on the mean flow from pressure gradient
365 48629142 : pguall(ii,k,m) = -pgu(i,k)
366 48629142 : pgdall(ii,k,m) = -pgd(i,k)
367 48629142 : icwu(ii,k,m) = conu(i,k)
368 67056534 : icwd(ii,k,m) = cond(i,k)
369 : end do
370 : end do
371 :
372 : ! Calculate momentum flux in units of mb*m/s2
373 :
374 1822032 : do k = ktm,pver
375 14304946 : do i = il1g,il2g
376 12482914 : ii = ideep(i)
377 12482914 : mflux(i,k,m) = &
378 12482914 : -mu(i,k)* (conu(i,k)-chat(i,k)) &
379 26589716 : -md(i,k)* (cond(i,k)-chat(i,k))
380 : end do
381 : end do
382 :
383 :
384 : ! Calculate winds at the end of the time step
385 :
386 1822032 : do k = ktm,pver
387 14304946 : do i = il1g,il2g
388 12482914 : ii = ideep(i)
389 12482914 : km1 = max(1,k-1)
390 12482914 : kp1 = k+1
391 14106802 : windf(i,k,m) = const(i,k) - (mflux(i,kp1,m) - mflux(i,k,m)) * dt /dp(i,k)
392 :
393 : end do
394 : end do
395 :
396 : end if ! for domomtran
397 : end do
398 :
399 : ! Need to add an energy fix to account for the dissipation of kinetic energy
400 : ! Formulation follows from Boville and Bretherton (2003)
401 : ! formulation by PJR
402 :
403 911016 : do k = ktm,pver
404 811944 : km1 = max(1,k-1)
405 811944 : kp1 = min(pver,k+1)
406 7152473 : do i = il1g,il2g
407 :
408 6241457 : ii = ideep(i)
409 :
410 : ! calculate the KE fluxes at top and bot of layer
411 : ! based on a discrete approximation to b&b eq(35) F_KE = u*F_u + v*F_v at interface
412 6241457 : utop = (wind0(i,k,1)+wind0(i,km1,1))/2._kind_phys
413 6241457 : vtop = (wind0(i,k,2)+wind0(i,km1,2))/2._kind_phys
414 6241457 : ubot = (wind0(i,kp1,1)+wind0(i,k,1))/2._kind_phys
415 6241457 : vbot = (wind0(i,kp1,2)+wind0(i,k,2))/2._kind_phys
416 6241457 : fket = utop*mflux(i,k,1) + vtop*mflux(i,k,2) ! top of layer
417 6241457 : fkeb = ubot*mflux(i,k+1,1) + vbot*mflux(i,k+1,2) ! bot of layer
418 :
419 : ! divergence of these fluxes should give a conservative redistribution of KE
420 6241457 : ketend_cons = (fket-fkeb)/dp(i,k)
421 :
422 : ! tendency in kinetic energy resulting from the momentum transport
423 6241457 : ketend = ((windf(i,k,1)**2 + windf(i,k,2)**2) - (wind0(i,k,1)**2 + wind0(i,k,2)**2))/dt
424 :
425 : ! the difference should be the dissipation
426 6241457 : gset2 = ketend_cons - ketend
427 7053401 : gseten(i,k) = gset2
428 :
429 : end do
430 :
431 : end do
432 :
433 : ! Scatter dry static energy to full array
434 9312768 : do k = 1,pver
435 33627339 : do i = il1g,il2g
436 24314571 : ii = ideep(i)
437 33528267 : seten(ii,k) = gseten(i,k)
438 :
439 : end do
440 : end do
441 :
442 : ! Split out the wind tendencies
443 154045440 : windu_tend(:,:) = wind_tends(:,:,1)
444 154045440 : windv_tend(:,:) = wind_tends(:,:,2)
445 :
446 154045440 : pguallu(:,:) = pguall(:,:,1)
447 154045440 : pguallv(:,:) = pguall(:,:,2)
448 154045440 : pgdallu(:,:) = pgdall(:,:,1)
449 154045440 : pgdallv(:,:) = pgdall(:,:,2)
450 154045440 : icwuu(:ncol,:) = icwu(:,:,1)
451 154045440 : icwuv(:ncol,:) = icwu(:,:,2)
452 154045440 : icwdu(:ncol,:) = icwd(:,:,1)
453 154045440 : icwdv(:ncol,:) = icwd(:,:,2)
454 :
455 99072 : return
456 : end subroutine zm_conv_momtran_run
457 :
458 :
459 : end module zm_conv_momtran
|