CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mynnedmf_wrapper.F90
1
4
9
10 contains
11
15 subroutine mynnedmf_wrapper_init ( &
16 & con_cp, con_grav, con_rd, con_rv, &
17 & con_cpv, con_cliq, con_cice, con_rcp, &
18 & con_XLV, con_XLF, con_p608, con_ep2, &
19 & con_karman, con_t0c, &
20 & do_mynnedmf, &
21 & errmsg, errflg )
22
23 use machine, only : kind_phys
25
26 implicit none
27
28 logical, intent(in) :: do_mynnedmf
29 character(len=*),intent(out):: errmsg
30 integer, intent(out) :: errflg
31
32 real(kind_phys),intent(in) :: con_xlv
33 real(kind_phys),intent(in) :: con_xlf
34 real(kind_phys),intent(in) :: con_rv
35 real(kind_phys),intent(in) :: con_rd
36 real(kind_phys),intent(in) :: con_ep2
37 real(kind_phys),intent(in) :: con_grav
38 real(kind_phys),intent(in) :: con_cp
39 real(kind_phys),intent(in) :: con_cpv
40 real(kind_phys),intent(in) :: con_rcp
41 real(kind_phys),intent(in) :: con_p608
42 real(kind_phys),intent(in) :: con_cliq
43 real(kind_phys),intent(in) :: con_cice
44 real(kind_phys),intent(in) :: con_karman
45 real(kind_phys),intent(in) :: con_t0c
46
47 ! Initialize CCPP error handling variables
48 errmsg = ''
49 errflg = 0
50
51 xlv = con_xlv
52 xlf = con_xlf
53 r_v = con_rv
54 r_d = con_rd
55 ep_2 = con_ep2
56 grav = con_grav
57 cp = con_cp
58 cpv = con_cpv
59 rcp = con_rcp
60 p608 = con_p608
61 cliq = con_cliq
62 cice = con_cice
63 karman = con_karman
64 t0c = con_t0c
65
66 xls = xlv+xlf != 2.85E6 (J/kg) sublimation
67 rvovrd = r_v/r_d != 1.608
68 ep_3 = 1.-ep_2 != 0.378
69 gtr = grav/tref
70 rk = cp/r_d
71 tv0 = p608*tref
72 tv1 = (1.+p608)*tref
73 xlscp = (xlv+xlf)/cp
74 xlvcp = xlv/cp
75 g_inv = 1./grav
76
77 ! Consistency checks
78 if (.not. do_mynnedmf) then
79 errmsg = 'Logic error: do_mynnedmf = .false.'
80 errflg = 1
81 return
82 end if
83
84 end subroutine mynnedmf_wrapper_init
85
91SUBROUTINE mynnedmf_wrapper_run( &
92 & im,levs, &
93 & flag_init,flag_restart, &
94 & lssav, ldiag3d, qdiag3d, &
95 & lsidea, cplflx, &
96 & delt,dtf,dx,zorl, &
97 & phii,u,v,omega,t3d, &
98 & qgrs_water_vapor, &
99 & qgrs_liquid_cloud, &
100 & qgrs_ice, &
101 & qgrs_snow, &
102 & qgrs_cloud_droplet_num_conc, &
103 & qgrs_cloud_ice_num_conc, &
104 & qgrs_ozone, &
105 & qgrs_water_aer_num_conc, &
106 & qgrs_ice_aer_num_conc, &
107 & qgrs_cccn, &
108 & prsl,prsi,exner, &
109 & slmsk,tsurf,qsfc,ps, &
110 & ust,ch,hflx,qflx,wspd,rb, &
111 & dtsfc1,dqsfc1, &
112 & dusfc1,dvsfc1, &
113 & dusfci_diag,dvsfci_diag, &
114 & dtsfci_diag,dqsfci_diag, &
115 & dusfc_diag,dvsfc_diag, &
116 & dtsfc_diag,dqsfc_diag, &
117 & dusfc_cice,dvsfc_cice, &
118 & dtsfc_cice,dqsfc_cice, &
119 & hflx_wat,qflx_wat,stress_wat, &
120 & oceanfrac,fice,wet,icy,dry, &
121 & dusfci_cpl,dvsfci_cpl, &
122 & dtsfci_cpl,dqsfci_cpl, &
123 & dusfc_cpl,dvsfc_cpl, &
124 & dtsfc_cpl,dqsfc_cpl, &
125 & recmol, &
126 & qke,qke_adv,Tsq,Qsq,Cov, &
127 & el_pbl,sh3d,sm3d,exch_h,exch_m, &
128 & dqke,qwt,qshear,qbuoy,qdiss, &
129 & Pblh,kpbl, &
130 & qc_bl,qi_bl,cldfra_bl, &
131 & edmf_a,edmf_w,edmf_qt, &
132 & edmf_thl,edmf_ent,edmf_qc, &
133 & sub_thl,sub_sqv,det_thl,det_sqv,&
134 & maxwidth,maxMF,ztop_plume, &
135 & ktop_plume, &
136 & dudt, dvdt, dtdt, &
137 & dqdt_water_vapor, dqdt_liquid_cloud, & ! <=== ntqv, ntcw
138 & dqdt_ice, dqdt_snow, & ! <=== ntiw, ntsw
139 & dqdt_ozone, & ! <=== ntoz
140 & dqdt_cloud_droplet_num_conc, dqdt_ice_num_conc, & ! <=== ntlnc, ntinc
141 & dqdt_water_aer_num_conc, dqdt_ice_aer_num_conc,& ! <=== ntwa, ntia
142 & dqdt_cccn, & ! <=== ntccn
143 & tmf, flag_for_pbl_generic_tend, &
144 & dtend, dtidx, index_of_temperature, &
145 & index_of_x_wind, index_of_y_wind, ntke, &
146 & ntqv, ntcw, ntiw, ntsw, &
147 & ntoz, ntlnc, ntinc, ntwa, ntia, &
148 & index_of_process_pbl, htrsw, htrlw, xmu, &
149 & tke_budget, bl_mynn_tkeadvect, &
150 & bl_mynn_cloudpdf, bl_mynn_mixlength, &
151 & bl_mynn_edmf, &
152 & bl_mynn_edmf_mom, bl_mynn_edmf_tke, &
153 & bl_mynn_cloudmix, bl_mynn_mixqt, &
154 & bl_mynn_output, bl_mynn_closure, &
155 & icloud_bl, do_mynnsfclay, &
156 & imp_physics, imp_physics_gfdl, &
157 & imp_physics_thompson, imp_physics_wsm6, &
158 & imp_physics_fa, imfdeepcnv, imfdeepcnv_c3, &
159 & imfdeepcnv_samf, &
160 & chem3d, frp, mix_chem, rrfs_sd, enh_mix, &
161 & nchem, ndvel, vdep, smoke_dbg, &
162 & imp_physics_nssl, nssl_ccn_on, &
163 & ltaerosol, mraerosol, spp_wts_pbl, spp_pbl, &
164 & lprnt, huge, errmsg, errflg )
165
166! should be moved to inside the mynn:
167 use machine, only: kind_phys
168 use bl_mynn_common, only: cp, r_d, grav, g_inv, zero, &
169 xlv, xlvcp, xlscp, p608
171
172!-------------------------------------------------------------------
173 implicit none
174!-------------------------------------------------------------------
175
176 real(kind_phys), intent(in) :: huge
177 character(len=*), intent(out) :: errmsg
178 integer, intent(out) :: errflg
179
180 logical, intent(in) :: lssav, ldiag3d, lsidea, qdiag3d
181 logical, intent(in) :: cplflx
182
183 !smoke/chem
184 integer, intent(in) :: nchem, ndvel
185 integer, parameter :: kdvel=1
186 logical, intent(in) :: smoke_dbg
187
188! NAMELIST OPTIONS (INPUT):
189 logical, intent(in) :: &
190 & bl_mynn_tkeadvect, &
191 & ltaerosol, &
192 & mraerosol, &
193 & lprnt, &
194 & do_mynnsfclay, &
195 & flag_for_pbl_generic_tend, &
196 & nssl_ccn_on
197 integer, intent(in) :: &
198 & bl_mynn_cloudpdf, &
199 & bl_mynn_mixlength, &
200 & icloud_bl, &
201 & bl_mynn_edmf, &
202 & bl_mynn_edmf_mom, &
203 & bl_mynn_edmf_tke, &
204 & bl_mynn_cloudmix, &
205 & bl_mynn_mixqt, &
206 & bl_mynn_output, &
207 & imp_physics, imp_physics_wsm6, &
208 & imp_physics_thompson, imp_physics_gfdl, &
209 & imp_physics_nssl, imp_physics_fa, imfdeepcnv, &
210 & imfdeepcnv_c3, imfdeepcnv_samf, &
211 & spp_pbl, &
212 & tke_budget
213 real(kind_phys), intent(in) :: &
214 & bl_mynn_closure
215
216!TENDENCY DIAGNOSTICS
217 real(kind_phys), intent(inout), optional :: dtend(:,:,:)
218 integer, intent(in) :: dtidx(:,:)
219 integer, intent(in) :: index_of_temperature, index_of_x_wind
220 integer, intent(in) :: index_of_y_wind, index_of_process_pbl
221 integer, intent(in) :: ntoz, ntqv, ntcw, ntiw, ntsw, ntlnc
222 integer, intent(in) :: ntinc, ntwa, ntia, ntke
223
224!MISC CONFIGURATION OPTIONS
225 INTEGER, PARAMETER :: &
226 & bl_mynn_mixscalars=1
227 LOGICAL :: &
228 & FLAG_QI, FLAG_QNI, FLAG_QC, FLAG_QS, FLAG_QNC, &
229 & FLAG_QNWFA, FLAG_QNIFA, FLAG_QNBCA, FLAG_OZONE
230 ! Define locally until needed from CCPP
231 LOGICAL, PARAMETER :: cycling = .false.
232
233!MYNN-1D
234 REAL(kind_phys), intent(in) :: delt, dtf
235 INTEGER, intent(in) :: im, levs
236 LOGICAL, intent(in) :: flag_init, flag_restart
237 INTEGER :: initflag, k, i
238 INTEGER :: IDS,IDE,JDS,JDE,KDS,KDE, &
239 & IMS,IME,JMS,JME,KMS,KME, &
240 & ITS,ITE,JTS,JTE,KTS,KTE
241
242 REAL(kind_phys) :: tem
243
244!MYNN-3D
245 real(kind_phys), dimension(:,:), intent(in) :: phii
246 real(kind_phys), dimension(:,:), intent(inout) :: &
247 & dtdt, dudt, dvdt, &
248 & dqdt_water_vapor, dqdt_liquid_cloud, dqdt_ice, &
249 & dqdt_snow, dqdt_ice_num_conc, dqdt_ozone
250 real(kind_phys), dimension(:,:), intent(inout), optional :: &
251 & dqdt_cloud_droplet_num_conc, dqdt_water_aer_num_conc, &
252 & dqdt_ice_aer_num_conc
253 real(kind_phys), dimension(:,:), intent(inout) :: qke, &
254 & EL_PBL, Sh3D, Sm3D, qc_bl, qi_bl, cldfra_bl
255 real(kind_phys), dimension(:,:), intent(inout), optional :: &
256 & dqdt_cccn
257 real(kind_phys), dimension(:,:), intent(inout) :: &
258 & qke_adv
259 real(kind_phys), dimension(:,:,:), intent(out) :: tmf
260 !These 10 arrays are only allocated when bl_mynn_output > 0
261 real(kind_phys), dimension(:,:), intent(inout), optional :: &
262 & edmf_a,edmf_w,edmf_qt, &
263 & edmf_thl,edmf_ent,edmf_qc, &
264 & sub_thl,sub_sqv,det_thl,det_sqv
265 real(kind_phys), dimension(:,:), intent(inout) :: &
266 & t3d,qgrs_water_vapor,qgrs_liquid_cloud,qgrs_ice, &
267 & qgrs_snow
268 real(kind_phys), dimension(:,:), intent(in) :: &
269 & qgrs_cloud_ice_num_conc, &
270 & u,v,omega, &
271 & exner,prsl,prsi, &
272 & qgrs_ozone
273 real(kind_phys), dimension(:,:), intent(in), optional :: &
274 & qgrs_water_aer_num_conc, &
275 & qgrs_cloud_droplet_num_conc, &
276 & qgrs_ice_aer_num_conc
277 real(kind_phys), dimension(:,:), intent(in), optional :: qgrs_cccn
278 real(kind_phys), dimension(:,:), intent(out) :: &
279 & Tsq, Qsq, Cov, exch_h, exch_m
280 real(kind_phys), dimension(:,:), intent(out), optional :: &
281 & dqke, qWT, qSHEAR, qBUOY, qDISS
282 real(kind_phys), dimension(:), intent(in) :: xmu
283 real(kind_phys), dimension(:,:), intent(in) :: htrsw, htrlw
284 ! spp_wts_pbl only allocated if spp_pbl == 1
285 real(kind_phys), dimension(:,:), intent(in), optional :: spp_wts_pbl
286
287 !LOCAL
288 real(kind_phys), dimension(im,levs) :: &
289 & sqv,sqc,sqi,sqs,qnc,qni,ozone,qnwfa,qnifa,qnbca, &
290 & dz, w, p, rho, th, qv, delp, &
291 & RUBLTEN, RVBLTEN, RTHBLTEN, RQVBLTEN, &
292 & RQCBLTEN, RQNCBLTEN, RQIBLTEN, RQNIBLTEN, RQSBLTEN, &
293 & RQNWFABLTEN, RQNIFABLTEN, RQNBCABLTEN
294 real(kind_phys), allocatable :: old_ozone(:,:)
295
296!smoke/chem arrays
297 real(kind_phys), dimension(:), intent(inout), optional :: frp
298 logical, intent(in) :: mix_chem, enh_mix, rrfs_sd
299 real(kind_phys), dimension(:,:,:), intent(inout), optional :: chem3d
300 real(kind_phys), dimension(:,: ), intent(in), optional :: vdep
301 real(kind_phys), dimension(im) :: emis_ant_no
302
303!MYNN-2D
304 real(kind_phys), dimension(:), intent(in) :: &
305 & dx,zorl,slmsk,tsurf,qsfc,ps, &
306 & hflx,qflx,ust,wspd,rb,recmol
307
308 real(kind_phys), dimension(:), intent(in), optional :: &
309 & dusfc_cice,dvsfc_cice,dtsfc_cice,dqsfc_cice
310 real(kind_phys), dimension(:), intent(in) :: &
311 & stress_wat,hflx_wat,qflx_wat, &
312 & oceanfrac,fice
313
314 logical, dimension(:), intent(in) :: &
315 & wet, dry, icy
316
317 real(kind_phys), dimension(:), intent(inout) :: &
318 & pblh,dusfc_diag,dvsfc_diag,dtsfc_diag,dqsfc_diag
319 real(kind_phys), dimension(:), intent(out) :: &
320 & ch,dtsfc1,dqsfc1,dusfc1,dvsfc1, &
321 & dtsfci_diag,dqsfci_diag,dusfci_diag,dvsfci_diag
322 real(kind_phys), dimension(:), intent(out) :: &
323 & maxMF,maxwidth,ztop_plume
324 integer, dimension(:), intent(inout) :: &
325 & kpbl
326 integer, dimension(:), intent(inout) :: &
327 & ktop_plume
328
329 real(kind_phys), dimension(:), intent(inout), optional :: &
330 & dusfc_cpl,dvsfc_cpl,dtsfc_cpl,dqsfc_cpl
331 real(kind_phys), dimension(:), intent(inout), optional :: &
332 & dusfci_cpl,dvsfci_cpl,dtsfci_cpl,dqsfci_cpl
333
334 !LOCAL
335 real(kind_phys), dimension(im) :: &
336 & hfx,qfx,rmol,xland,uoce,voce,znt,ts
337 integer :: idtend
338 real(kind_phys), dimension(im) :: dusfci1,dvsfci1,dtsfci1,dqsfci1
339 real(kind_phys), allocatable :: save_qke_adv(:,:)
340 real(kind_phys), dimension(levs) :: kzero
341
342 ! Initialize CCPP error handling variables
343 errmsg = ''
344 errflg = 0
345
346 if (lprnt) then
347 write(0,*)"=============================================="
348 write(0,*)"in mynn wrapper..."
349 write(0,*)"flag_init=",flag_init
350 write(0,*)"flag_restart=",flag_restart
351 endif
352
353 if (.not. flag_for_pbl_generic_tend .and. ldiag3d) then
354 idtend = dtidx(ntke+100,index_of_process_pbl)
355 if (idtend>=1) then
356 allocate(save_qke_adv(im,levs))
357 save_qke_adv=qke_adv
358 endif
359 endif
360
361 ! DH* TODO: Use flag_restart to distinguish which fields need
362 ! to be initialized and which are read from restart files
363 if (flag_init) then
364 initflag=1
365 !print*,"in MYNN, initflag=",initflag
366 else
367 initflag=0
368 !print*,"in MYNN, initflag=",initflag
369 endif
370
371 kzero = zero !generic zero array
372 !initialize arrays for test
373 emis_ant_no = 0.
374
375 flag_ozone = ntoz>0
376
377 ! Assign variables for each microphysics scheme
378 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa) then
379 ! WSM6 or Ferrier-Aligo
380 flag_qi = .true.
381 flag_qni= .false.
382 flag_qc = .true.
383 flag_qnc= .false.
384 flag_qs = .false.
385 flag_qnwfa= .false.
386 flag_qnifa= .false.
387 flag_qnbca= .false.
388 do k=1,levs
389 do i=1,im
390 sqv(i,k) = qgrs_water_vapor(i,k)
391 sqc(i,k) = qgrs_liquid_cloud(i,k)
392 sqi(i,k) = qgrs_ice(i,k)
393 sqs(i,k) = 0.
394 ozone(i,k) = qgrs_ozone(i,k)
395 qnc(i,k) = 0.
396 qni(i,k) = 0.
397 qnwfa(i,k) = 0.
398 qnifa(i,k) = 0.
399 qnbca(i,k) = 0.
400 enddo
401 enddo
402 elseif (imp_physics == imp_physics_nssl ) then
403 ! NSSL
404 flag_qi = .true.
405 flag_qni= .true.
406 flag_qc = .true.
407 flag_qnc= .true.
408 flag_qs = .true.
409 flag_qnwfa= nssl_ccn_on ! ERM: Perhaps could use this field for CCN field?
410 flag_qnifa= .false.
411 flag_qnbca= .false.
412 do k=1,levs
413 do i=1,im
414 sqv(i,k) = qgrs_water_vapor(i,k)
415 sqc(i,k) = qgrs_liquid_cloud(i,k)
416 sqi(i,k) = qgrs_ice(i,k)
417 sqs(i,k) = qgrs_snow(i,k)
418 ozone(i,k) = qgrs_ozone(i,k)
419 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
420 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
421 qnwfa(i,k) = 0.
422 IF ( nssl_ccn_on ) THEN
423 qnwfa(i,k) = qgrs_cccn(i,k)
424 ENDIF
425 qnifa(i,k) = 0.
426 qnbca(i,k) = 0.
427 enddo
428 enddo
429 elseif (imp_physics == imp_physics_thompson) then
430 ! Thompson
431 if(ltaerosol) then
432 flag_qi = .true.
433 flag_qni= .true.
434 flag_qc = .true.
435 flag_qs = .true. !pipe it in, but do not mix
436 flag_qnc= .true.
437 flag_qnwfa= .true.
438 flag_qnifa= .true.
439 flag_qnbca= .false.
440 do k=1,levs
441 do i=1,im
442 sqv(i,k) = qgrs_water_vapor(i,k)
443 sqc(i,k) = qgrs_liquid_cloud(i,k)
444 sqi(i,k) = qgrs_ice(i,k)
445 sqs(i,k) = qgrs_snow(i,k)
446 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
447 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
448 ozone(i,k) = qgrs_ozone(i,k)
449 qnwfa(i,k) = qgrs_water_aer_num_conc(i,k)
450 qnifa(i,k) = qgrs_ice_aer_num_conc(i,k)
451 qnbca(i,k) = 0.
452 enddo
453 enddo
454 else if(mraerosol) then
455 flag_qi = .true.
456 flag_qni= .true.
457 flag_qc = .true.
458 flag_qs = .true.
459 flag_qnc= .true.
460 flag_qnwfa= .false.
461 flag_qnifa= .false.
462 flag_qnbca= .false.
463 do k=1,levs
464 do i=1,im
465 sqv(i,k) = qgrs_water_vapor(i,k)
466 sqc(i,k) = qgrs_liquid_cloud(i,k)
467 sqi(i,k) = qgrs_ice(i,k)
468 sqs(i,k) = qgrs_snow(i,k)
469 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
470 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
471 ozone(i,k) = qgrs_ozone(i,k)
472 qnwfa(i,k) = 0.
473 qnifa(i,k) = 0.
474 qnbca(i,k) = 0.
475 enddo
476 enddo
477 else
478 flag_qi = .true.
479 flag_qni= .true.
480 flag_qc = .true.
481 flag_qs = .true.
482 flag_qnc= .false.
483 flag_qnwfa= .false.
484 flag_qnifa= .false.
485 flag_qnbca= .false.
486 do k=1,levs
487 do i=1,im
488 sqv(i,k) = qgrs_water_vapor(i,k)
489 sqc(i,k) = qgrs_liquid_cloud(i,k)
490 sqi(i,k) = qgrs_ice(i,k)
491 sqs(i,k) = qgrs_snow(i,k)
492 qnc(i,k) = 0.
493 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
494 ozone(i,k) = qgrs_ozone(i,k)
495 qnwfa(i,k) = 0.
496 qnifa(i,k) = 0.
497 qnbca(i,k) = 0.
498 enddo
499 enddo
500 endif
501 elseif (imp_physics == imp_physics_gfdl) then
502 ! GFDL MP
503 flag_qi = .true.
504 flag_qni= .false.
505 flag_qc = .true.
506 flag_qnc= .false.
507 flag_qs = .false.
508 flag_qnwfa= .false.
509 flag_qnifa= .false.
510 flag_qnbca= .false.
511 do k=1,levs
512 do i=1,im
513 sqv(i,k) = qgrs_water_vapor(i,k)
514 sqc(i,k) = qgrs_liquid_cloud(i,k)
515 sqi(i,k) = qgrs_ice(i,k)
516 qnc(i,k) = 0.
517 qni(i,k) = 0.
518 sqs(i,k) = 0.
519 qnwfa(i,k) = 0.
520 qnifa(i,k) = 0.
521 qnbca(i,k) = 0.
522 ozone(i,k) = qgrs_ozone(i,k)
523 enddo
524 enddo
525 else
526 print*,"In MYNN wrapper. Unknown microphysics scheme, imp_physics=",imp_physics
527 print*,"Defaulting to qc and qv species only..."
528 flag_qi = .false.
529 flag_qni= .false.
530 flag_qc = .true.
531 flag_qnc= .false.
532 flag_qs = .false.
533 flag_qnwfa= .false.
534 flag_qnifa= .false.
535 flag_qnbca= .false.
536 do k=1,levs
537 do i=1,im
538 sqv(i,k) = qgrs_water_vapor(i,k)
539 sqc(i,k) = qgrs_liquid_cloud(i,k)
540 sqi(i,k) = 0.
541 sqs(i,k) = 0.
542 qnc(i,k) = 0.
543 qni(i,k) = 0.
544 qnwfa(i,k) = 0.
545 qnifa(i,k) = 0.
546 qnbca(i,k) = 0.
547 ozone(i,k) = qgrs_ozone(i,k)
548 enddo
549 enddo
550 endif
551 if(ldiag3d .and. dtidx(100+ntoz,index_of_process_pbl)>1) then
552 allocate(old_ozone(im,levs))
553 old_ozone = ozone
554 endif
555
556 do k=1,levs
557 do i=1,im
558 th(i,k)=t3d(i,k)/exner(i,k)
559 rho(i,k)=prsl(i,k)/(r_d*t3d(i,k)*(1.+p608*max(sqv(i,k),1e-8)))
560 w(i,k) = -omega(i,k)/(rho(i,k)*grav)
561 enddo
562 enddo
563
564 do k=1,levs
565 do i=1,im
566 tmf(i,k,1)=0.
567 enddo
568 enddo
569
570
571 ! Check incoming moist species to ensure non-negative values
572 ! First, create height difference (dz)
573 do k=1,levs
574 do i=1,im
575 dz(i,k)=(phii(i,k+1) - phii(i,k))*g_inv
576 enddo
577 enddo
578
579 do i=1,im
580 do k=1,levs
581 delp(i,k) = prsi(i,k) - prsi(i,k+1)
582 enddo
583 enddo
584
585 do i=1,im
586 call moisture_check2(levs, delt, &
587 delp(i,:), exner(i,:), &
588 sqv(i,:), sqc(i,:), &
589 sqi(i,:), kzero(:), &
590 t3d(i,:) )
591 enddo
592
593 !intialize more variables
594 do i=1,im
595 if (slmsk(i)==1. .or. slmsk(i)==2.) then !sea/land/ice mask (=0/1/2) in FV3
596 xland(i)=1.0 !but land/water = (1/2) in SFCLAY_mynn
597 else
598 xland(i)=2.0
599 endif
600 uoce(i)=0.0
601 voce(i)=0.0
602 !ust(i) = sqrt(stress(i))
603 ch(i)=0.0
604 hfx(i)=hflx(i)*rho(i,1)*cp
605 qfx(i)=qflx(i)*rho(i,1)
606 !filter bad incoming fluxes
607 if (hfx(i) > 1200.)hfx(i) = 1200.
608 if (hfx(i) < -500.)hfx(i) = -500.
609 if (qfx(i) > .0005)qfx(i) = 0.0005
610 if (qfx(i) < -.0002)qfx(i) = -0.0002
611
612 dtsfc1(i) = hfx(i)
613 dqsfc1(i) = qfx(i)*xlv
614 dusfc1(i) = -1.*rho(i,1)*ust(i)*ust(i)*u(i,1)/wspd(i)
615 dvsfc1(i) = -1.*rho(i,1)*ust(i)*ust(i)*v(i,1)/wspd(i)
616
617 !BWG: diagnostic surface fluxes for scalars & momentum
618 dtsfci_diag(i) = dtsfc1(i)
619 dqsfci_diag(i) = dqsfc1(i)
620 dtsfc_diag(i) = dtsfc_diag(i) + dtsfc1(i)*delt
621 dqsfc_diag(i) = dqsfc_diag(i) + dqsfc1(i)*delt
622 dusfci_diag(i) = dusfc1(i)
623 dvsfci_diag(i) = dvsfc1(i)
624 dusfc_diag(i) = dusfc_diag(i) + dusfci_diag(i)*delt
625 dvsfc_diag(i) = dvsfc_diag(i) + dvsfci_diag(i)*delt
626
627 znt(i)=zorl(i)*0.01 !cm -> m?
628 if (do_mynnsfclay) then
629 rmol(i)=recmol(i)
630 else
631 if (hfx(i) .ge. 0.)then
632 rmol(i)=-hfx(i)/(200.*dz(i,1)*0.5)
633 else
634 rmol(i)=abs(rb(i))*1./(dz(i,1)*0.5)
635 endif
636 endif
637 ts(i)=tsurf(i)/exner(i,1) !theta
638! qsfc(i)=qss(i)
639! ps(i)=pgr(i)
640! wspd(i)=wind(i)
641 enddo
642
643 ! BWG: Coupling insertion
644 if (cplflx) then
645 do i=1,im
646 if (oceanfrac(i) > zero) then ! Ocean only, NO LAKES
647 if ( .not. wet(i)) then ! no open water, use results from CICE
648 dusfci_cpl(i) = dusfc_cice(i)
649 dvsfci_cpl(i) = dvsfc_cice(i)
650 dtsfci_cpl(i) = dtsfc_cice(i)
651 dqsfci_cpl(i) = dqsfc_cice(i)
652 elseif (icy(i) .or. dry(i)) then ! use stress_ocean for opw component at mixed point
653 if (wspd(i) > zero) then
654 dusfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*u(i,1)/wspd(i) ! U-momentum flux
655 dvsfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*v(i,1)/wspd(i) ! V-momentum flux
656 else
657 dusfci_cpl(i) = zero
658 dvsfci_cpl(i) = zero
659 endif
660 dtsfci_cpl(i) = cp*rho(i,1)*hflx_wat(i) ! sensible heat flux over open ocean
661 dqsfci_cpl(i) = xlv*rho(i,1)*qflx_wat(i) ! latent heat flux over open ocean
662 else ! use results from this scheme for 100% open ocean
663 dusfci_cpl(i) = dusfci_diag(i)
664 dvsfci_cpl(i) = dvsfci_diag(i)
665 dtsfci_cpl(i) = dtsfci_diag(i)
666 dqsfci_cpl(i) = dqsfci_diag(i)
667 endif
668!
669 dusfc_cpl(i) = dusfc_cpl(i) + dusfci_cpl(i) * delt
670 dvsfc_cpl(i) = dvsfc_cpl(i) + dvsfci_cpl(i) * delt
671 dtsfc_cpl(i) = dtsfc_cpl(i) + dtsfci_cpl(i) * delt
672 dqsfc_cpl(i) = dqsfc_cpl(i) + dqsfci_cpl(i) * delt
673 else ! If no ocean
674 dusfc_cpl(i) = huge
675 dvsfc_cpl(i) = huge
676 dtsfc_cpl(i) = huge
677 dqsfc_cpl(i) = huge
678 endif ! Ocean only, NO LAKES
679 enddo
680 endif ! End coupling insertion
681
682 if (lprnt) then
683 print*
684 write(0,*)"===CALLING mynn_bl_driver; input:"
685 print*,"tke_budget=",tke_budget," bl_mynn_tkeadvect=",bl_mynn_tkeadvect
686 print*,"bl_mynn_cloudpdf=",bl_mynn_cloudpdf," bl_mynn_mixlength=",bl_mynn_mixlength
687 print*,"bl_mynn_edmf=",bl_mynn_edmf," bl_mynn_edmf_mom=",bl_mynn_edmf_mom
688 print*,"bl_mynn_edmf_tke=",bl_mynn_edmf_tke
689 print*,"bl_mynn_cloudmix=",bl_mynn_cloudmix," bl_mynn_mixqt=",bl_mynn_mixqt
690 print*,"icloud_bl=",icloud_bl
691 print*,"T:",t3d(1,1),t3d(1,2),t3d(1,levs)
692 print*,"TH:",th(1,1),th(1,2),th(1,levs)
693 print*,"rho:",rho(1,1),rho(1,2),rho(1,levs)
694 print*,"exner:",exner(1,1),exner(1,2),exner(1,levs)
695 print*,"prsl:",prsl(1,1),prsl(1,2),prsl(1,levs)
696 print*,"dz:",dz(1,1),dz(1,2),dz(1,levs)
697 print*,"u:",u(1,1),u(1,2),u(1,levs)
698 print*,"v:",v(1,1),v(1,2),v(1,levs)
699 print*,"sqv:",sqv(1,1),sqv(1,2),sqv(1,levs)
700 print*,"sqc:",sqc(1,1),sqc(1,2),sqc(1,levs)
701 print*,"sqi:",sqi(1,1),sqi(1,2),sqi(1,levs)
702 print*,"rmol:",rmol(1)," ust:",ust(1)
703 print*," dx=",dx(1),"initflag=",initflag
704 print*,"Tsurf:",tsurf(1)," Thetasurf:",ts(1)
705 print*,"HFX:",hfx(1)," qfx",qfx(1)
706 print*,"qsfc:",qsfc(1)," ps:",ps(1)
707 print*,"wspd:",wspd(1)," rb=",rb(1)
708 print*,"znt:",znt(1)," delt=",delt
709 print*,"im=",im," levs=",levs
710 print*,"PBLH=",pblh(1)," KPBL=",kpbl(1)," xland=",xland(1)
711 print*,"ch=",ch(1)
712 !print*,"TKE:",TKE_PBL(1,1),TKE_PBL(1,2),TKE_PBL(1,levs)
713 print*,"qke:",qke(1,1),qke(1,2),qke(1,levs)
714 print*,"el_pbl:",el_pbl(1,1),el_pbl(1,2),el_pbl(1,levs)
715 print*,"Sh3d:",sh3d(1,1),sh3d(1,2),sh3d(1,levs)
716 !print*,"exch_h:",exch_h(1,1),exch_h(1,2),exch_h(1,levs) ! - intent(out)
717 !print*,"exch_m:",exch_m(1,1),exch_m(1,2),exch_m(1,levs) ! - intent(out)
718 print*,"max cf_bl:",maxval(cldfra_bl(1,:))
719 endif
720
721
722 CALL mynn_bl_driver( &
723 & initflag=initflag,restart=flag_restart, &
724 & cycling=cycling, &
725 & delt=delt,dz=dz,dx=dx,znt=znt, &
726 & u=u,v=v,w=w,th=th,sqv3d=sqv,sqc3d=sqc, &
727 & sqi3d=sqi,sqs3d=sqs,qnc=qnc,qni=qni, &
728 & qnwfa=qnwfa,qnifa=qnifa,qnbca=qnbca,ozone=ozone, &
729 & p=prsl,exner=exner,rho=rho,t3d=t3d, &
730 & xland=xland,ts=ts,qsfc=qsfc,ps=ps, &
731 & ust=ust,ch=ch,hfx=hfx,qfx=qfx,rmol=rmol, &
732 & wspd=wspd,uoce=uoce,voce=voce, & !input
733 & qke=qke,qke_adv=qke_adv, & !output
734 & sh3d=sh3d,sm3d=sm3d, &
735!chem/smoke
736 & nchem=nchem,kdvel=kdvel,ndvel=ndvel, &
737 & chem3d=chem3d,vdep=vdep,smoke_dbg=smoke_dbg, &
738 & frp=frp,emis_ant_no=emis_ant_no, &
739 & mix_chem=mix_chem,enh_mix=enh_mix, &
740 & rrfs_sd=rrfs_sd, &
741!-----
742 & tsq=tsq,qsq=qsq,cov=cov, & !output
743 & rublten=rublten,rvblten=rvblten,rthblten=rthblten, & !output
744 & rqvblten=rqvblten,rqcblten=rqcblten, &
745 & rqiblten=rqiblten,rqncblten=rqncblten, & !output
746 & rqsblten=rqsblten, & !output
747 & rqniblten=rqniblten,rqnwfablten=rqnwfablten, & !output
748 & rqnifablten=rqnifablten,rqnbcablten=rqnbcablten, & !output
749 & dozone=dqdt_ozone, & !output
750 & exch_h=exch_h,exch_m=exch_m, & !output
751 & pblh=pblh,kpbl=kpbl, & !output
752 & el_pbl=el_pbl, & !output
753 & dqke=dqke, & !output
754 & qwt=qwt,qshear=qshear,qbuoy=qbuoy,qdiss=qdiss, & !output
755 & bl_mynn_tkeadvect=bl_mynn_tkeadvect, &
756 & tke_budget=tke_budget, & !input parameter
757 & bl_mynn_cloudpdf=bl_mynn_cloudpdf, & !input parameter
758 & bl_mynn_mixlength=bl_mynn_mixlength, & !input parameter
759 & icloud_bl=icloud_bl, & !input parameter
760 & qc_bl=qc_bl,qi_bl=qi_bl,cldfra_bl=cldfra_bl, & !output
761 & closure=bl_mynn_closure,bl_mynn_edmf=bl_mynn_edmf, & !input parameter
762 & bl_mynn_edmf_mom=bl_mynn_edmf_mom, & !input parameter
763 & bl_mynn_edmf_tke=bl_mynn_edmf_tke, & !input parameter
764 & bl_mynn_mixscalars=bl_mynn_mixscalars, & !input parameter
765 & bl_mynn_output=bl_mynn_output, & !input parameter
766 & bl_mynn_cloudmix=bl_mynn_cloudmix, & !input parameter
767 & bl_mynn_mixqt=bl_mynn_mixqt, & !input parameter
768 & edmf_a=edmf_a,edmf_w=edmf_w,edmf_qt=edmf_qt, & !output
769 & edmf_thl=edmf_thl,edmf_ent=edmf_ent,edmf_qc=edmf_qc,& !output
770 & sub_thl3d=sub_thl,sub_sqv3d=sub_sqv, &
771 & det_thl3d=det_thl,det_sqv3d=det_sqv, &
772 & maxwidth=maxwidth,maxmf=maxmf,ztop_plume=ztop_plume,& !output
773 & ktop_plume=ktop_plume, & !output
774 & spp_pbl=spp_pbl,pattern_spp_pbl=spp_wts_pbl, & !input
775 & rthraten=htrlw, & !input
776 & flag_qi=flag_qi,flag_qni=flag_qni, & !input
777 & flag_qc=flag_qc,flag_qnc=flag_qnc,flag_qs=flag_qs, & !input
778 & flag_qnwfa=flag_qnwfa,flag_qnifa=flag_qnifa, & !input
779 & flag_qnbca=flag_qnbca,flag_ozone=flag_ozone, & !input
780 & ids=1,ide=im,jds=1,jde=1,kds=1,kde=levs, & !input
781 & ims=1,ime=im,jms=1,jme=1,kms=1,kme=levs, & !input
782 & its=1,ite=im,jts=1,jte=1,kts=1,kte=levs ) !input
783
784
785 ! POST MYNN (INTERSTITIAL) WORK:
786 !update/save MYNN-only variables
787 !do k=1,levs
788 ! do i=1,im
789 ! gq0(i,k,4)=qke(i,k,1) !tke*2
790 ! enddo
791 !enddo
792 !For MYNN, convert TH-tend to T-tend
793 do k = 1, levs
794 do i = 1, im
795 dtdt(i,k) = dtdt(i,k) + rthblten(i,k)*exner(i,k)
796 dudt(i,k) = dudt(i,k) + rublten(i,k)
797 dvdt(i,k) = dvdt(i,k) + rvblten(i,k)
798 enddo
799 enddo
800 accum_duvt3dt: if(ldiag3d .or. lsidea) then
801 call dtend_helper(index_of_x_wind,rublten)
802 call dtend_helper(index_of_y_wind,rvblten)
803 call dtend_helper(index_of_temperature,rthblten,exner)
804 if(ldiag3d) then
805 call dtend_helper(100+ntoz,dqdt_ozone)
806 ! idtend = dtidx(100+ntoz,index_of_process_pbl)
807 ! if(idtend>=1) then
808 ! dtend(:,:,idtend) = dtend(:,:,idtend) + (ozone-old_ozone)
809 ! deallocate(old_ozone)
810 ! endif
811 endif
812 endif accum_duvt3dt
813 !Update T, U and V:
814 !do k = 1, levs
815 ! do i = 1, im
816 ! T3D(i,k) = T3D(i,k) + RTHBLTEN(i,k)*exner(i,k)*delt
817 ! u(i,k) = u(i,k) + RUBLTEN(i,k)*delt
818 ! v(i,k) = v(i,k) + RVBLTEN(i,k)*delt
819 ! enddo
820 !enddo
821
822 !DO moist/scalar/tracer tendencies:
823 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa) then
824 ! WSM6
825 do k=1,levs
826 do i=1,im
827 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
828 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
829 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
830 dqdt_snow(i,k) = rqsblten(i,k) !/(1.0 + qv(i,k))
831 !dqdt_ozone(i,k) = 0.0
832 enddo
833 enddo
834 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
835 call dtend_helper(100+ntqv,rqvblten)
836 call dtend_helper(100+ntcw,rqcblten)
837 call dtend_helper(100+ntiw,rqiblten)
838 endif
839 !Update moist species:
840 !do k=1,levs
841 ! do i=1,im
842 ! qgrs_water_vapor(i,k) = qgrs_water_vapor(i,k) + (RQVBLTEN(i,k)/(1.0+RQVBLTEN(i,k)))*delt
843 ! qgrs_liquid_cloud(i,k) = qgrs_liquid_cloud(i,k) + RQCBLTEN(i,k)*delt
844 ! qgrs_ice(i,k) = qgrs_ice(i,k) + RQIBLTEN(i,k)*delt
845 ! !dqdt_ozone(i,k) = 0.0
846 ! enddo
847 !enddo
848 elseif (imp_physics == imp_physics_thompson) then
849 ! Thompson-Aerosol
850 if(ltaerosol) then
851 do k=1,levs
852 do i=1,im
853 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
854 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
855 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
856 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
857 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
858 dqdt_snow(i,k) = rqsblten(i,k) !/(1.0 + qv(i,k))
859 !dqdt_ozone(i,k) = 0.0
860 dqdt_water_aer_num_conc(i,k) = rqnwfablten(i,k)
861 dqdt_ice_aer_num_conc(i,k) = rqnifablten(i,k)
862 enddo
863 enddo
864 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
865 call dtend_helper(100+ntqv,rqvblten)
866 call dtend_helper(100+ntcw,rqcblten)
867 call dtend_helper(100+ntlnc,rqncblten)
868 call dtend_helper(100+ntiw,rqiblten)
869 call dtend_helper(100+ntinc,rqniblten)
870 call dtend_helper(100+ntwa,rqnwfablten)
871 call dtend_helper(100+ntia,rqnifablten)
872 endif
873 !do k=1,levs
874 ! do i=1,im
875 ! qgrs_water_vapor(i,k) = qgrs_water_vapor(i,k) + (RQVBLTEN(i,k)/(1.0+RQVBLTEN(i,k)))*delt
876 ! qgrs_liquid_cloud(i,k) = qgrs_liquid_cloud(i,k) + RQCBLTEN(i,k)*delt
877 ! qgrs_ice(i,k) = qgrs_ice(i,k) + RQIBLTEN(i,k)*delt
878 ! qgrs_cloud_droplet_num_conc(i,k) = qgrs_cloud_droplet_num_conc(i,k) + RQNCBLTEN(i,k)*delt
879 ! qgrs_cloud_ice_num_conc(i,k) = qgrs_cloud_ice_num_conc(i,k) + RQNIBLTEN(i,k)*delt
880 ! !dqdt_ozone(i,k) = 0.0
881 ! !qgrs_water_aer_num_conc(i,k) = qgrs_water_aer_num_conc(i,k) + RQNWFABLTEN(i,k)*delt
882 ! !qgrs_ice_aer_num_conc(i,k) = qgrs_ice_aer_num_conc(i,k) + RQNIFABLTEN(i,k)*delt
883 ! enddo
884 !enddo
885 else if(mraerosol) then
886 do k=1,levs
887 do i=1,im
888 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
889 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
890 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
891 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
892 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
893 dqdt_snow(i,k) = rqsblten(i,k) !/(1.0 + qv(i,k))
894 enddo
895 enddo
896 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
897 call dtend_helper(100+ntqv,rqvblten)
898 call dtend_helper(100+ntcw,rqcblten)
899 call dtend_helper(100+ntlnc,rqncblten)
900 call dtend_helper(100+ntiw,rqiblten)
901 call dtend_helper(100+ntinc,rqniblten)
902 endif
903 else
904 !Thompson (2008)
905 do k=1,levs
906 do i=1,im
907 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
908 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
909 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
910 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
911 dqdt_snow(i,k) = rqsblten(i,k) !/(1.0 + qv(i,k))
912 !dqdt_ozone(i,k) = 0.0
913 enddo
914 enddo
915 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
916 call dtend_helper(100+ntqv,rqvblten)
917 call dtend_helper(100+ntcw,rqcblten)
918 call dtend_helper(100+ntiw,rqiblten)
919 call dtend_helper(100+ntinc,rqniblten)
920 !call dtend_helper(100+ntsw,RQSBLTEN)
921 endif
922 !do k=1,levs
923 ! do i=1,im
924 ! qgrs_water_vapor(i,k) = qgrs_water_vapor(i,k) + (RQVBLTEN(i,k)/(1.0+RQVBLTEN(i,k)))*delt
925 ! qgrs_liquid_cloud(i,k) = qgrs_liquid_cloud(i,k) + RQCBLTEN(i,k)*delt
926 ! qgrs_ice(i,k) = qgrs_ice(i,k) + RQIBLTEN(i,k)*delt
927 ! qgrs_cloud_ice_num_conc(i,k) = qgrs_cloud_ice_num_conc(i,k) + RQNIBLTEN(i,k)*delt
928 ! !dqdt_ozone(i,k) = 0.0
929 ! enddo
930 !enddo
931 endif !end thompson choice
932 elseif (imp_physics == imp_physics_nssl) then
933 ! NSSL
934 do k=1,levs
935 do i=1,im
936 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
937 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
938 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
939 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
940 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
941 dqdt_snow(i,k) = rqsblten(i,k) !/(1.0 + qv(i,k))
942 IF ( nssl_ccn_on ) THEN !
943 dqdt_cccn(i,k) = rqnwfablten(i,k)
944 ENDIF
945 enddo
946 enddo
947
948 elseif (imp_physics == imp_physics_gfdl) then
949 ! GFDL MP
950 do k=1,levs
951 do i=1,im
952 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
953 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
954 dqdt_ice(i,k) = rqiblten(i,k) !/(1.0 + qv(i,k))
955 !dqdt_rain(i,k) = 0.0
956 !dqdt_snow(i,k) = 0.0
957 !dqdt_graupel(i,k) = 0.0
958 !dqdt_ozone(i,k) = 0.0
959 enddo
960 enddo
961 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
962 call dtend_helper(100+ntqv,rqvblten)
963 call dtend_helper(100+ntcw,rqcblten)
964 call dtend_helper(100+ntiw,rqiblten)
965 endif
966 !do k=1,levs
967 ! do i=1,im
968 ! qgrs_water_vapor(i,k) = qgrs_water_vapor(i,k) + (RQVBLTEN(i,k)/(1.0+RQVBLTEN(i,k)))*delt
969 ! qgrs_liquid_cloud(i,k) = qgrs_liquid_cloud(i,k) + RQCBLTEN(i,k)*delt
970 ! qgrs_ice(i,k) = qgrs_ice(i,k) + RQIBLTEN(i,k)*delt
971 ! !dqdt_ozone(i,k) = 0.0
972 ! enddo
973 !enddo
974 else
975! print*,"In MYNN wrapper. Unknown microphysics scheme, imp_physics=",imp_physics
976 do k=1,levs
977 do i=1,im
978 dqdt_water_vapor(i,k) = rqvblten(i,k) !/(1.0 + qv(i,k))
979 dqdt_liquid_cloud(i,k) = rqcblten(i,k) !/(1.0 + qv(i,k))
980 dqdt_ice(i,k) = 0.0
981 !dqdt_rain(i,k) = 0.0
982 !dqdt_snow(i,k) = 0.0
983 !dqdt_graupel(i,k) = 0.0
984 !dqdt_ozone(i,k) = 0.0
985 enddo
986 enddo
987 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
988 call dtend_helper(100+ntqv,rqvblten)
989 call dtend_helper(100+ntcw,rqcblten)
990 call dtend_helper(100+ntiw,rqiblten)
991 endif
992 endif
993
994 if (lprnt) then
995 print*
996 print*,"===Finished with mynn_bl_driver; output:"
997 print*,"T:",t3d(1,1),t3d(1,2),t3d(1,levs)
998 print*,"TH:",th(1,1),th(1,2),th(1,levs)
999 print*,"rho:",rho(1,1),rho(1,2),rho(1,levs)
1000 print*,"exner:",exner(1,1),exner(1,2),exner(1,levs)
1001 print*,"prsl:",prsl(1,1),prsl(1,2),prsl(1,levs)
1002 print*,"dz:",dz(1,1),dz(1,2),dz(1,levs)
1003 print*,"u:",u(1,1),u(1,2),u(1,levs)
1004 print*,"v:",v(1,1),v(1,2),v(1,levs)
1005 print*,"sqv:",sqv(1,1),sqv(1,2),sqv(1,levs)
1006 print*,"sqc:",sqc(1,1),sqc(1,2),sqc(1,levs)
1007 print*,"sqi:",sqi(1,1),sqi(1,2),sqi(1,levs)
1008 print*,"rmol:",rmol(1)," ust:",ust(1)
1009 print*,"dx(1)=",dx(1),"initflag=",initflag
1010 print*,"Tsurf:",tsurf(1)," Thetasurf:",ts(1)
1011 print*,"HFX:",hfx(1)," qfx",qfx(1)
1012 print*,"qsfc:",qsfc(1)," ps:",ps(1)
1013 print*,"wspd:",wspd(1)," rb=",rb(1)
1014 print*,"znt:",znt(1)," delt=",delt
1015 print*,"im=",im," levs=",levs
1016 print*,"PBLH=",pblh(1)," KPBL=",kpbl(1)," xland=",xland(1)
1017 print*,"ch=",ch(1)
1018 print*,"qke:",qke(1,1),qke(1,2),qke(1,levs)
1019 print*,"el_pbl:",el_pbl(1,1),el_pbl(1,2),el_pbl(1,levs)
1020 print*,"Sh3d:",sh3d(1,1),sh3d(1,2),sh3d(1,levs)
1021 print*,"exch_h:",exch_h(1,1),exch_h(1,2),exch_h(1,levs)
1022 print*,"exch_m:",exch_m(1,1),exch_m(1,2),exch_m(1,levs)
1023 print*,"max cf_bl:",maxval(cldfra_bl(1,:))
1024 print*,"max qc_bl:",maxval(qc_bl(1,:))
1025 print*,"dtdt:",dtdt(1,1),dtdt(1,2),dtdt(1,levs)
1026 print*,"dudt:",dudt(1,1),dudt(1,2),dudt(1,levs)
1027 print*,"dvdt:",dvdt(1,1),dvdt(1,2),dvdt(1,levs)
1028 print*,"dqdt:",dqdt_water_vapor(1,1),dqdt_water_vapor(1,2),dqdt_water_vapor(1,levs)
1029 print*,"ztop_plume:",ztop_plume(1)," maxmf:",maxmf(1)
1030 print*,"maxwidth:",maxwidth(1)
1031 print*
1032 endif
1033
1034 if(allocated(save_qke_adv)) then
1035 if(ldiag3d .and. .not. flag_for_pbl_generic_tend) then
1036 idtend = dtidx(100+ntke,index_of_process_pbl)
1037 if(idtend>=1) then
1038 dtend(:,:,idtend) = dtend(:,:,idtend) + qke_adv-save_qke_adv
1039 endif
1040 endif
1041 deallocate(save_qke_adv)
1042 endif
1043
1044 if(imfdeepcnv == imfdeepcnv_c3 .or. imfdeepcnv == imfdeepcnv_samf)then
1045 !LB: save PBL q-tendency for use in prognostic closure
1046 do k=1,levs
1047 do i=1,im
1048 tmf(i,k,1)=dqdt_water_vapor(i,k)
1049 enddo
1050 enddo
1051 endif
1052
1053 CONTAINS
1054
1055 SUBROUTINE dtend_helper(itracer,field,mult)
1056 real(kind_phys), intent(in) :: field(im,levs)
1057 real(kind_phys), intent(in), optional :: mult(im,levs)
1058 integer, intent(in) :: itracer
1059 integer :: idtend
1060
1061 idtend=dtidx(itracer,index_of_process_pbl)
1062 if(idtend>=1) then
1063 if(present(mult)) then
1064 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf*mult
1065 else
1066 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf
1067 endif
1068 endif
1069 END SUBROUTINE dtend_helper
1070
1071! ==================================================================
1072 SUBROUTINE moisture_check2(kte, delt, dp, exner, &
1073 qv, qc, qi, qs, th )
1074 !
1075 ! If qc < qcmin, qi < qimin, or qv < qvmin happens in any layer,
1076 ! force them to be larger than minimum value by (1) condensating
1077 ! water vapor into liquid or ice, and (2) by transporting water vapor
1078 ! from the very lower layer.
1079 !
1080 ! We then update the final state variables and tendencies associated
1081 ! with this correction. If any condensation happens, update theta/temperature too.
1082 ! Note that (qv,qc,qi,th) are the final state variables after
1083 ! applying corresponding input tendencies and corrective tendencies.
1084
1085 implicit none
1086 integer, intent(in) :: kte
1087 real(kind_phys), intent(in) :: delt
1088 real(kind_phys), dimension(kte), intent(in) :: dp, exner
1089 real(kind_phys), dimension(kte), intent(inout) :: qv, qc, qi, qs, th
1090 integer k
1091 real :: dqc2, dqi2, dqs2, dqv2, sum, aa, dum
1092 real, parameter :: qvmin1= 1e-8, & !min at k=1
1093 qvmin = 1e-20, & !min above k=1
1094 qcmin = 0.0, &
1095 qimin = 0.0
1096
1097 do k = kte, 1, -1 ! From the top to the surface
1098 dqc2 = max(0.0, qcmin-qc(k)) !qc deficit (>=0)
1099 dqi2 = max(0.0, qimin-qi(k)) !qi deficit (>=0)
1100 dqs2 = max(0.0, qimin-qs(k)) !qs deficit (>=0)
1101
1102 !update species
1103 qc(k) = qc(k) + dqc2
1104 qi(k) = qi(k) + dqi2
1105 qs(k) = qs(k) + dqs2
1106 qv(k) = qv(k) - dqc2 - dqi2 - dqs2
1107 !for theta
1108 !th(k) = th(k) + xlvcp/exner(k)*dqc2 + &
1109 ! xlscp/exner(k)*dqi2
1110 !for temperature
1111 th(k) = th(k) + xlvcp*dqc2 + &
1112 xlscp*(dqi2+dqs2)
1113
1114 !then fix qv if lending qv made it negative
1115 if (k .eq. 1) then
1116 dqv2 = max(0.0, qvmin1-qv(k)) !qv deficit (>=0)
1117 qv(k) = qv(k) + dqv2
1118 qv(k) = max(qv(k),qvmin1)
1119 dqv2 = 0.0
1120 else
1121 dqv2 = max(0.0, qvmin-qv(k)) !qv deficit (>=0)
1122 qv(k) = qv(k) + dqv2
1123 qv(k-1)= qv(k-1) - dqv2*dp(k)/dp(k-1)
1124 qv(k) = max(qv(k),qvmin)
1125 endif
1126 qc(k) = max(qc(k),qcmin)
1127 qi(k) = max(qi(k),qimin)
1128 qs(k) = max(qs(k),qimin)
1129 end do
1130
1131 ! Extra moisture used to satisfy 'qv(1)>=qvmin' is proportionally
1132 ! extracted from all the layers that has 'qv > 2*qvmin'. This fully
1133 ! preserves column moisture.
1134 if( dqv2 .gt. 1.e-20 ) then
1135 sum = 0.0
1136 do k = 1, kte
1137 if( qv(k) .gt. 2.0*qvmin ) sum = sum + qv(k)*dp(k)
1138 enddo
1139 aa = dqv2*dp(1)/max(1.e-20,sum)
1140 if( aa .lt. 0.5 ) then
1141 do k = 1, kte
1142 if( qv(k) .gt. 2.0*qvmin ) then
1143 dum = aa*qv(k)
1144 qv(k) = qv(k) - dum
1145 endif
1146 enddo
1147 else
1148 ! For testing purposes only (not yet found in any output):
1149 ! write(*,*) 'Full moisture conservation is impossible'
1150 endif
1151 endif
1152
1153 return
1154
1155 END SUBROUTINE moisture_check2
1156
1157 END SUBROUTINE mynnedmf_wrapper_run
1158
1159!###=================================================================
1160
1161END MODULE mynnedmf_wrapper
subroutine mynn_bl_driver(initflag, restart, cycling, delt, dz, dx, znt, u, v, w, th, sqv3d, sqc3d, sqi3d, sqs3d, qnc, qni, qnwfa, qnifa, qnbca, ozone, p, exner, rho, t3d, xland, ts, qsfc, ps, ust, ch, hfx, qfx, rmol, wspd, uoce, voce, qke, qke_adv, sh3d, sm3d, nchem, kdvel, ndvel, chem3d, vdep, smoke_dbg, frp, emis_ant_no, mix_chem, enh_mix, rrfs_sd, tsq, qsq, cov, rublten, rvblten, rthblten, rqvblten, rqcblten, rqiblten, rqncblten, rqniblten, rqsblten, rqnwfablten, rqnifablten, rqnbcablten, dozone, exch_h, exch_m, pblh, kpbl, el_pbl, dqke, qwt, qshear, qbuoy, qdiss, qc_bl, qi_bl, cldfra_bl, bl_mynn_tkeadvect, tke_budget, bl_mynn_cloudpdf, bl_mynn_mixlength, icloud_bl, closure, bl_mynn_edmf, bl_mynn_edmf_mom, bl_mynn_edmf_tke, bl_mynn_mixscalars, bl_mynn_output, bl_mynn_cloudmix, bl_mynn_mixqt, edmf_a, edmf_w, edmf_qt, edmf_thl, edmf_ent, edmf_qc, sub_thl3d, sub_sqv3d, det_thl3d, det_sqv3d, maxwidth, maxmf, ztop_plume, ktop_plume, spp_pbl, pattern_spp_pbl, rthraten, flag_qc, flag_qi, flag_qnc, flag_qni, flag_qs, flag_qnwfa, flag_qnifa, flag_qnbca, flag_ozone, ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, its, ite, jts, jte, kts, kte)
This subroutine is the MYNN-EDNF PBL driver routine,which encompassed the majority of the subroutines...
This module defines model-specific constants/parameters.
This module contains the entity of MYNN-EDMF PBL scheme.
The following references best describe the code within Olson et al. (2019, NOAA Technical Memorandum)...