1 !> \file module_bl_mynn_wrapper.F90
2 !! This serves as the interface between the WRF PBL driver and the MYNN
3 !! eddy-diffusivity mass-flux scheme in module_bl_mynn.F.
5 !>\ingroup gsd_mynn_edmf
6 !> The following references best describe the code within
7 !! Olson et al. (2019, NOAA Technical Memorandum)
8 !! Nakanishi and Niino (2009) \cite NAKANISHI_2009
9 MODULE module_bl_mynn_wrapper
11 use module_bl_mynn_common
15 !> \section arg_table_mynnedmf_wrapper_init Argument Table
16 !! \htmlinclude mynnedmf_wrapper_init.html
18 subroutine mynnedmf_wrapper_init ( &
19 & RUBLTEN,RVBLTEN,RTHBLTEN,RQVBLTEN,RQCBLTEN,&
21 & restart,allowed_to_read, &
22 & P_QC,P_QI,PARAM_FIRST_SCALAR, &
23 & IDS,IDE,JDS,JDE,KDS,KDE, &
24 & IMS,IME,JMS,JME,KMS,KME, &
25 & ITS,ITE,JTS,JTE,KTS,KTE )
29 LOGICAL,INTENT(IN) :: ALLOWED_TO_READ,RESTART
31 INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE, &
32 & IMS,IME,JMS,JME,KMS,KME, &
33 & ITS,ITE,JTS,JTE,KTS,KTE
36 REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: &
37 &RUBLTEN,RVBLTEN,RTHBLTEN,RQVBLTEN, &
38 &RQCBLTEN,RQIBLTEN,QKE
40 INTEGER, intent(in) :: P_QC,P_QI,PARAM_FIRST_SCALAR
42 INTEGER :: I,J,K,ITF,JTF,KTF
48 IF (.NOT.RESTART) THEN
56 if( p_qc >= param_first_scalar ) RQCBLTEN(i,k,j)=0.
57 if( p_qi >= param_first_scalar ) RQIBLTEN(i,k,j)=0.
63 end subroutine mynnedmf_wrapper_init
65 subroutine mynnedmf_wrapper_finalize ()
66 end subroutine mynnedmf_wrapper_finalize
68 ! \brief This scheme (1) performs pre-mynnedmf work, (2) runs the mynnedmf, and (3) performs post-mynnedmf work
69 !> \section arg_table_mynnedmf_wrapper_run Argument Table
70 !! \htmlinclude mynnedmf_wrapper_run.html
72 SUBROUTINE mynnedmf_wrapper_run( &
73 & initflag,restart,cycling, &
76 & qv,qc,qi,qs,qnc,qni,qnwfa,qnifa,qnbca, &
80 & ust,ch,hfx,qfx,rmol,wspd, &
82 & qke,qke_adv,sh3d,sm3d, &
85 & mix_chem,chem3d,vd3d,nchem,kdvel, &
86 & ndvel,num_vert_mix, &
87 ! & frp_mean,emis_ant_no,enh_mix, & !to be included soon
91 & rublten,rvblten,rthblten, &
92 & rqvblten,rqcblten,rqiblten,rqsblten, &
93 & rqncblten,rqniblten, &
94 & rqnwfablten,rqnifablten,rqnbcablten, &
96 & exch_h,exch_m,pblh,kpbl,el_pbl, &
97 & dqke,qwt,qshear,qbuoy,qdiss, &
98 & qc_bl,qi_bl,cldfra_bl, &
99 & edmf_a,edmf_w,edmf_qt, &
100 & edmf_thl,edmf_ent,edmf_qc, &
101 & sub_thl3d,sub_sqv3d, &
102 & det_thl3d,det_sqv3d, &
103 & maxwidth,maxMF,ztop_plume,ktop_plume, &
105 & tke_budget, bl_mynn_tkeadvect, &
106 & bl_mynn_cloudpdf, bl_mynn_mixlength, &
107 & icloud_bl, bl_mynn_edmf, &
108 & bl_mynn_edmf_mom, bl_mynn_edmf_tke, &
109 & bl_mynn_cloudmix, bl_mynn_mixqt, &
110 & bl_mynn_output, bl_mynn_closure, &
111 & bl_mynn_mixscalars, &
112 & spp_pbl,pattern_spp_pbl, &
113 & flag_qc,flag_qi,flag_qs, &
114 & flag_qnc,flag_qni, &
115 & flag_qnwfa,flag_qnifa,flag_qnbca, &
116 & ids,ide,jds,jde,kds,kde, &
117 & ims,ime,jms,jme,kms,kme, &
118 & its,ite,jts,jte,kts,kte )
120 use module_bl_mynn, only: mynn_bl_driver
122 !-------------------------------------------------------------------
124 !-------------------------------------------------------------------
126 !smoke/chem: disclaimer: all smoke-related variables are still
127 !considered under development in CCPP. Until that work is
128 !completed, these flags/arrays must be kept hard-coded as is.
130 logical, intent(in) :: mix_chem
131 integer, intent(in) :: nchem, ndvel, kdvel, num_vert_mix
132 logical, parameter :: &
133 & rrfs_sd =.false., &
134 & smoke_dbg =.false., &
137 logical, parameter :: &
138 & mix_chem =.false., &
139 & enh_mix =.false., &
140 & rrfs_sd =.false., &
142 integer, parameter :: nchem=2, ndvel=2, kdvel=1, &
146 ! NAMELIST OPTIONS (INPUT):
147 logical, intent(in) :: &
148 & bl_mynn_tkeadvect, &
150 integer, intent(in) :: &
151 & bl_mynn_cloudpdf, &
152 & bl_mynn_mixlength, &
155 & bl_mynn_edmf_mom, &
156 & bl_mynn_edmf_tke, &
157 & bl_mynn_cloudmix, &
160 & bl_mynn_mixscalars, &
163 real(kind_phys), intent(in) :: &
166 logical, intent(in) :: &
167 & FLAG_QI, FLAG_QNI, FLAG_QC, FLAG_QNC, &
168 & FLAG_QS, FLAG_QNWFA, FLAG_QNIFA, FLAG_QNBCA
169 logical, parameter :: FLAG_OZONE = .false.
172 REAL(kind_phys), intent(in) :: delt, dxc
173 LOGICAL, intent(in) :: restart
174 INTEGER :: i, j, k, itf, jtf, ktf, n
175 INTEGER, intent(in) :: initflag, &
176 & IDS,IDE,JDS,JDE,KDS,KDE, &
177 & IMS,IME,JMS,JME,KMS,KME, &
178 & ITS,ITE,JTS,JTE,KTS,KTE
181 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), intent(in) :: &
182 & u,v,w,t3d,th,rho,exner,p,dz
183 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), intent(inout) :: &
184 & rublten,rvblten,rthblten, &
185 & rqvblten,rqcblten,rqiblten,rqsblten, &
186 & rqncblten,rqniblten, &
187 & rqnwfablten,rqnifablten,rqnbcablten !,ro3blten
188 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), intent(inout) :: &
189 & qke, qke_adv, el_pbl, sh3d, sm3d
190 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), intent(inout) :: &
191 & Tsq, Qsq, Cov, exch_h, exch_m
192 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), intent(in) :: rthraten
195 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), optional, intent(in) :: &
197 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), optional, intent(inout) :: &
198 & qc_bl, qi_bl, cldfra_bl
199 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), optional, intent(inout) :: &
200 & edmf_a,edmf_w,edmf_qt, &
201 & edmf_thl,edmf_ent,edmf_qc, &
202 & sub_thl3d,sub_sqv3d,det_thl3d,det_sqv3d
203 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), optional, intent(inout) :: &
204 & dqke,qWT,qSHEAR,qBUOY,qDISS
205 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme), optional, intent(inout) :: &
206 & qv,qc,qi,qs,qnc,qni,qnwfa,qnifa,qnbca!,o3
208 !optional 2D arrays for passing into module_bl_myn.F
209 real(kind_phys), allocatable, dimension(:,:) :: &
210 & qc_bl2d, qi_bl2d, cldfra_bl2d, pattern_spp_pbl2d
211 real(kind_phys), allocatable, dimension(:,:) :: &
212 & edmf_a2d,edmf_w2d,edmf_qt2d, &
213 & edmf_thl2d,edmf_ent2d,edmf_qc2d, &
214 & sub_thl2d,sub_sqv2d,det_thl2d,det_sqv2d
215 real(kind_phys), allocatable, dimension(:,:) :: &
216 & dqke2d,qWT2d,qSHEAR2d,qBUOY2d,qDISS2d
217 real(kind_phys), allocatable, dimension(:,:) :: &
218 & qc2d,qi2d,qs2d,qnc2d,qni2d,qnwfa2d,qnifa2d,qnbca2d!,o32d
220 !smoke/chem arrays - no if-defs in module_bl_mynn.F, so must define arrays
222 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme,nchem), intent(in) :: chem3d
223 real(kind_phys), dimension(ims:ime,kdvel,jms:jme, ndvel), intent(in) :: vd3d
224 real(kind_phys), dimension(ims:ime,kms:kme,nchem) :: chem
225 real(kind_phys), dimension(ims:ime,ndvel) :: vd
226 real(kind_phys), dimension(ims:ime) :: frp_mean, emis_ant_no
228 real(kind_phys), dimension(ims:ime,kms:kme,nchem) :: chem
229 real(kind_phys), dimension(ims:ime,ndvel) :: vd
230 real(kind_phys), dimension(ims:ime) :: frp_mean, emis_ant_no
234 real(kind_phys), dimension(ims:ime,jms:jme), intent(in) :: &
235 & xland,ts,qsfc,ps,ch
236 real(kind_phys), dimension(ims:ime,jms:jme), intent(inout) :: &
237 & znt,pblh,maxwidth,maxmf,ztop_plume,rmol,hfx,qfx,ust,wspd, &
239 integer, dimension(ims:ime,jms:jme), intent(inout) :: &
243 real(kind_phys), dimension(ims:ime,kms:kme) :: delp,sqv,sqc,sqi,sqs,ikzero
244 real(kind_phys), dimension(ims:ime) :: dx
245 logical, parameter :: debug = .false.
246 real(kind_phys), dimension(ims:ime,kms:kme,jms:jme) :: ozone,rO3blten
248 !write(0,*)"=============================================="
249 !write(0,*)"in mynn wrapper..."
250 !write(0,*)"initflag=",initflag
251 !write(0,*)"restart =",restart
257 !For now, initialized bogus array
262 !Allocate any arrays being used
263 if (icloud_bl > 0) then
264 allocate(qc_bl2d(ims:ime,kms:kme))
265 allocate(qi_bl2d(ims:ime,kms:kme))
266 allocate(cldfra_bl2d(ims:ime,kms:kme))
271 if (spp_pbl > 0) then
272 allocate(pattern_spp_pbl2d(ims:ime,kms:kme))
274 if (bl_mynn_output > 0) then
275 allocate(edmf_a2d(ims:ime,kms:kme))
276 allocate(edmf_w2d(ims:ime,kms:kme))
277 allocate(edmf_qt2d(ims:ime,kms:kme))
278 allocate(edmf_thl2d(ims:ime,kms:kme))
279 allocate(edmf_ent2d(ims:ime,kms:kme))
280 allocate(edmf_qc2d(ims:ime,kms:kme))
281 allocate(sub_thl2d(ims:ime,kms:kme))
282 allocate(sub_sqv2d(ims:ime,kms:kme))
283 allocate(det_thl2d(ims:ime,kms:kme))
284 allocate(det_sqv2d(ims:ime,kms:kme))
286 if (tke_budget .eq. 1) then
287 allocate(dqke2d(ims:ime,kms:kme))
288 allocate(qWT2d(ims:ime,kms:kme))
289 allocate(qSHEAR2d(ims:ime,kms:kme))
290 allocate(qBUOY2d(ims:ime,kms:kme))
291 allocate(qDISS2d(ims:ime,kms:kme))
299 allocate(qc2d(ims:ime,kms:kme))
303 allocate(qi2d(ims:ime,kms:kme))
307 allocate(qs2d(ims:ime,kms:kme))
311 allocate(qnc2d(ims:ime,kms:kme))
315 allocate(qni2d(ims:ime,kms:kme))
319 allocate(qnwfa2d(ims:ime,kms:kme))
323 allocate(qnifa2d(ims:ime,kms:kme))
327 allocate(qnbca2d(ims:ime,kms:kme))
330 !---------------------------------
331 !Begin looping in the j-direction
332 !---------------------------------
335 !need sgs cloud info input for diagnostic-decay
336 if (icloud_bl > 0) then
339 qc_bl2d(i,k) = qc_bl(i,k,j)
340 qi_bl2d(i,k) = qi_bl(i,k,j)
341 cldfra_bl2d(i,k) = cldfra_bl(i,k,j)
347 if (spp_pbl > 0) then
350 pattern_spp_pbl2d(i,k) = pattern_spp_pbl(i,k,j)
355 !intialize moist species
359 qc2d(i,k) = qc(i,k,j)
366 qi2d(i,k) = qi(i,k,j)
373 qs2d(i,k) = qs(i,k,j)
380 qnc2d(i,k) = qnc(i,k,j)
387 qni2d(i,k) = qni(i,k,j)
394 qnwfa2d(i,k) = qnwfa(i,k,j)
401 qnifa2d(i,k) = qnifa(i,k,j)
408 qnbca2d(i,k) = qnbca(i,k,j)
418 chem(i,k,n)=chem3d(i,k,j,n)
426 vd(i,n) =vd3d(i,1,j,n)
439 ! Check incoming moist species to ensure non-negative values
440 ! First, create pressure differences (delp) across model layers
446 ! call moisture_check2(kte, delt, &
447 ! delp(i,:), exner(i,:,j), &
448 ! qv(i,:,j), qc(i,:,j), &
449 ! qi(i,:,j), t3d(i,:,j) )
452 !In WRF, mixing ratio is incoming. Convert to specific humidity:
455 sqv(i,k)=qv(i,k,j)/(1.0 + qv(i,k,j))
456 sqc(i,k)=qc2d(i,k)/(1.0 + qv(i,k,j))
462 sqi(i,k)=qi2d(i,k)/(1.0 + qv(i,k,j))
471 sqs(i,k)=qs2d(i,k)/(1.0 + qv(i,k,j))
480 write(0,*)"===CALLING mynn_bl_driver; input:"
481 print*,"tke_budget=",tke_budget
482 print*,"bl_mynn_tkeadvect=",bl_mynn_tkeadvect
483 print*,"bl_mynn_cloudpdf=",bl_mynn_cloudpdf
484 print*,"bl_mynn_mixlength=",bl_mynn_mixlength
485 print*,"bl_mynn_edmf=",bl_mynn_edmf
486 print*,"bl_mynn_edmf_mom=",bl_mynn_edmf_mom
487 print*,"bl_mynn_edmf_tke=",bl_mynn_edmf_tke
488 print*,"bl_mynn_cloudmix=",bl_mynn_cloudmix
489 print*,"bl_mynn_mixqt=",bl_mynn_mixqt
490 print*,"icloud_bl=",icloud_bl
491 print*,"T:",t3d(its,1,j),t3d(its,2,j),t3d(its,kte,j)
492 print*,"TH:",th(its,1,j),th(its,2,j),th(its,kte,j)
493 print*,"rho:",rho(its,1,j),rho(its,2,j),rho(its,kte,j)
494 print*,"exner:",exner(its,1,j),exner(its,2,j),exner(its,kte,j)
495 print*,"p:",p(its,1,j),p(its,2,j),p(its,kte,j)
496 print*,"dz:",dz(its,1,j),dz(its,2,j),dz(its,kte,j)
497 print*,"u:",u(its,1,j),u(its,2,j),u(its,kte,j)
498 print*,"v:",v(its,1,j),v(its,2,j),v(its,kte,j)
499 print*,"sqv:",sqv(its,1),sqv(its,2),sqv(its,kte)
500 print*,"sqc:",sqc(its,1),sqc(its,2),sqc(its,kte)
501 print*,"sqi:",sqi(its,1),sqi(its,2),sqi(its,kte)
502 print*,"rmol:",rmol(its,j)," ust:",ust(its,j)
503 print*,"dx=",dx(its),"initflag=",initflag
504 print*,"Thetasurf:",ts(its,j)
505 print*,"HFX:",hfx(its,j)," qfx",qfx(its,j)
506 print*,"qsfc:",qsfc(its,j)," ps:",ps(its,j)
507 print*,"wspd:",wspd(its,j)
508 print*,"znt:",znt(its,j)," delt=",delt
509 print*,"ite=",ite," kte=",kte
510 print*,"PBLH=",pblh(its,j)," KPBL=",KPBL(its,j)," xland=",xland(its,j)
511 print*," ch=",ch(its,j)
512 print*,"qke:",qke(its,1,j),qke(its,2,j),qke(its,kte,j)
513 print*,"el_pbl:",el_pbl(its,1,j),el_pbl(its,2,j),el_pbl(its,kte,j)
514 print*,"Sh3d:",Sh3d(its,1,j),sh3d(its,2,j),sh3d(its,kte,j)
515 print*,"max cf_bl:",maxval(cldfra_bl(its,:,j))
518 !print*,"In mynn wrapper, calling mynn_bl_driver"
519 CALL mynn_bl_driver( &
520 & initflag=initflag,restart=restart,cycling=cycling, &
521 & delt=delt,dz=dz(:,:,j),dx=dx,znt=znt(:,j), &
522 & u=u(:,:,j),v=v(:,:,j),w=w(:,:,j), &
523 & th=th(:,:,j),sqv3D=sqv,sqc3D=sqc, &
524 & sqi3D=sqi,sqs3D=sqs,qnc=qnc2d,qni=qni2d, &
525 & qnwfa=qnwfa2d,qnifa=qnifa2d,qnbca=qnbca2d, &
526 & ozone=ozone(:,:,j), &
527 & p=p(:,:,j),exner=exner(:,:,j),rho=rho(:,:,j), &
528 & T3D=t3d(:,:,j),xland=xland(:,j), &
529 & ts=ts(:,j),qsfc=qsfc(:,j),ps=ps(:,j), &
530 & ust=ust(:,j),ch=ch(:,j),hfx=hfx(:,j),qfx=qfx(:,j), &
531 & rmol=rmol(:,j),wspd=wspd(:,j), &
532 & uoce=uoce(:,j),voce=voce(:,j), & !input
533 & qke=QKE(:,:,j),qke_adv=qke_adv(:,:,j), & !output
534 & sh3d=Sh3d(:,:,j),sm3d=Sm3d(:,:,j), & !output
535 & nchem=nchem,kdvel=kdvel,ndvel=ndvel, & !chem/smoke
536 & Chem3d=chem,Vdep=vd, &
537 & FRP=frp_mean,EMIS_ANT_NO=emis_ant_no, &
538 & mix_chem=mix_chem,enh_mix=enh_mix, &
539 & rrfs_sd=rrfs_sd,smoke_dbg=smoke_dbg, & !end chem/smoke
540 & tsq=tsq(:,:,j),qsq=qsq(:,:,j),cov=cov(:,:,j), & !output
541 & RUBLTEN=RUBLTEN(:,:,j),RVBLTEN=RVBLTEN(:,:,j), & !output
542 & RTHBLTEN=RTHBLTEN(:,:,j),RQVBLTEN=RQVBLTEN(:,:,j), & !output
543 & RQCBLTEN=rqcblten(:,:,j),RQIBLTEN=rqiblten(:,:,j), & !output
544 & RQNCBLTEN=rqncblten(:,:,j),RQNIBLTEN=rqniblten(:,:,j), & !output
545 & RQSBLTEN=ikzero, & !there is no RQSBLTEN, so use dummy arary
546 & RQNWFABLTEN=RQNWFABLTEN(:,:,j), & !output
547 & RQNIFABLTEN=RQNIFABLTEN(:,:,j), & !output
548 & RQNBCABLTEN=RQNBCABLTEN(:,:,j), & !output
549 & dozone=rO3blten(:,:,j), & !output
550 & EXCH_H=exch_h(:,:,j),EXCH_M=exch_m(:,:,j), & !output
551 & pblh=pblh(:,j),KPBL=KPBL(:,j), & !output
552 & el_pbl=el_pbl(:,:,j), & !output
553 & dqke=dqke2d,qWT=qWT2d,qSHEAR=qSHEAR2d, & !output
554 & qBUOY=qBUOY2d,qDISS=qDISS2d, & !output
555 & qc_bl=qc_bl2d,qi_bl=qi_bl2d,cldfra_bl=cldfra_bl2d, & !output
556 & bl_mynn_tkeadvect=bl_mynn_tkeadvect, & !input parameter
557 & tke_budget=tke_budget, & !input parameter
558 & bl_mynn_cloudpdf=bl_mynn_cloudpdf, & !input parameter
559 & bl_mynn_mixlength=bl_mynn_mixlength, & !input parameter
560 & icloud_bl=icloud_bl, & !input parameter
561 & closure=bl_mynn_closure,bl_mynn_edmf=bl_mynn_edmf, & !input parameter
562 & bl_mynn_edmf_mom=bl_mynn_edmf_mom, & !input parameter
563 & bl_mynn_edmf_tke=bl_mynn_edmf_tke, & !input parameter
564 & bl_mynn_mixscalars=bl_mynn_mixscalars, & !input parameter
565 & bl_mynn_output=bl_mynn_output, & !input parameter
566 & bl_mynn_cloudmix=bl_mynn_cloudmix, & !input parameter
567 & bl_mynn_mixqt=bl_mynn_mixqt, & !input parameter
568 & edmf_a=edmf_a2d,edmf_w=edmf_w2d, & !output
569 & edmf_qt=edmf_qt2d,edmf_thl=edmf_thl2d, & !output
570 & edmf_ent=edmf_ent2d,edmf_qc=edmf_qc2d, & !output
571 & sub_thl3D=sub_thl2d,sub_sqv3D=sub_sqv2d, & !output
572 & det_thl3D=det_thl2d,det_sqv3D=det_sqv2d, & !output
573 & maxwidth=maxwidth(:,j),maxMF=maxMF(:,j), & !output
574 & ztop_plume=ztop_plume(:,j),ktop_plume=ktop_plume(:,j), & !output
575 & spp_pbl=spp_pbl,pattern_spp_pbl=pattern_spp_pbl2d, & !input
576 & RTHRATEN=rthraten(:,:,j), & !input
577 & FLAG_QI=flag_qi,FLAG_QNI=flag_qni,FLAG_QS=flag_qs, & !input
578 & FLAG_QC=flag_qc,FLAG_QNC=flag_qnc, & !input
579 & FLAG_QNWFA=FLAG_QNWFA,FLAG_QNIFA=FLAG_QNIFA, & !input
580 & FLAG_QNBCA=FLAG_QNBCA,FLAG_OZONE=flag_ozone, & !input
581 & IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde, & !input
582 & IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme, & !input
583 & ITS=its,ITE=itf,JTS=jts,JTE=jtf,KTS=kts,KTE=kte) !input
584 !print*,"In mynn wrapper, after bl_mynn_driver"
586 !- Convert spec hum to mixing ratio:
589 RQVBLTEN(i,k,j) = RQVBLTEN(i,k,j)/(1.0 - sqv(i,k))
590 RQCBLTEN(i,k,j) = RQCBLTEN(i,k,j)/(1.0 - sqv(i,k))
591 RQIBLTEN(i,k,j) = RQIBLTEN(i,k,j)/(1.0 - sqv(i,k))
594 if (.false.) then !as of now, there is no RQSBLTEN in WRF
597 RQSBLTEN(i,k,j) = RQSBLTEN(i,k,j)/(1.0 - sqv(i,k))
603 if (icloud_bl > 0) then
606 qc_bl(i,k,j) = qc_bl2d(i,k)/(1.0 - sqv(i,k))
607 qi_bl(i,k,j) = qi_bl2d(i,k)/(1.0 - sqv(i,k))
608 cldfra_bl(i,k,j) = cldfra_bl2d(i,k)
613 if (tke_budget .eq. 1) then
616 dqke(i,k,j) = dqke2d(i,k)
617 qwt(i,k,j) = qwt2d(i,k)
618 qshear(i,k,j) = qshear2d(i,k)
619 qbuoy(i,k,j) = qbuoy2d(i,k)
620 qdiss(i,k,j) = qdiss2d(i,k)
625 if (bl_mynn_output > 0) then
628 edmf_a(i,k,j) = edmf_a2d(i,k)
629 edmf_w(i,k,j) = edmf_w2d(i,k)
630 edmf_qt(i,k,j) = edmf_qt2d(i,k)
631 edmf_thl(i,k,j) = edmf_thl2d(i,k)
632 edmf_ent(i,k,j) = edmf_ent2d(i,k)
633 edmf_qc(i,k,j) = edmf_qc2d(i,k)
634 sub_thl3d(i,k,j) = sub_thl2d(i,k)
635 sub_sqv3d(i,k,j) = sub_sqv2d(i,k)
636 det_thl3d(i,k,j) = det_thl2d(i,k)
637 det_sqv3d(i,k,j) = det_sqv2d(i,k)
644 print*,"===Finished with mynn_bl_driver; output:"
645 print*,"T:",t3d(its,1,j),t3d(its,2,j),t3d(its,kte,j)
646 print*,"TH:",th(its,1,j),th(its,2,j),th(its,kte,j)
647 print*,"rho:",rho(its,1,j),rho(its,2,j),rho(its,kte,j)
648 print*,"exner:",exner(its,1,j),exner(its,2,j),exner(its,kte,j)
649 print*,"p:",p(its,1,j),p(its,2,j),p(its,kte,j)
650 print*,"dz:",dz(its,1,j),dz(its,2,j),dz(its,kte,j)
651 print*,"u:",u(its,1,j),u(its,2,j),u(its,kte,j)
652 print*,"v:",v(its,1,j),v(its,2,j),v(its,kte,j)
653 print*,"sqv:",sqv(its,1),sqv(its,2),sqv(its,kte)
654 print*,"sqc:",sqc(its,1),sqc(its,2),sqc(its,kte)
655 print*,"sqi:",sqi(its,1),sqi(its,2),sqi(its,kte)
656 print*,"rmol:",rmol(its,j)," ust:",ust(its,j)
657 print*,"dx(its,j)=",dx(its),"initflag=",initflag
658 print*,"Thetasurf:",ts(its,j)
659 print*,"HFX:",hfx(its,j)," qfx",qfx(its,j)
660 print*,"qsfc:",qsfc(its,j)," ps:",ps(its,j)
661 print*,"wspd:",wspd(its,j)
662 print*,"znt:",znt(its,j)," delt=",delt
663 print*,"im=",ite," kte=",kte
664 print*,"PBLH=",pblh(its,j)," KPBL=",KPBL(its,j)," xland=",xland(its,j)
665 print*,"ch=",ch(its,j)
666 print*,"qke:",qke(its,1,j),qke(its,2,j),qke(its,kte,j)
667 print*,"el_pbl:",el_pbl(its,1,j),el_pbl(its,2,j),el_pbl(its,kte,j)
668 print*,"Sh3d:",Sh3d(its,1,j),sh3d(its,2,j),sh3d(its,kte,j)
669 print*,"exch_h:",exch_h(its,1,j),exch_h(its,2,j),exch_h(its,kte,j)
670 print*,"exch_m:",exch_m(its,1,j),exch_m(its,2,j),exch_m(its,kte,j)
671 print*,"max cf_bl:",maxval(cldfra_bl(its,:,j))
672 print*,"max qc_bl:",maxval(qc_bl(its,:,j))
673 print*,"dtdt:",rthblten(its,1,j),rthblten(its,2,j),rthblten(its,kte,j)
674 print*,"dudt:",rublten(its,1,j),rublten(its,2,j),rublten(its,kte,j)
675 print*,"dvdt:",rvblten(its,1,j),rvblten(its,2,j),rvblten(its,kte,j)
676 print*,"dqdt:",rqvblten(its,1,j),rqvblten(its,2,j),rqvblten(its,kte,j)
677 print*,"ztop_plume:",ztop_plume(its,j)," maxmf:",maxmf(its,j)
683 !Deallocate all temporary interface arrays
684 if (bl_mynn_output > 0) then
687 deallocate(edmf_qt2d)
688 deallocate(edmf_thl2d)
689 deallocate(edmf_ent2d)
690 deallocate(edmf_qc2d)
691 deallocate(sub_thl2d)
692 deallocate(sub_sqv2d)
693 deallocate(det_thl2d)
694 deallocate(det_sqv2d)
696 if (tke_budget .eq. 1) then
703 if (icloud_bl > 0) then
706 deallocate(cldfra_bl2d)
708 if (flag_qc) deallocate(qc2d)
709 if (flag_qi) deallocate(qi2d)
710 if (flag_qs) deallocate(qs2d)
711 if (flag_qnc) deallocate(qnc2d)
712 if (flag_qni) deallocate(qni2d)
713 if (flag_qnwfa)deallocate(qnwfa2d)
714 if (flag_qnifa)deallocate(qnifa2d)
715 if (flag_qnbca)deallocate(qnbca2d)
716 if (spp_pbl > 0) then
717 deallocate(pattern_spp_pbl2d)
720 !print*,"In mynn wrapper, at end"
724 ! ==================================================================
725 SUBROUTINE moisture_check2(kte, delt, dp, exner, &
728 ! If qc < qcmin, qi < qimin, or qv < qvmin happens in any layer,
729 ! force them to be larger than minimum value by (1) condensating
730 ! water vapor into liquid or ice, and (2) by transporting water vapor
731 ! from the very lower layer.
733 ! We then update the final state variables and tendencies associated
734 ! with this correction. If any condensation happens, update theta/temperature too.
735 ! Note that (qv,qc,qi,th) are the final state variables after
736 ! applying corresponding input tendencies and corrective tendencies.
739 integer, intent(in) :: kte
740 real, intent(in) :: delt
741 real, dimension(kte), intent(in) :: dp
742 real, dimension(kte), intent(in) :: exner
743 real, dimension(kte), intent(inout) :: qv, qc, qi, th
745 real :: dqc2, dqi2, dqv2, sum, aa, dum
746 real, parameter :: qvmin1= 1e-8, & !min at k=1
747 qvmin = 1e-20, & !min above k=1
751 do k = kte, 1, -1 ! From the top to the surface
752 dqc2 = max(0.0, qcmin-qc(k)) !qc deficit (>=0)
753 dqi2 = max(0.0, qimin-qi(k)) !qi deficit (>=0)
758 qv(k) = qv(k) - dqc2 - dqi2
760 !th(k) = th(k) + xlvcp/exner(k)*dqc2 + &
761 ! xlscp/exner(k)*dqi2
763 th(k) = th(k) + xlvcp*dqc2 + &
766 !then fix qv if lending qv made it negative
768 dqv2 = max(0.0, qvmin1-qv(k)) !qv deficit (>=0)
770 qv(k) = max(qv(k),qvmin1)
773 dqv2 = max(0.0, qvmin-qv(k)) !qv deficit (>=0)
775 qv(k-1)= qv(k-1) - dqv2*dp(k)/dp(k-1)
776 qv(k) = max(qv(k),qvmin)
778 qc(k) = max(qc(k),qcmin)
779 qi(k) = max(qi(k),qimin)
782 ! Extra moisture used to satisfy 'qv(1)>=qvmin' is proportionally
783 ! extracted from all the layers that has 'qv > 2*qvmin'. This fully
784 ! preserves column moisture.
785 if( dqv2 .gt. 1.e-20 ) then
788 if( qv(k) .gt. 2.0*qvmin ) sum = sum + qv(k)*dp(k)
790 aa = dqv2*dp(1)/max(1.e-20,sum)
791 if( aa .lt. 0.5 ) then
793 if( qv(k) .gt. 2.0*qvmin ) then
799 ! For testing purposes only (not yet found in any output):
800 ! write(*,*) 'Full moisture conservation is impossible'
806 END SUBROUTINE moisture_check2
808 END SUBROUTINE mynnedmf_wrapper_run
810 !###=================================================================
812 END MODULE module_bl_mynn_wrapper