91SUBROUTINE mynnedmf_wrapper_run( &
93 & flag_init,flag_restart, &
94 & lssav, ldiag3d, qdiag3d, &
97 & phii,u,v,omega,t3d, &
99 & qgrs_liquid_cloud, &
102 & qgrs_cloud_droplet_num_conc, &
103 & qgrs_cloud_ice_num_conc, &
105 & qgrs_water_aer_num_conc, &
106 & qgrs_ice_aer_num_conc, &
109 & slmsk,tsurf,qsfc,ps, &
110 & ust,ch,hflx,qflx,wspd,rb, &
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, &
126 & qke,qke_adv,Tsq,Qsq,Cov, &
127 & el_pbl,sh3d,sm3d,exch_h,exch_m, &
128 & dqke,qwt,qshear,qbuoy,qdiss, &
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, &
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, &
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, &
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 )
167 use machine,
only: kind_phys
169 xlv, xlvcp, xlscp, p608
176 real(kind_phys),
intent(in) :: huge
177 character(len=*),
intent(out) :: errmsg
178 integer,
intent(out) :: errflg
180 logical,
intent(in) :: lssav, ldiag3d, lsidea, qdiag3d
181 logical,
intent(in) :: cplflx
184 integer,
intent(in) :: nchem, ndvel
185 integer,
parameter :: kdvel=1
186 logical,
intent(in) :: smoke_dbg
189 logical,
intent(in) :: &
190 & bl_mynn_tkeadvect, &
195 & flag_for_pbl_generic_tend, &
197 integer,
intent(in) :: &
198 & bl_mynn_cloudpdf, &
199 & bl_mynn_mixlength, &
202 & bl_mynn_edmf_mom, &
203 & bl_mynn_edmf_tke, &
204 & bl_mynn_cloudmix, &
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, &
213 real(kind_phys),
intent(in) :: &
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
225 INTEGER,
PARAMETER :: &
226 & bl_mynn_mixscalars=1
228 & FLAG_QI, FLAG_QNI, FLAG_QC, FLAG_QS, FLAG_QNC, &
229 & FLAG_QNWFA, FLAG_QNIFA, FLAG_QNBCA, FLAG_OZONE
231 LOGICAL,
PARAMETER :: cycling = .false.
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
242 REAL(kind_phys) :: tem
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 :: &
257 real(kind_phys),
dimension(:,:),
intent(inout) :: &
259 real(kind_phys),
dimension(:,:,:),
intent(out) :: tmf
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, &
268 real(kind_phys),
dimension(:,:),
intent(in) :: &
269 & qgrs_cloud_ice_num_conc, &
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
285 real(kind_phys),
dimension(:,:),
intent(in),
optional :: spp_wts_pbl
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(:,:)
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
304 real(kind_phys),
dimension(:),
intent(in) :: &
305 & dx,zorl,slmsk,tsurf,qsfc,ps, &
306 & hflx,qflx,ust,wspd,rb,recmol
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, &
314 logical,
dimension(:),
intent(in) :: &
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) :: &
326 integer,
dimension(:),
intent(inout) :: &
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
335 real(kind_phys),
dimension(im) :: &
336 & hfx,qfx,rmol,xland,uoce,voce,znt,ts
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
347 write(0,*)
"=============================================="
348 write(0,*)
"in mynn wrapper..."
349 write(0,*)
"flag_init=",flag_init
350 write(0,*)
"flag_restart=",flag_restart
353 if (.not. flag_for_pbl_generic_tend .and. ldiag3d)
then
354 idtend = dtidx(ntke+100,index_of_process_pbl)
356 allocate(save_qke_adv(im,levs))
378 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa)
then
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)
394 ozone(i,k) = qgrs_ozone(i,k)
402 elseif (imp_physics == imp_physics_nssl )
then
409 flag_qnwfa= nssl_ccn_on
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)
422 IF ( nssl_ccn_on )
THEN
423 qnwfa(i,k) = qgrs_cccn(i,k)
429 elseif (imp_physics == imp_physics_thompson)
then
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)
454 else if(mraerosol)
then
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)
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)
493 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
494 ozone(i,k) = qgrs_ozone(i,k)
501 elseif (imp_physics == imp_physics_gfdl)
then
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)
522 ozone(i,k) = qgrs_ozone(i,k)
526 print*,
"In MYNN wrapper. Unknown microphysics scheme, imp_physics=",imp_physics
527 print*,
"Defaulting to qc and qv species only..."
538 sqv(i,k) = qgrs_water_vapor(i,k)
539 sqc(i,k) = qgrs_liquid_cloud(i,k)
547 ozone(i,k) = qgrs_ozone(i,k)
551 if(ldiag3d .and. dtidx(100+ntoz,index_of_process_pbl)>1)
then
552 allocate(old_ozone(im,levs))
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)
575 dz(i,k)=(phii(i,k+1) - phii(i,k))*g_inv
581 delp(i,k) = prsi(i,k) - prsi(i,k+1)
586 call moisture_check2(levs, delt, &
587 delp(i,:), exner(i,:), &
588 sqv(i,:), sqc(i,:), &
589 sqi(i,:), kzero(:), &
595 if (slmsk(i)==1. .or. slmsk(i)==2.)
then
604 hfx(i)=hflx(i)*rho(i,1)*cp
605 qfx(i)=qflx(i)*rho(i,1)
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
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)
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
628 if (do_mynnsfclay)
then
631 if (hfx(i) .ge. 0.)
then
632 rmol(i)=-hfx(i)/(200.*dz(i,1)*0.5)
634 rmol(i)=abs(rb(i))*1./(dz(i,1)*0.5)
637 ts(i)=tsurf(i)/exner(i,1)
646 if (oceanfrac(i) > zero)
then
647 if ( .not. wet(i))
then
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
653 if (wspd(i) > zero)
then
654 dusfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*u(i,1)/wspd(i)
655 dvsfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*v(i,1)/wspd(i)
660 dtsfci_cpl(i) = cp*rho(i,1)*hflx_wat(i)
661 dqsfci_cpl(i) = xlv*rho(i,1)*qflx_wat(i)
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)
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
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)
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)
718 print*,
"max cf_bl:",maxval(cldfra_bl(1,:))
723 & initflag=initflag,restart=flag_restart, &
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, &
733 & qke=qke,qke_adv=qke_adv, &
734 & sh3d=sh3d,sm3d=sm3d, &
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, &
742 & tsq=tsq,qsq=qsq,cov=cov, &
743 & rublten=rublten,rvblten=rvblten,rthblten=rthblten, &
744 & rqvblten=rqvblten,rqcblten=rqcblten, &
745 & rqiblten=rqiblten,rqncblten=rqncblten, &
746 & rqsblten=rqsblten, &
747 & rqniblten=rqniblten,rqnwfablten=rqnwfablten, &
748 & rqnifablten=rqnifablten,rqnbcablten=rqnbcablten, &
749 & dozone=dqdt_ozone, &
750 & exch_h=exch_h,exch_m=exch_m, &
751 & pblh=pblh,kpbl=kpbl, &
754 & qwt=qwt,qshear=qshear,qbuoy=qbuoy,qdiss=qdiss, &
755 & bl_mynn_tkeadvect=bl_mynn_tkeadvect, &
756 & tke_budget=tke_budget, &
757 & bl_mynn_cloudpdf=bl_mynn_cloudpdf, &
758 & bl_mynn_mixlength=bl_mynn_mixlength, &
759 & icloud_bl=icloud_bl, &
760 & qc_bl=qc_bl,qi_bl=qi_bl,cldfra_bl=cldfra_bl, &
761 & closure=bl_mynn_closure,bl_mynn_edmf=bl_mynn_edmf, &
762 & bl_mynn_edmf_mom=bl_mynn_edmf_mom, &
763 & bl_mynn_edmf_tke=bl_mynn_edmf_tke, &
764 & bl_mynn_mixscalars=bl_mynn_mixscalars, &
765 & bl_mynn_output=bl_mynn_output, &
766 & bl_mynn_cloudmix=bl_mynn_cloudmix, &
767 & bl_mynn_mixqt=bl_mynn_mixqt, &
768 & edmf_a=edmf_a,edmf_w=edmf_w,edmf_qt=edmf_qt, &
769 & edmf_thl=edmf_thl,edmf_ent=edmf_ent,edmf_qc=edmf_qc,&
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,&
773 & ktop_plume=ktop_plume, &
774 & spp_pbl=spp_pbl,pattern_spp_pbl=spp_wts_pbl, &
776 & flag_qi=flag_qi,flag_qni=flag_qni, &
777 & flag_qc=flag_qc,flag_qnc=flag_qnc,flag_qs=flag_qs, &
778 & flag_qnwfa=flag_qnwfa,flag_qnifa=flag_qnifa, &
779 & flag_qnbca=flag_qnbca,flag_ozone=flag_ozone, &
780 & ids=1,ide=im,jds=1,jde=1,kds=1,kde=levs, &
781 & ims=1,ime=im,jms=1,jme=1,kms=1,kme=levs, &
782 & its=1,ite=im,jts=1,jte=1,kts=1,kte=levs )
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)
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)
805 call dtend_helper(100+ntoz,dqdt_ozone)
823 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa)
then
827 dqdt_water_vapor(i,k) = rqvblten(i,k)
828 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
829 dqdt_ice(i,k) = rqiblten(i,k)
830 dqdt_snow(i,k) = rqsblten(i,k)
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)
848 elseif (imp_physics == imp_physics_thompson)
then
853 dqdt_water_vapor(i,k) = rqvblten(i,k)
854 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
855 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
856 dqdt_ice(i,k) = rqiblten(i,k)
857 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
858 dqdt_snow(i,k) = rqsblten(i,k)
860 dqdt_water_aer_num_conc(i,k) = rqnwfablten(i,k)
861 dqdt_ice_aer_num_conc(i,k) = rqnifablten(i,k)
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)
885 else if(mraerosol)
then
888 dqdt_water_vapor(i,k) = rqvblten(i,k)
889 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
890 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
891 dqdt_ice(i,k) = rqiblten(i,k)
892 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
893 dqdt_snow(i,k) = rqsblten(i,k)
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)
907 dqdt_water_vapor(i,k) = rqvblten(i,k)
908 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
909 dqdt_ice(i,k) = rqiblten(i,k)
910 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
911 dqdt_snow(i,k) = rqsblten(i,k)
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)
932 elseif (imp_physics == imp_physics_nssl)
then
936 dqdt_water_vapor(i,k) = rqvblten(i,k)
937 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
938 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
939 dqdt_ice(i,k) = rqiblten(i,k)
940 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
941 dqdt_snow(i,k) = rqsblten(i,k)
942 IF ( nssl_ccn_on )
THEN
943 dqdt_cccn(i,k) = rqnwfablten(i,k)
948 elseif (imp_physics == imp_physics_gfdl)
then
952 dqdt_water_vapor(i,k) = rqvblten(i,k)
953 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
954 dqdt_ice(i,k) = rqiblten(i,k)
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)
978 dqdt_water_vapor(i,k) = rqvblten(i,k)
979 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
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)
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)
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)
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)
1038 dtend(:,:,idtend) = dtend(:,:,idtend) + qke_adv-save_qke_adv
1041 deallocate(save_qke_adv)
1044 if(imfdeepcnv == imfdeepcnv_c3 .or. imfdeepcnv == imfdeepcnv_samf)
then
1048 tmf(i,k,1)=dqdt_water_vapor(i,k)
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
1061 idtend=dtidx(itracer,index_of_process_pbl)
1063 if(
present(mult))
then
1064 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf*mult
1066 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf
1069 END SUBROUTINE dtend_helper
1072 SUBROUTINE moisture_check2(kte, delt, dp, exner, &
1073 qv, qc, qi, qs, th )
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
1091 real :: dqc2, dqi2, dqs2, dqv2, sum, aa, dum
1092 real,
parameter :: qvmin1= 1e-8, &
1098 dqc2 = max(0.0, qcmin-qc(k))
1099 dqi2 = max(0.0, qimin-qi(k))
1100 dqs2 = max(0.0, qimin-qs(k))
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
1111 th(k) = th(k) + xlvcp*dqc2 + &
1116 dqv2 = max(0.0, qvmin1-qv(k))
1117 qv(k) = qv(k) + dqv2
1118 qv(k) = max(qv(k),qvmin1)
1121 dqv2 = max(0.0, qvmin-qv(k))
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)
1126 qc(k) = max(qc(k),qcmin)
1127 qi(k) = max(qi(k),qimin)
1128 qs(k) = max(qs(k),qimin)
1134 if( dqv2 .gt. 1.e-20 )
then
1137 if( qv(k) .gt. 2.0*qvmin ) sum = sum + qv(k)*dp(k)
1139 aa = dqv2*dp(1)/max(1.e-20,sum)
1140 if( aa .lt. 0.5 )
then
1142 if( qv(k) .gt. 2.0*qvmin )
then
1155 END SUBROUTINE moisture_check2