118 subroutine gfdl_cloud_microphys_v3_run(fast_mp_consv, &
119 levs, im, rainmin, con_g, con_fvirt, con_rd, con_eps, garea, slmsk, snowd, &
120 gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, gq0_ntsw, gq0_ntgl, gq0_ntclamt, aerfld, &
121 gt0, gu0, gv0, vvl, prsl, phii, del, &
122 rain0, ice0, snow0, graupel0, prcp0, sr, oro, &
123 dtp, hydrostatic, lradar, refl_10cm, &
124 reset, effr_in, rew, rei, rer, res, reg, &
125 cplchm, pfi_lsan, pfl_lsan, con_one, con_p001, con_secinday, errmsg, errflg)
127 use machine,
only: kind_phys, kind_dyn, kind_dbl_prec
132 integer,
intent(in ) :: levs, im
133 real(kind=kind_phys),
intent(in ) :: con_g, con_fvirt, con_rd, con_eps, rainmin
134 real(kind=kind_phys),
intent(in ) :: con_one, con_p001, con_secinday
135 real(kind=kind_phys),
intent(in ),
dimension(:) :: garea, slmsk, snowd, oro
136 real(kind=kind_phys),
intent(inout),
dimension(:,:) :: gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, &
137 gq0_ntsw, gq0_ntgl, gq0_ntclamt
138 real(kind_phys),
intent(in ),
dimension(:,:,:) :: aerfld
139 real(kind=kind_phys),
intent(inout),
dimension(:,:) :: gt0, gu0, gv0
140 real(kind=kind_phys),
intent(in ),
dimension(:,:) :: vvl, prsl, del
141 real(kind=kind_phys),
intent(in ),
dimension(:,:) :: phii
145 real(kind_phys),
intent(out ),
dimension(:),
optional :: rain0
146 real(kind_phys),
intent(out ),
dimension(:),
optional :: snow0
147 real(kind_phys),
intent(out ),
dimension(:),
optional :: ice0
148 real(kind_phys),
intent(out ),
dimension(:),
optional :: graupel0
149 real(kind_phys),
intent(out ),
dimension(:) :: prcp0
150 real(kind_phys),
intent(out ),
dimension(:) :: sr
152 real(kind_phys),
intent(in) :: dtp
153 logical,
intent (in) :: hydrostatic, fast_mp_consv
155 logical,
intent (in) :: lradar
156 real(kind=kind_phys),
intent(inout),
dimension(:,:) :: refl_10cm
157 logical,
intent (in) :: reset, effr_in
158 real(kind=kind_phys),
intent(inout),
dimension(:,:),
optional :: rew, rei, rer, res, reg
159 logical,
intent (in) :: cplchm
161 real(kind=kind_phys),
intent(inout),
dimension(:,:),
optional :: pfi_lsan, pfl_lsan
163 character(len=*),
intent(out) :: errmsg
164 integer,
intent(out) :: errflg
167 integer :: iis, iie, jjs, jje, kks, kke, kbot, ktop
169 real(kind=kind_phys),
dimension(1:im,1:levs) :: delp, dz, uin, vin, pt, qv1, ql1, qi1, qr1, qs1, qg1, &
170 qa1, qnl, qni, pt_dt, qa_dt, u_dt, v_dt, w, qv_dt, ql_dt,&
171 qr_dt, qi_dt, qs_dt, qg_dt, p123, refl
172 real(kind=kind_phys),
dimension(1:im,1:levs) :: q_con, cappa
173 real(kind=kind_phys),
dimension(1:im,1,1:levs) :: pfils, pflls
174 real(kind=kind_phys),
dimension(1:im,1,1:levs) :: adj_vmr, te
175 real(kind=kind_phys),
dimension(1:im,1:levs) :: prefluxw, prefluxr, prefluxi, prefluxs, prefluxg
176 real(kind=kind_phys),
dimension(1:im) :: hs, gsize
177 real(kind=kind_dbl_prec),
dimension(1:im) :: dte
179 real(kind=kind_phys),
dimension(1:im) :: water0
180 real(kind=kind_phys) :: onebg
181 real(kind=kind_phys) :: tem
182 logical last_step, do_inline_mp
198 onebg = con_one/con_g
213 qnl(i,k) = aerfld(i,kk,11)
224 ql1(i,k) = gq0_ntcw(i,kk)
225 qr1(i,k) = gq0_ntrw(i,kk)
226 qi1(i,k) = gq0_ntiw(i,kk)
227 qs1(i,k) = gq0_ntsw(i,kk)
228 qg1(i,k) = gq0_ntgl(i,kk)
229 qa1(i,k) = gq0_ntclamt(i,kk)
231 w(i,k) = -vvl(i,kk) * (con_one+con_fvirt * gq0(i,kk)) &
232 * gt0(i,kk) / prsl(i,kk) * (con_rd*onebg)
235 delp(i,k) = del(i,kk)
236 dz(i,k) = (phii(i,kk)-phii(i,kk+1))*onebg
237 p123(i,k) = prsl(i,kk)
253 do_inline_mp = .false.
255 gsize = sqrt(garea(:))
257 call gfdl_cloud_microphys_v3_mod_driver( qv1, ql1, qr1, qi1, qs1, qg1, qa1, qnl, qni, pt, w,&
258 uin, vin, dz, delp, gsize, dtp, hs, water0, rain0, &
259 ice0, snow0, graupel0, hydrostatic, iis, iie, kks, kke, q_con, cappa, &
260 fast_mp_consv, adj_vmr, te, dte, prefluxw, prefluxr, prefluxi, prefluxs, &
261 prefluxg, last_step, do_inline_mp )
262 tem = dtp*con_p001/con_secinday
270 if(water0(i)*tem < rainmin)
then
273 if(rain0(i)*tem < rainmin)
then
276 if(ice0(i)*tem < rainmin)
then
279 if(snow0(i)*tem < rainmin)
then
282 if(graupel0(i)*tem < rainmin)
then
290 prcp0(i) = (rain0(i)+snow0(i)+ice0(i)+graupel0(i)) * tem
291 if ( prcp0(i) > rainmin )
then
292 sr(i) = (snow0(i) + ice0(i) + graupel0(i)) &
293 / (rain0(i) + snow0(i) + ice0(i) + graupel0(i))
304 graupel0 = graupel0*tem
311 gq0_ntcw(i,k) = ql1(i,kk)
312 gq0_ntrw(i,k) = qr1(i,kk)
313 gq0_ntiw(i,k) = qi1(i,kk)
314 gq0_ntsw(i,k) = qs1(i,kk)
315 gq0_ntgl(i,k) = qg1(i,kk)
316 gq0_ntclamt(i,k) = qa1(i,kk)
320 refl_10cm(i,k) = refl(i,kk)
329 pfi_lsan(i,k) = prefluxi(i,kk) + prefluxs(i,kk) + prefluxg(i,kk)
330 pfl_lsan(i,k) = prefluxr(i,kk)
336 call cld_eff_rad (1, im, 1, levs, slmsk(1:im), &
337 prsl(1:im,1:levs), del(1:im,1:levs), &
338 gt0(1:im,1:levs), gq0(1:im,1:levs), &
339 gq0_ntcw(1:im,1:levs), gq0_ntiw(1:im,1:levs), &
340 gq0_ntrw(1:im,1:levs), gq0_ntsw(1:im,1:levs), &
341 gq0_ntgl(1:im,1:levs), gq0_ntclamt(1:im,1:levs), &
342 rew(1:im,1:levs), rei(1:im,1:levs), rer(1:im,1:levs),&
343 res(1:im,1:levs), reg(1:im,1:levs),snowd(1:im))
347 call rad_ref (1, im, 1, 1, qv1(1:im,1:levs), qr1(1:im,1:levs), &
348 qs1(1:im,1:levs),qg1(1:im,1:levs),pt(1:im,1:levs), &
349 delp(1:im,1:levs), dz(1:im,1:levs), refl(1:im,1:levs), levs, hydrostatic, &
355 refl_10cm(i,k) = max(-35.,refl(i,kk))