2 use ccpp_kind_types,only: kind_phys
4 use mp_wsm6,only: mp_wsm6_run,refl10cm_wsm6
5 use mp_wsm6_effectrad,only: mp_wsm6_effectRad_run
16 !=================================================================================================================
17 subroutine wsm6(th,q,qc,qr,qi,qs,qg,den,pii,p,delz, &
18 delt,g,cpd,cpv,rd,rv,t0c,ep1,ep2,qmin, &
19 xls,xlv0,xlf0,den0,denr,cliq,cice,psat, &
20 rain,rainncv,snow,snowncv,graupel,graupelncv,sr, &
21 refl_10cm,diagflag,do_radar_ref, &
22 has_reqc,has_reqi,has_reqs, &
23 re_qc_bg,re_qi_bg,re_qs_bg, &
24 re_qc_max,re_qi_max,re_qs_max, &
25 re_cloud,re_ice,re_snow, &
26 ids,ide,jds,jde,kds,kde, &
27 ims,ime,jms,jme,kms,kme, &
28 its,ite,jts,jte,kts,kte, &
31 ,wetscav_on,evapprod,rainprod &
34 !=================================================================================================================
37 logical,intent(in),optional:: diagflag
39 integer,intent(in):: ids,ide,jds,jde,kds,kde, &
40 ims,ime,jms,jme,kms,kme, &
41 its,ite,jts,jte,kts,kte
43 integer,intent(in):: has_reqc,has_reqi,has_reqs
44 integer,intent(in),optional:: do_radar_ref
46 real(kind=kind_phys),intent(in):: &
47 delt,g,rd,rv,t0c,den0,cpd,cpv,ep1,ep2,qmin,xls,xlv0,xlf0, &
50 real(kind=kind_phys),intent(in):: &
51 re_qc_bg,re_qi_bg,re_qs_bg,re_qc_max,re_qi_max,re_qs_max
53 real(kind=kind_phys),intent(in),dimension(ims:ime,kms:kme,jms:jme ):: &
60 real(kind=kind_phys),intent(inout),dimension(ims:ime,jms:jme):: &
63 real(kind=kind_phys),intent(inout),dimension(ims:ime,jms:jme),optional:: &
66 real(kind=kind_phys),intent(inout),dimension(ims:ime,jms:jme),optional:: &
69 real(kind=kind_phys),intent(inout),dimension(ims:ime,kms:kme,jms:jme):: &
78 real(kind=kind_phys),intent(inout),dimension(ims:ime,kms:kme,jms:jme):: &
83 real(kind=kind_phys),intent(inout),dimension(ims:ime,kms:kme,jms:jme),optional:: &
87 logical,intent(in):: wetscav_on
88 real(kind=kind_phys),intent(inout),dimension(ims:ime,kms:kme,jms:jme ):: &
90 real(kind=kind_phys),dimension(its:ite,kts:kte):: rainprod_hv,evapprod_hv
94 character(len=*),intent(out):: errmsg
95 integer,intent(out):: errflg
97 !local variables and arrays:
98 logical:: do_microp_re
101 real(kind=kind_phys),dimension(kts:kte):: qv1d,t1d,p1d,qr1d,qs1d,qg1d,dBZ
102 real(kind=kind_phys),dimension(kts:kte):: den1d,qc1d,qi1d,re_qc,re_qi,re_qs
104 real(kind=kind_phys),dimension(its:ite):: rainncv_hv,rain_hv,sr_hv
105 real(kind=kind_phys),dimension(its:ite):: snowncv_hv,snow_hv
106 real(kind=kind_phys),dimension(its:ite):: graupelncv_hv,graupel_hv
107 real(kind=kind_phys),dimension(its:ite,kts:kte):: t_hv,den_hv,p_hv,delz_hv
108 real(kind=kind_phys),dimension(its:ite,kts:kte):: qv_hv,qc_hv,qi_hv,qr_hv,qs_hv,qg_hv
109 real(kind=kind_phys),dimension(its:ite,kts:kte):: re_qc_hv,re_qi_hv,re_qs_hv
111 !-----------------------------------------------------------------------------------------------------------------
118 den_hv(i,k) = den(i,k,j)
120 delz_hv(i,k) = delz(i,k,j)
124 rain_hv(i) = rain(i,j)
127 t_hv(i,k) = th(i,k,j)*pii(i,k,j)
128 qv_hv(i,k) = q(i,k,j)
129 qc_hv(i,k) = qc(i,k,j)
130 qi_hv(i,k) = qi(i,k,j)
131 qr_hv(i,k) = qr(i,k,j)
132 qs_hv(i,k) = qs(i,k,j)
133 qg_hv(i,k) = qg(i,k,j)
137 if(present(snow) .and. present(snowncv)) then
139 snow_hv(i) = snow(i,j)
142 if(present(graupel) .and. present(graupelncv)) then
144 graupel_hv(i) = graupel(i,j)
148 !--- call to cloud microphysics scheme:
149 call mp_wsm6_run(t=t_hv,q=qv_hv,qc=qc_hv,qi=qi_hv,qr=qr_hv,qs=qs_hv,qg=qg_hv, &
150 den=den_hv,p=p_hv,delz=delz_hv,delt=delt,g=g,cpd=cpd,cpv=cpv, &
151 rd=rd,rv=rv,t0c=t0c,ep1=ep1,ep2=ep2,qmin=qmin,xls=xls,xlv0=xlv0, &
152 xlf0=xlf0,den0=den0,denr=denr,cliq=cliq,cice=cice,psat=psat, &
153 rain=rain_hv,rainncv=rainncv_hv,sr=sr_hv,snow=snow_hv, &
154 snowncv=snowncv_hv,graupel=graupel_hv,graupelncv=graupelncv_hv, &
155 its=its,ite=ite,kts=kts,kte=kte,errmsg=errmsg,errflg=errflg &
157 ,rainprod2d=rainprod_hv,evapprod2d=evapprod_hv &
163 rain(i,j) = rain_hv(i)
164 rainncv(i,j) = rainncv_hv(i)
168 th(i,k,j) = t_hv(i,k)/pii(i,k,j)
169 q(i,k,j) = qv_hv(i,k)
170 qc(i,k,j) = qc_hv(i,k)
171 qi(i,k,j) = qi_hv(i,k)
172 qr(i,k,j) = qr_hv(i,k)
173 qs(i,k,j) = qs_hv(i,k)
174 qg(i,k,j) = qg_hv(i,k)
178 if(present(snow) .and. present(snowncv)) then
180 snow(i,j) = snow_hv(i)
181 snowncv(i,j) = snowncv_hv(i)
184 if(present(graupel) .and. present(graupelncv)) then
186 graupel(i,j) = graupel_hv(i)
187 graupelncv(i,j) = graupelncv_hv(i)
195 rainprod(i,k,j) = rainprod_hv(i,k)
196 evapprod(i,k,j) = evapprod_hv(i,k)
209 !--- call to computation of effective radii for cloud water, cloud ice, and snow:
210 do_microp_re = .false.
211 if(has_reqc == 1 .and. has_reqi == 1 .and. has_reqs == 1) do_microp_re = .true.
215 t_hv(i,k) = th(i,k,j)*pii(i,k,j)
216 re_qc_hv(i,k) = re_cloud(i,k,j)
217 re_qi_hv(i,k) = re_ice(i,k,j)
218 re_qs_hv(i,k) = re_snow(i,k,j)
222 call mp_wsm6_effectRad_run(do_microp_re,t_hv,qc_hv,qi_hv,qs_hv,den_hv,qmin,t0c, &
223 re_qc_bg,re_qi_bg,re_qs_bg,re_qc_max,re_qi_max,re_qs_max,re_qc_hv, &
224 re_qi_hv,re_qs_hv,its,ite,kts,kte,errmsg,errflg)
228 re_cloud(i,k,j) = re_qc_hv(i,k)
229 re_ice(i,k,j) = re_qi_hv(i,k)
230 re_snow(i,k,j) = re_qs_hv(i,k)
234 if (diagflag .and. do_radar_ref == 1) then
237 t1d(k)=th(i,k,j)*pii(i,k,j)
244 call refl10cm_wsm6 (qv1d,qr1d,qs1d,qg1d,t1d,p1d,dBZ,kts,kte)
246 refl_10cm(i,k,j) = MAX(-35., dBZ(k))
255 !=================================================================================================================
256 end module module_mp_wsm6
257 !=================================================================================================================