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 & 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, &
159 & chem3d, frp, mix_chem, rrfs_sd, enh_mix, &
160 & nchem, ndvel, vdep, smoke_dbg, &
161 & imp_physics_nssl, nssl_ccn_on, &
162 & ltaerosol, mraerosol, spp_wts_pbl, spp_pbl, &
163 & lprnt, huge, errmsg, errflg )
166 use machine,
only: kind_phys
168 xlv, xlvcp, xlscp, p608
175 real(kind_phys),
intent(in) :: huge
176 character(len=*),
intent(out) :: errmsg
177 integer,
intent(out) :: errflg
179 logical,
intent(in) :: lssav, ldiag3d, lsidea, qdiag3d
180 logical,
intent(in) :: cplflx
183 integer,
intent(in) :: nchem, ndvel
184 integer,
parameter :: kdvel=1
185 logical,
intent(in) :: smoke_dbg
188 logical,
intent(in) :: &
189 & bl_mynn_tkeadvect, &
194 & flag_for_pbl_generic_tend, &
196 integer,
intent(in) :: &
197 & bl_mynn_cloudpdf, &
198 & bl_mynn_mixlength, &
201 & bl_mynn_edmf_mom, &
202 & bl_mynn_edmf_tke, &
203 & bl_mynn_cloudmix, &
206 & imp_physics, imp_physics_wsm6, &
207 & imp_physics_thompson, imp_physics_gfdl, &
208 & imp_physics_nssl, imp_physics_fa, &
211 real(kind_phys),
intent(in) :: &
215 real(kind_phys),
intent(inout),
optional :: dtend(:,:,:)
216 integer,
intent(in) :: dtidx(:,:)
217 integer,
intent(in) :: index_of_temperature, index_of_x_wind
218 integer,
intent(in) :: index_of_y_wind, index_of_process_pbl
219 integer,
intent(in) :: ntoz, ntqv, ntcw, ntiw, ntsw, ntlnc
220 integer,
intent(in) :: ntinc, ntwa, ntia, ntke
223 INTEGER,
PARAMETER :: &
224 & bl_mynn_mixscalars=1
226 & FLAG_QI, FLAG_QNI, FLAG_QC, FLAG_QS, FLAG_QNC, &
227 & FLAG_QNWFA, FLAG_QNIFA, FLAG_QNBCA, FLAG_OZONE
229 LOGICAL,
PARAMETER :: cycling = .false.
232 REAL(kind_phys),
intent(in) :: delt, dtf
233 INTEGER,
intent(in) :: im, levs
234 LOGICAL,
intent(in) :: flag_init, flag_restart
235 INTEGER :: initflag, k, i
236 INTEGER :: IDS,IDE,JDS,JDE,KDS,KDE, &
237 & IMS,IME,JMS,JME,KMS,KME, &
238 & ITS,ITE,JTS,JTE,KTS,KTE
240 REAL(kind_phys) :: tem
243 real(kind_phys),
dimension(:,:),
intent(in) :: phii
244 real(kind_phys),
dimension(:,:),
intent(inout) :: &
245 & dtdt, dudt, dvdt, &
246 & dqdt_water_vapor, dqdt_liquid_cloud, dqdt_ice, &
247 & dqdt_snow, dqdt_ice_num_conc, dqdt_ozone
248 real(kind_phys),
dimension(:,:),
intent(inout),
optional :: &
249 & dqdt_cloud_droplet_num_conc, dqdt_water_aer_num_conc, &
250 & dqdt_ice_aer_num_conc
251 real(kind_phys),
dimension(:,:),
intent(inout),
optional :: qke, &
252 & EL_PBL, Sh3D, Sm3D, qc_bl, qi_bl, cldfra_bl, dqdt_cccn
253 real(kind_phys),
dimension(:,:),
intent(inout) :: &
256 real(kind_phys),
dimension(:,:),
intent(inout),
optional :: &
257 & edmf_a,edmf_w,edmf_qt, &
258 & edmf_thl,edmf_ent,edmf_qc, &
259 & sub_thl,sub_sqv,det_thl,det_sqv
260 real(kind_phys),
dimension(:,:),
intent(inout) :: &
261 & t3d,qgrs_water_vapor,qgrs_liquid_cloud,qgrs_ice, &
263 real(kind_phys),
dimension(:,:),
intent(in) :: &
264 & qgrs_cloud_ice_num_conc, &
268 real(kind_phys),
dimension(:,:),
intent(in),
optional :: &
269 & qgrs_water_aer_num_conc, &
270 & qgrs_cloud_droplet_num_conc, &
271 & qgrs_ice_aer_num_conc
272 real(kind_phys),
dimension(:,:),
intent(in),
optional :: qgrs_cccn
273 real(kind_phys),
dimension(:,:),
intent(out),
optional :: &
274 & Tsq, Qsq, Cov, exch_h, exch_m, dqke, qWT, qSHEAR, qBUOY, &
276 real(kind_phys),
dimension(:),
intent(in) :: xmu
277 real(kind_phys),
dimension(:,:),
intent(in) :: htrsw, htrlw
279 real(kind_phys),
dimension(:,:),
intent(in),
optional :: spp_wts_pbl
282 real(kind_phys),
dimension(im,levs) :: &
283 & sqv,sqc,sqi,sqs,qnc,qni,ozone,qnwfa,qnifa,qnbca, &
284 & dz, w, p, rho, th, qv, delp, &
285 & RUBLTEN, RVBLTEN, RTHBLTEN, RQVBLTEN, &
286 & RQCBLTEN, RQNCBLTEN, RQIBLTEN, RQNIBLTEN, RQSBLTEN, &
287 & RQNWFABLTEN, RQNIFABLTEN, RQNBCABLTEN
288 real(kind_phys),
allocatable :: old_ozone(:,:)
291 real(kind_phys),
dimension(:),
intent(inout),
optional :: frp
292 logical,
intent(in) :: mix_chem, enh_mix, rrfs_sd
293 real(kind_phys),
dimension(:,:,:),
intent(inout),
optional :: chem3d
294 real(kind_phys),
dimension(:,: ),
intent(in),
optional :: vdep
295 real(kind_phys),
dimension(im) :: emis_ant_no
298 real(kind_phys),
dimension(:),
intent(in) :: &
299 & dx,zorl,slmsk,tsurf,qsfc,ps, &
300 & hflx,qflx,ust,wspd,rb,recmol
302 real(kind_phys),
dimension(:),
intent(in),
optional :: &
303 & dusfc_cice,dvsfc_cice,dtsfc_cice,dqsfc_cice
304 real(kind_phys),
dimension(:),
intent(in) :: &
305 & stress_wat,hflx_wat,qflx_wat, &
308 logical,
dimension(:),
intent(in) :: &
311 real(kind_phys),
dimension(:),
intent(inout) :: &
312 & pblh,dusfc_diag,dvsfc_diag,dtsfc_diag,dqsfc_diag
313 real(kind_phys),
dimension(:),
intent(out) :: &
314 & ch,dtsfc1,dqsfc1,dusfc1,dvsfc1, &
315 & dtsfci_diag,dqsfci_diag,dusfci_diag,dvsfci_diag
316 real(kind_phys),
dimension(:),
intent(out),
optional :: &
317 & maxMF,maxwidth,ztop_plume
318 integer,
dimension(:),
intent(inout) :: &
320 integer,
dimension(:),
intent(inout),
optional :: &
323 real(kind_phys),
dimension(:),
intent(inout),
optional :: &
324 & dusfc_cpl,dvsfc_cpl,dtsfc_cpl,dqsfc_cpl
325 real(kind_phys),
dimension(:),
intent(inout),
optional :: &
326 & dusfci_cpl,dvsfci_cpl,dtsfci_cpl,dqsfci_cpl
329 real(kind_phys),
dimension(im) :: &
330 & hfx,qfx,rmol,xland,uoce,voce,znt,ts
332 real(kind_phys),
dimension(im) :: dusfci1,dvsfci1,dtsfci1,dqsfci1
333 real(kind_phys),
allocatable :: save_qke_adv(:,:)
334 real(kind_phys),
dimension(levs) :: kzero
341 write(0,*)
"=============================================="
342 write(0,*)
"in mynn wrapper..."
343 write(0,*)
"flag_init=",flag_init
344 write(0,*)
"flag_restart=",flag_restart
347 if (.not. flag_for_pbl_generic_tend .and. ldiag3d)
then
348 idtend = dtidx(ntke+100,index_of_process_pbl)
350 allocate(save_qke_adv(im,levs))
372 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa)
then
384 sqv(i,k) = qgrs_water_vapor(i,k)
385 sqc(i,k) = qgrs_liquid_cloud(i,k)
386 sqi(i,k) = qgrs_ice(i,k)
388 ozone(i,k) = qgrs_ozone(i,k)
396 elseif (imp_physics == imp_physics_nssl )
then
403 flag_qnwfa= nssl_ccn_on
408 sqv(i,k) = qgrs_water_vapor(i,k)
409 sqc(i,k) = qgrs_liquid_cloud(i,k)
410 sqi(i,k) = qgrs_ice(i,k)
411 sqs(i,k) = qgrs_snow(i,k)
412 ozone(i,k) = qgrs_ozone(i,k)
413 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
414 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
416 IF ( nssl_ccn_on )
THEN
417 qnwfa(i,k) = qgrs_cccn(i,k)
423 elseif (imp_physics == imp_physics_thompson)
then
436 sqv(i,k) = qgrs_water_vapor(i,k)
437 sqc(i,k) = qgrs_liquid_cloud(i,k)
438 sqi(i,k) = qgrs_ice(i,k)
439 sqs(i,k) = qgrs_snow(i,k)
440 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
441 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
442 ozone(i,k) = qgrs_ozone(i,k)
443 qnwfa(i,k) = qgrs_water_aer_num_conc(i,k)
444 qnifa(i,k) = qgrs_ice_aer_num_conc(i,k)
448 else if(mraerosol)
then
459 sqv(i,k) = qgrs_water_vapor(i,k)
460 sqc(i,k) = qgrs_liquid_cloud(i,k)
461 sqi(i,k) = qgrs_ice(i,k)
462 sqs(i,k) = qgrs_snow(i,k)
463 qnc(i,k) = qgrs_cloud_droplet_num_conc(i,k)
464 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
465 ozone(i,k) = qgrs_ozone(i,k)
482 sqv(i,k) = qgrs_water_vapor(i,k)
483 sqc(i,k) = qgrs_liquid_cloud(i,k)
484 sqi(i,k) = qgrs_ice(i,k)
485 sqs(i,k) = qgrs_snow(i,k)
487 qni(i,k) = qgrs_cloud_ice_num_conc(i,k)
488 ozone(i,k) = qgrs_ozone(i,k)
495 elseif (imp_physics == imp_physics_gfdl)
then
507 sqv(i,k) = qgrs_water_vapor(i,k)
508 sqc(i,k) = qgrs_liquid_cloud(i,k)
509 sqi(i,k) = qgrs_ice(i,k)
516 ozone(i,k) = qgrs_ozone(i,k)
520 print*,
"In MYNN wrapper. Unknown microphysics scheme, imp_physics=",imp_physics
521 print*,
"Defaulting to qc and qv species only..."
532 sqv(i,k) = qgrs_water_vapor(i,k)
533 sqc(i,k) = qgrs_liquid_cloud(i,k)
541 ozone(i,k) = qgrs_ozone(i,k)
545 if(ldiag3d .and. dtidx(100+ntoz,index_of_process_pbl)>1)
then
546 allocate(old_ozone(im,levs))
552 th(i,k)=t3d(i,k)/exner(i,k)
553 rho(i,k)=prsl(i,k)/(r_d*t3d(i,k)*(1.+p608*max(sqv(i,k),1e-8)))
554 w(i,k) = -omega(i,k)/(rho(i,k)*grav)
562 dz(i,k)=(phii(i,k+1) - phii(i,k))*g_inv
568 delp(i,k) = prsi(i,k) - prsi(i,k+1)
573 call moisture_check2(levs, delt, &
574 delp(i,:), exner(i,:), &
575 sqv(i,:), sqc(i,:), &
576 sqi(i,:), kzero(:), &
582 if (slmsk(i)==1. .or. slmsk(i)==2.)
then
591 hfx(i)=hflx(i)*rho(i,1)*cp
592 qfx(i)=qflx(i)*rho(i,1)
594 if (hfx(i) > 1200.)hfx(i) = 1200.
595 if (hfx(i) < -500.)hfx(i) = -500.
596 if (qfx(i) > .0005)qfx(i) = 0.0005
597 if (qfx(i) < -.0002)qfx(i) = -0.0002
600 dqsfc1(i) = qfx(i)*xlv
601 dusfc1(i) = -1.*rho(i,1)*ust(i)*ust(i)*u(i,1)/wspd(i)
602 dvsfc1(i) = -1.*rho(i,1)*ust(i)*ust(i)*v(i,1)/wspd(i)
605 dtsfci_diag(i) = dtsfc1(i)
606 dqsfci_diag(i) = dqsfc1(i)
607 dtsfc_diag(i) = dtsfc_diag(i) + dtsfc1(i)*delt
608 dqsfc_diag(i) = dqsfc_diag(i) + dqsfc1(i)*delt
609 dusfci_diag(i) = dusfc1(i)
610 dvsfci_diag(i) = dvsfc1(i)
611 dusfc_diag(i) = dusfc_diag(i) + dusfci_diag(i)*delt
612 dvsfc_diag(i) = dvsfc_diag(i) + dvsfci_diag(i)*delt
615 if (do_mynnsfclay)
then
618 if (hfx(i) .ge. 0.)
then
619 rmol(i)=-hfx(i)/(200.*dz(i,1)*0.5)
621 rmol(i)=abs(rb(i))*1./(dz(i,1)*0.5)
624 ts(i)=tsurf(i)/exner(i,1)
633 if (oceanfrac(i) > zero)
then
634 if ( .not. wet(i))
then
635 dusfci_cpl(i) = dusfc_cice(i)
636 dvsfci_cpl(i) = dvsfc_cice(i)
637 dtsfci_cpl(i) = dtsfc_cice(i)
638 dqsfci_cpl(i) = dqsfc_cice(i)
639 elseif (icy(i) .or. dry(i))
then
640 if (wspd(i) > zero)
then
641 dusfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*u(i,1)/wspd(i)
642 dvsfci_cpl(i) = -1.*rho(i,1)*stress_wat(i)*v(i,1)/wspd(i)
647 dtsfci_cpl(i) = cp*rho(i,1)*hflx_wat(i)
648 dqsfci_cpl(i) = xlv*rho(i,1)*qflx_wat(i)
650 dusfci_cpl(i) = dusfci_diag(i)
651 dvsfci_cpl(i) = dvsfci_diag(i)
652 dtsfci_cpl(i) = dtsfci_diag(i)
653 dqsfci_cpl(i) = dqsfci_diag(i)
656 dusfc_cpl(i) = dusfc_cpl(i) + dusfci_cpl(i) * delt
657 dvsfc_cpl(i) = dvsfc_cpl(i) + dvsfci_cpl(i) * delt
658 dtsfc_cpl(i) = dtsfc_cpl(i) + dtsfci_cpl(i) * delt
659 dqsfc_cpl(i) = dqsfc_cpl(i) + dqsfci_cpl(i) * delt
671 write(0,*)
"===CALLING mynn_bl_driver; input:"
672 print*,
"tke_budget=",tke_budget,
" bl_mynn_tkeadvect=",bl_mynn_tkeadvect
673 print*,
"bl_mynn_cloudpdf=",bl_mynn_cloudpdf,
" bl_mynn_mixlength=",bl_mynn_mixlength
674 print*,
"bl_mynn_edmf=",bl_mynn_edmf,
" bl_mynn_edmf_mom=",bl_mynn_edmf_mom
675 print*,
"bl_mynn_edmf_tke=",bl_mynn_edmf_tke
676 print*,
"bl_mynn_cloudmix=",bl_mynn_cloudmix,
" bl_mynn_mixqt=",bl_mynn_mixqt
677 print*,
"icloud_bl=",icloud_bl
678 print*,
"T:",t3d(1,1),t3d(1,2),t3d(1,levs)
679 print*,
"TH:",th(1,1),th(1,2),th(1,levs)
680 print*,
"rho:",rho(1,1),rho(1,2),rho(1,levs)
681 print*,
"exner:",exner(1,1),exner(1,2),exner(1,levs)
682 print*,
"prsl:",prsl(1,1),prsl(1,2),prsl(1,levs)
683 print*,
"dz:",dz(1,1),dz(1,2),dz(1,levs)
684 print*,
"u:",u(1,1),u(1,2),u(1,levs)
685 print*,
"v:",v(1,1),v(1,2),v(1,levs)
686 print*,
"sqv:",sqv(1,1),sqv(1,2),sqv(1,levs)
687 print*,
"sqc:",sqc(1,1),sqc(1,2),sqc(1,levs)
688 print*,
"sqi:",sqi(1,1),sqi(1,2),sqi(1,levs)
689 print*,
"rmol:",rmol(1),
" ust:",ust(1)
690 print*,
" dx=",dx(1),
"initflag=",initflag
691 print*,
"Tsurf:",tsurf(1),
" Thetasurf:",ts(1)
692 print*,
"HFX:",hfx(1),
" qfx",qfx(1)
693 print*,
"qsfc:",qsfc(1),
" ps:",ps(1)
694 print*,
"wspd:",wspd(1),
" rb=",rb(1)
695 print*,
"znt:",znt(1),
" delt=",delt
696 print*,
"im=",im,
" levs=",levs
697 print*,
"PBLH=",pblh(1),
" KPBL=",kpbl(1),
" xland=",xland(1)
700 print*,
"qke:",qke(1,1),qke(1,2),qke(1,levs)
701 print*,
"el_pbl:",el_pbl(1,1),el_pbl(1,2),el_pbl(1,levs)
702 print*,
"Sh3d:",sh3d(1,1),sh3d(1,2),sh3d(1,levs)
705 print*,
"max cf_bl:",maxval(cldfra_bl(1,:))
710 & initflag=initflag,restart=flag_restart, &
712 & delt=delt,dz=dz,dx=dx,znt=znt, &
713 & u=u,v=v,w=w,th=th,sqv3d=sqv,sqc3d=sqc, &
714 & sqi3d=sqi,sqs3d=sqs,qnc=qnc,qni=qni, &
715 & qnwfa=qnwfa,qnifa=qnifa,qnbca=qnbca,ozone=ozone, &
716 & p=prsl,exner=exner,rho=rho,t3d=t3d, &
717 & xland=xland,ts=ts,qsfc=qsfc,ps=ps, &
718 & ust=ust,ch=ch,hfx=hfx,qfx=qfx,rmol=rmol, &
719 & wspd=wspd,uoce=uoce,voce=voce, &
720 & qke=qke,qke_adv=qke_adv, &
721 & sh3d=sh3d,sm3d=sm3d, &
723 & nchem=nchem,kdvel=kdvel,ndvel=ndvel, &
724 & chem3d=chem3d,vdep=vdep,smoke_dbg=smoke_dbg, &
725 & frp=frp,emis_ant_no=emis_ant_no, &
726 & mix_chem=mix_chem,enh_mix=enh_mix, &
729 & tsq=tsq,qsq=qsq,cov=cov, &
730 & rublten=rublten,rvblten=rvblten,rthblten=rthblten, &
731 & rqvblten=rqvblten,rqcblten=rqcblten, &
732 & rqiblten=rqiblten,rqncblten=rqncblten, &
733 & rqsblten=rqsblten, &
734 & rqniblten=rqniblten,rqnwfablten=rqnwfablten, &
735 & rqnifablten=rqnifablten,rqnbcablten=rqnbcablten, &
736 & dozone=dqdt_ozone, &
737 & exch_h=exch_h,exch_m=exch_m, &
738 & pblh=pblh,kpbl=kpbl, &
741 & qwt=qwt,qshear=qshear,qbuoy=qbuoy,qdiss=qdiss, &
742 & bl_mynn_tkeadvect=bl_mynn_tkeadvect, &
743 & tke_budget=tke_budget, &
744 & bl_mynn_cloudpdf=bl_mynn_cloudpdf, &
745 & bl_mynn_mixlength=bl_mynn_mixlength, &
746 & icloud_bl=icloud_bl, &
747 & qc_bl=qc_bl,qi_bl=qi_bl,cldfra_bl=cldfra_bl, &
748 & closure=bl_mynn_closure,bl_mynn_edmf=bl_mynn_edmf, &
749 & bl_mynn_edmf_mom=bl_mynn_edmf_mom, &
750 & bl_mynn_edmf_tke=bl_mynn_edmf_tke, &
751 & bl_mynn_mixscalars=bl_mynn_mixscalars, &
752 & bl_mynn_output=bl_mynn_output, &
753 & bl_mynn_cloudmix=bl_mynn_cloudmix, &
754 & bl_mynn_mixqt=bl_mynn_mixqt, &
755 & edmf_a=edmf_a,edmf_w=edmf_w,edmf_qt=edmf_qt, &
756 & edmf_thl=edmf_thl,edmf_ent=edmf_ent,edmf_qc=edmf_qc,&
757 & sub_thl3d=sub_thl,sub_sqv3d=sub_sqv, &
758 & det_thl3d=det_thl,det_sqv3d=det_sqv, &
759 & maxwidth=maxwidth,maxmf=maxmf,ztop_plume=ztop_plume,&
760 & ktop_plume=ktop_plume, &
761 & spp_pbl=spp_pbl,pattern_spp_pbl=spp_wts_pbl, &
763 & flag_qi=flag_qi,flag_qni=flag_qni, &
764 & flag_qc=flag_qc,flag_qnc=flag_qnc,flag_qs=flag_qs, &
765 & flag_qnwfa=flag_qnwfa,flag_qnifa=flag_qnifa, &
766 & flag_qnbca=flag_qnbca,flag_ozone=flag_ozone, &
767 & ids=1,ide=im,jds=1,jde=1,kds=1,kde=levs, &
768 & ims=1,ime=im,jms=1,jme=1,kms=1,kme=levs, &
769 & its=1,ite=im,jts=1,jte=1,kts=1,kte=levs )
782 dtdt(i,k) = dtdt(i,k) + rthblten(i,k)*exner(i,k)
783 dudt(i,k) = dudt(i,k) + rublten(i,k)
784 dvdt(i,k) = dvdt(i,k) + rvblten(i,k)
787 accum_duvt3dt:
if(ldiag3d .or. lsidea)
then
788 call dtend_helper(index_of_x_wind,rublten)
789 call dtend_helper(index_of_y_wind,rvblten)
790 call dtend_helper(index_of_temperature,rthblten,exner)
792 call dtend_helper(100+ntoz,dqdt_ozone)
810 if (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_fa)
then
814 dqdt_water_vapor(i,k) = rqvblten(i,k)
815 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
816 dqdt_ice(i,k) = rqiblten(i,k)
817 dqdt_snow(i,k) = rqsblten(i,k)
821 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
822 call dtend_helper(100+ntqv,rqvblten)
823 call dtend_helper(100+ntcw,rqcblten)
824 call dtend_helper(100+ntiw,rqiblten)
835 elseif (imp_physics == imp_physics_thompson)
then
840 dqdt_water_vapor(i,k) = rqvblten(i,k)
841 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
842 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
843 dqdt_ice(i,k) = rqiblten(i,k)
844 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
845 dqdt_snow(i,k) = rqsblten(i,k)
847 dqdt_water_aer_num_conc(i,k) = rqnwfablten(i,k)
848 dqdt_ice_aer_num_conc(i,k) = rqnifablten(i,k)
851 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
852 call dtend_helper(100+ntqv,rqvblten)
853 call dtend_helper(100+ntcw,rqcblten)
854 call dtend_helper(100+ntlnc,rqncblten)
855 call dtend_helper(100+ntiw,rqiblten)
856 call dtend_helper(100+ntinc,rqniblten)
857 call dtend_helper(100+ntwa,rqnwfablten)
858 call dtend_helper(100+ntia,rqnifablten)
872 else if(mraerosol)
then
875 dqdt_water_vapor(i,k) = rqvblten(i,k)
876 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
877 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
878 dqdt_ice(i,k) = rqiblten(i,k)
879 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
880 dqdt_snow(i,k) = rqsblten(i,k)
883 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
884 call dtend_helper(100+ntqv,rqvblten)
885 call dtend_helper(100+ntcw,rqcblten)
886 call dtend_helper(100+ntlnc,rqncblten)
887 call dtend_helper(100+ntiw,rqiblten)
888 call dtend_helper(100+ntinc,rqniblten)
894 dqdt_water_vapor(i,k) = rqvblten(i,k)
895 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
896 dqdt_ice(i,k) = rqiblten(i,k)
897 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
898 dqdt_snow(i,k) = rqsblten(i,k)
902 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
903 call dtend_helper(100+ntqv,rqvblten)
904 call dtend_helper(100+ntcw,rqcblten)
905 call dtend_helper(100+ntiw,rqiblten)
906 call dtend_helper(100+ntinc,rqniblten)
919 elseif (imp_physics == imp_physics_nssl)
then
923 dqdt_water_vapor(i,k) = rqvblten(i,k)
924 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
925 dqdt_cloud_droplet_num_conc(i,k) = rqncblten(i,k)
926 dqdt_ice(i,k) = rqiblten(i,k)
927 dqdt_ice_num_conc(i,k) = rqniblten(i,k)
928 dqdt_snow(i,k) = rqsblten(i,k)
929 IF ( nssl_ccn_on )
THEN
930 dqdt_cccn(i,k) = rqnwfablten(i,k)
935 elseif (imp_physics == imp_physics_gfdl)
then
939 dqdt_water_vapor(i,k) = rqvblten(i,k)
940 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
941 dqdt_ice(i,k) = rqiblten(i,k)
948 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
949 call dtend_helper(100+ntqv,rqvblten)
950 call dtend_helper(100+ntcw,rqcblten)
951 call dtend_helper(100+ntiw,rqiblten)
965 dqdt_water_vapor(i,k) = rqvblten(i,k)
966 dqdt_liquid_cloud(i,k) = rqcblten(i,k)
974 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
975 call dtend_helper(100+ntqv,rqvblten)
976 call dtend_helper(100+ntcw,rqcblten)
977 call dtend_helper(100+ntiw,rqiblten)
983 print*,
"===Finished with mynn_bl_driver; output:"
984 print*,
"T:",t3d(1,1),t3d(1,2),t3d(1,levs)
985 print*,
"TH:",th(1,1),th(1,2),th(1,levs)
986 print*,
"rho:",rho(1,1),rho(1,2),rho(1,levs)
987 print*,
"exner:",exner(1,1),exner(1,2),exner(1,levs)
988 print*,
"prsl:",prsl(1,1),prsl(1,2),prsl(1,levs)
989 print*,
"dz:",dz(1,1),dz(1,2),dz(1,levs)
990 print*,
"u:",u(1,1),u(1,2),u(1,levs)
991 print*,
"v:",v(1,1),v(1,2),v(1,levs)
992 print*,
"sqv:",sqv(1,1),sqv(1,2),sqv(1,levs)
993 print*,
"sqc:",sqc(1,1),sqc(1,2),sqc(1,levs)
994 print*,
"sqi:",sqi(1,1),sqi(1,2),sqi(1,levs)
995 print*,
"rmol:",rmol(1),
" ust:",ust(1)
996 print*,
"dx(1)=",dx(1),
"initflag=",initflag
997 print*,
"Tsurf:",tsurf(1),
" Thetasurf:",ts(1)
998 print*,
"HFX:",hfx(1),
" qfx",qfx(1)
999 print*,
"qsfc:",qsfc(1),
" ps:",ps(1)
1000 print*,
"wspd:",wspd(1),
" rb=",rb(1)
1001 print*,
"znt:",znt(1),
" delt=",delt
1002 print*,
"im=",im,
" levs=",levs
1003 print*,
"PBLH=",pblh(1),
" KPBL=",kpbl(1),
" xland=",xland(1)
1005 print*,
"qke:",qke(1,1),qke(1,2),qke(1,levs)
1006 print*,
"el_pbl:",el_pbl(1,1),el_pbl(1,2),el_pbl(1,levs)
1007 print*,
"Sh3d:",sh3d(1,1),sh3d(1,2),sh3d(1,levs)
1008 print*,
"exch_h:",exch_h(1,1),exch_h(1,2),exch_h(1,levs)
1009 print*,
"exch_m:",exch_m(1,1),exch_m(1,2),exch_m(1,levs)
1010 print*,
"max cf_bl:",maxval(cldfra_bl(1,:))
1011 print*,
"max qc_bl:",maxval(qc_bl(1,:))
1012 print*,
"dtdt:",dtdt(1,1),dtdt(1,2),dtdt(1,levs)
1013 print*,
"dudt:",dudt(1,1),dudt(1,2),dudt(1,levs)
1014 print*,
"dvdt:",dvdt(1,1),dvdt(1,2),dvdt(1,levs)
1015 print*,
"dqdt:",dqdt_water_vapor(1,1),dqdt_water_vapor(1,2),dqdt_water_vapor(1,levs)
1016 print*,
"ztop_plume:",ztop_plume(1),
" maxmf:",maxmf(1)
1017 print*,
"maxwidth:",maxwidth(1)
1021 if(
allocated(save_qke_adv))
then
1022 if(ldiag3d .and. .not. flag_for_pbl_generic_tend)
then
1023 idtend = dtidx(100+ntke,index_of_process_pbl)
1025 dtend(:,:,idtend) = dtend(:,:,idtend) + qke_adv-save_qke_adv
1028 deallocate(save_qke_adv)
1033 SUBROUTINE dtend_helper(itracer,field,mult)
1034 real(kind_phys),
intent(in) :: field(im,levs)
1035 real(kind_phys),
intent(in),
optional :: mult(im,levs)
1036 integer,
intent(in) :: itracer
1039 idtend=dtidx(itracer,index_of_process_pbl)
1041 if(
present(mult))
then
1042 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf*mult
1044 dtend(:,:,idtend) = dtend(:,:,idtend) + field*dtf
1047 END SUBROUTINE dtend_helper
1050 SUBROUTINE moisture_check2(kte, delt, dp, exner, &
1051 qv, qc, qi, qs, th )
1064 integer,
intent(in) :: kte
1065 real(kind_phys),
intent(in) :: delt
1066 real(kind_phys),
dimension(kte),
intent(in) :: dp, exner
1067 real(kind_phys),
dimension(kte),
intent(inout) :: qv, qc, qi, qs, th
1069 real :: dqc2, dqi2, dqs2, dqv2, sum, aa, dum
1070 real,
parameter :: qvmin1= 1e-8, &
1076 dqc2 = max(0.0, qcmin-qc(k))
1077 dqi2 = max(0.0, qimin-qi(k))
1078 dqs2 = max(0.0, qimin-qs(k))
1081 qc(k) = qc(k) + dqc2
1082 qi(k) = qi(k) + dqi2
1083 qs(k) = qs(k) + dqs2
1084 qv(k) = qv(k) - dqc2 - dqi2 - dqs2
1089 th(k) = th(k) + xlvcp*dqc2 + &
1094 dqv2 = max(0.0, qvmin1-qv(k))
1095 qv(k) = qv(k) + dqv2
1096 qv(k) = max(qv(k),qvmin1)
1099 dqv2 = max(0.0, qvmin-qv(k))
1100 qv(k) = qv(k) + dqv2
1101 qv(k-1)= qv(k-1) - dqv2*dp(k)/dp(k-1)
1102 qv(k) = max(qv(k),qvmin)
1104 qc(k) = max(qc(k),qcmin)
1105 qi(k) = max(qi(k),qimin)
1106 qs(k) = max(qs(k),qimin)
1112 if( dqv2 .gt. 1.e-20 )
then
1115 if( qv(k) .gt. 2.0*qvmin ) sum = sum + qv(k)*dp(k)
1117 aa = dqv2*dp(1)/max(1.e-20,sum)
1118 if( aa .lt. 0.5 )
then
1120 if( qv(k) .gt. 2.0*qvmin )
then
1133 END SUBROUTINE moisture_check2