CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
cu_gf_driver.F90
1
3
6
7 ! DH* TODO: replace constants with arguments to cu_gf_driver_run
8 !use physcons , g => con_g, cp => con_cp, xlv => con_hvap, r_v => con_rv
9 use machine , only: kind_phys
11 use cu_gf_sh , only: cu_gf_sh_run
12
13 implicit none
14
15 private
16
18
19contains
20
31 subroutine cu_gf_driver_init(imfshalcnv, imfshalcnv_gf, imfdeepcnv, &
32 imfdeepcnv_gf,mpirank, mpiroot, errmsg, errflg)
33
34 implicit none
35
36 integer, intent(in) :: imfshalcnv, imfshalcnv_gf
37 integer, intent(in) :: imfdeepcnv, imfdeepcnv_gf
38 integer, intent(in) :: mpirank
39 integer, intent(in) :: mpiroot
40 character(len=*), intent( out) :: errmsg
41 integer, intent( out) :: errflg
42
43 ! initialize ccpp error handling variables
44 errmsg = ''
45 errflg = 0
46
47 end subroutine cu_gf_driver_init
48
49!
50! t2di is temp after advection, but before physics
51! t = current temp (t2di + physics up to now)
52!===================
53
59 subroutine cu_gf_driver_run(ntracer,garea,im,km,dt,flag_init,flag_restart,&
60 cactiv,cactiv_m,g,cp,xlv,r_v,forcet,forceqv_spechum,phil,raincv, &
61 qv_spechum,t,cld1d,us,vs,t2di,w,qv2di_spechum,p2di,psuri, &
62 hbot,htop,kcnv,xland,hfx2,qfx2,aod_gf,cliw,clcw, &
63 pbl,ud_mf,dd_mf,dt_mf,cnvw_moist,cnvc,imfshalcnv, &
64 flag_for_scnv_generic_tend,flag_for_dcnv_generic_tend, &
65 dtend,dtidx,ntqv,ntiw,ntcw,index_of_temperature,index_of_x_wind, &
66 index_of_y_wind,index_of_process_scnv,index_of_process_dcnv, &
67 fhour,fh_dfi_radar,ix_dfi_radar,num_dfi_radar,cap_suppress, &
68 dfi_radar_max_intervals,ldiag3d,qci_conv,do_cap_suppress, &
69 maxupmf,maxMF,do_mynnedmf,ichoice_in,ichoicem_in,ichoice_s_in, &
70 spp_cu_deep,spp_wts_cu_deep,nchem,chem3d,fscav,wetdpc_deep, &
71 do_smoke_transport,kdt,errmsg,errflg)
72!-------------------------------------------------------------
73 implicit none
74 integer, parameter :: maxiens=1
75 integer, parameter :: maxens=1
76 integer, parameter :: maxens2=1
77 integer, parameter :: maxens3=16
78 integer, parameter :: ensdim=16
79 integer :: imid_gf=1 ! gf congest conv.
80 integer, parameter :: ideep=1
81 integer :: ichoice=0 ! 0 2 5 13 8
82 integer :: ichoicem=13 ! 0 2 5 13
83 integer :: ichoice_s=3 ! 0 1 2 3
84 integer, intent(in) :: spp_cu_deep ! flag for using SPP perturbations
85 real(kind_phys), dimension(:,:), intent(in),optional :: &
86 & spp_wts_cu_deep
87 real(kind=kind_phys) :: spp_wts_cu_deep_tmp
88
89 logical, intent(in) :: do_cap_suppress, do_smoke_transport
90 real(kind=kind_phys), parameter :: aodc0=0.14
91 real(kind=kind_phys), parameter :: aodreturn=30.
92 real(kind=kind_phys) :: dts,fpi,fp
93 integer, parameter :: dicycle=0 ! diurnal cycle flag
94 integer, parameter :: dicycle_m=0 !- diurnal cycle flag
95 integer :: ishallow_g3 ! depend on imfshalcnv
96!-------------------------------------------------------------
97 integer :: its,ite, jts,jte, kts,kte
98 integer, intent(in ) :: im,km,ntracer,nchem,kdt
99 integer, intent(in ) :: ichoice_in,ichoicem_in,ichoice_s_in
100 logical, intent(in ) :: flag_init, flag_restart, do_mynnedmf
101 logical, intent(in ) :: flag_for_scnv_generic_tend,flag_for_dcnv_generic_tend
102 real (kind=kind_phys), intent(in) :: g,cp,xlv,r_v
103 logical, intent(in ) :: ldiag3d
104
105 real(kind=kind_phys), intent(inout), optional :: dtend(:,:,:)
106!$acc declare copy(dtend)
107 integer, intent(in) :: dtidx(:,:), &
108 index_of_x_wind, index_of_y_wind, index_of_temperature, &
109 index_of_process_scnv, index_of_process_dcnv, ntqv, ntcw, ntiw
110!$acc declare copyin(dtidx)
111 real(kind=kind_phys), dimension( : , : ), intent(in ), optional :: forcet,forceqv_spechum
112 real(kind=kind_phys), dimension( : , : ), intent(in ) :: w,phil
113 real(kind=kind_phys), dimension( : , : ), intent(inout ) :: t,us,vs
114 real(kind=kind_phys), dimension( : , : ), intent(inout ), optional :: qci_conv
115 real(kind=kind_phys), dimension( : , : ), intent(out ) :: cnvw_moist,cnvc
116 real(kind=kind_phys), dimension( : , : ), intent(inout ) :: cliw, clcw
117!$acc declare copyin(forcet,forceqv_spechum,w,phil)
118!$acc declare copy(t,us,vs,qci_conv,cliw, clcw)
119!$acc declare copyout(cnvw_moist,cnvc)
120
121 real(kind=kind_phys), allocatable :: clcw_save(:,:), cliw_save(:,:)
122
123 integer, intent(in) :: dfi_radar_max_intervals
124 real(kind=kind_phys), intent(in) :: fhour, fh_dfi_radar(:)
125 integer, intent(in) :: num_dfi_radar, ix_dfi_radar(:)
126 real(kind=kind_phys), intent(in), optional :: cap_suppress(:,:)
127!$acc declare copyin(fh_dfi_radar,ix_dfi_radar,cap_suppress)
128
129 integer, dimension (:), intent(out) :: hbot,htop,kcnv
130 integer, dimension (:), intent(in) :: xland
131 real(kind=kind_phys), dimension (:), intent(in) :: pbl
132 real(kind=kind_phys), dimension (:), intent(in), optional :: maxmf
133!$acc declare copyout(hbot,htop,kcnv)
134!$acc declare copyin(xland,pbl)
135 integer, dimension (im) :: tropics
136!$acc declare create(tropics)
137! ruc variable
138 real(kind=kind_phys), dimension (:), intent(in) :: hfx2,qfx2,psuri
139 real(kind=kind_phys), dimension (:,:), intent(out) :: dd_mf,dt_mf
140 real(kind=kind_phys), dimension (:,:), intent(out), optional :: ud_mf
141 real(kind=kind_phys), dimension (:), intent(out) :: raincv,cld1d
142 real(kind=kind_phys), dimension (:), intent(out), optional :: maxupmf
143 real(kind=kind_phys), dimension (:,:), intent(in) :: t2di,p2di
144!$acc declare copyin(hfx2,qfx2,psuri,t2di,p2di)
145!$acc declare copyout(ud_mf,dd_mf,dt_mf,raincv,cld1d)
146 ! Specific humidity from FV3
147 real(kind=kind_phys), dimension (:,:), intent(in) :: qv2di_spechum
148 real(kind=kind_phys), dimension (:,:), intent(inout) :: qv_spechum
149 real(kind=kind_phys), dimension (:), intent(inout), optional :: aod_gf
150!$acc declare copyin(qv2di_spechum) copy(qv_spechum,aod_gf)
151 ! Local water vapor mixing ratios and cloud water mixing ratios
152 real(kind=kind_phys), dimension (im,km) :: qv2di, qv, forceqv, cnvw
153!$acc declare create(qv2di, qv, forceqv, cnvw)
154 !
155 real(kind=kind_phys), dimension(:),intent(in) :: garea
156!$acc declare copyin(garea)
157 real(kind=kind_phys), intent(in ) :: dt
158
159 integer, intent(in ) :: imfshalcnv
160 integer, dimension(:), intent(inout), optional :: cactiv,cactiv_m
161 real(kind_phys), dimension(:), intent(in) :: fscav
162!$acc declare copyin(fscav)
163 real(kind_phys), dimension(:,:,:), intent(inout), optional :: chem3d
164 real(kind_phys), dimension(:,:), intent(inout), optional :: wetdpc_deep
165!$acc declare copy(cactiv,cactiv_m,chem3d,wetdpc_deep)
166
167 character(len=*), intent(out) :: errmsg
168 integer, intent(out) :: errflg
169
170! local variables
171 integer, dimension(im) :: k22_shallow,kbcon_shallow,ktop_shallow
172 real(kind=kind_phys), dimension (im) :: rand_mom,rand_vmas
173 real(kind=kind_phys), dimension (im,4) :: rand_clos
174 real(kind=kind_phys), dimension (im,km,11) :: gdc,gdc2
175 real(kind=kind_phys), dimension (im) :: ht
176 real(kind=kind_phys), dimension (im) :: ccn_gf,ccn_m
177 real(kind=kind_phys) :: ccnclean
178 real(kind=kind_phys), dimension (im) :: dx
179 real(kind=kind_phys), dimension (im) :: frhm,frhd
180 real(kind=kind_phys), dimension (im,km) :: outt,outq,outqc,phh,subm,cupclw,cupclws
181 real(kind=kind_phys), dimension (im,km) :: dhdt,zu,zus,zd,phf,zum,zdm,outum,outvm
182 real(kind=kind_phys), dimension (im,km) :: outts,outqs,outqcs,outu,outv,outus,outvs
183 real(kind=kind_phys), dimension (im,km) :: outtm,outqm,outqcm,submm,cupclwm
184 real(kind=kind_phys), dimension (im,km) :: cnvwt,cnvwts,cnvwtm
185 real(kind=kind_phys), dimension (im,km) :: hco,hcdo,zdo,zdd,hcom,hcdom,zdom
186 real(kind=kind_phys), dimension (km) :: zh
187 real(kind=kind_phys), dimension (im) :: tau_ecmwf,edt,edtm,edtd,ter11,aa0,xlandi
188 real(kind=kind_phys), dimension (im) :: pret,prets,pretm,hexec
189 real(kind=kind_phys), dimension (im,10) :: forcing,forcing2
190 real(kind=kind_phys), dimension (im,nchem) :: wetdpc_mid
191
192 integer, dimension (im) :: kbcon, ktop,ierr,ierrs,ierrm,kpbli
193 integer, dimension (im) :: k22s,kbcons,ktops,k22,jmin,jminm
194 integer, dimension (im) :: kbconm,ktopm,k22m
195!$acc declare create(k22_shallow,kbcon_shallow,ktop_shallow,rand_mom,rand_vmas, &
196!$acc rand_clos,gdc,gdc2,ht,ccn_gf,ccn_m,dx,frhm,frhd,wetdpc_mid, &
197!$acc outt,outq,outqc,phh,subm,cupclw,cupclws, &
198!$acc dhdt,zu,zus,zd,phf,zum,zdm,outum,outvm, &
199!$acc outts,outqs,outqcs,outu,outv,outus,outvs, &
200!$acc outtm,outqm,outqcm,submm,cupclwm, &
201!$acc cnvwt,cnvwts,cnvwtm,hco,hcdo,zdo,zdd,hcom,hcdom,zdom, &
202!$acc tau_ecmwf,edt,edtm,edtd,ter11,aa0,xlandi, &
203!$acc pret,prets,pretm,hexec,forcing,forcing2,wetdpc_mid, &
204!$acc kbcon, ktop,ierr,ierrs,ierrm,kpbli, &
205!$acc k22s,kbcons,ktops,k22,jmin,jminm,kbconm,ktopm,k22m)
206
207 integer :: iens,ibeg,iend,jbeg,jend,n
208 integer :: ibegh,iendh,jbegh,jendh
209 integer :: ibegc,iendc,jbegc,jendc,kstop
210 real(kind=kind_phys), dimension(im,km) :: rho_dryar
211!$acc declare create(rho_dryar)
212 real(kind=kind_phys) :: pten,pqen,paph,zrho,pahfs,pqhfl,zkhvfl,pgeoh
213 integer, parameter :: ipn = 0
214
215!
216! basic environmental input includes moisture convergence (mconv)
217! omega (omeg), windspeed (us,vs), and a flag (ierr) to turn off
218! convection for this call only and at that particular gridpoint
219!
220 real(kind=kind_phys), dimension (im,km) :: qcheck,zo,t2d,q2d,po,p2d,rhoi,clw_ten
221 real(kind=kind_phys), dimension (im,km) :: tn,qo,tshall,qshall,dz8w,omeg
222 real(kind=kind_phys), dimension (im) :: z1,psur,cuten,cutens,cutenm
223 real(kind=kind_phys), dimension (im) :: umean,vmean,pmean
224 real(kind=kind_phys), dimension (im) :: xmbs,xmbs2,xmb,xmbm,xmb_dumm,mconv
225!$acc declare create(qcheck,zo,t2d,q2d,po,p2d,rhoi,clw_ten,tn,qo,tshall,qshall,dz8w,omeg, &
226!$acc z1,psur,cuten,cutens,cutenm,umean,vmean,pmean, &
227!$acc xmbs,xmbs2,xmb,xmbm,xmb_dumm,mconv)
228
229 integer :: i,j,k,icldck,ipr,jpr,jpr_deep,ipr_deep,uidx,vidx,tidx,qidx
230 integer :: itf,jtf,ktf,iss,jss,nbegin,nend,cliw_idx,clcw_idx
231 integer :: high_resolution
232 real(kind=kind_phys) :: clwtot,clwtot1,excess,tcrit,tscl_kf,dp,dq,sub_spread,subcenter
233 real(kind=kind_phys) :: dsubclw,dsubclws,dsubclwm,dtime_max,ztm,ztq,hfm,qfm,rkbcon,rktop
234 real(kind=kind_phys), dimension(km) :: massflx,trcflx_in1,clw_in1,po_cup
235! real(kind=kind_phys), dimension(km) :: trcflx_in2,clw_in2,clw_ten2
236 real(kind=kind_phys), dimension (im) :: flux_tun,tun_rad_mid,tun_rad_shall,tun_rad_deep
237!$acc declare create(flux_tun,tun_rad_mid,tun_rad_shall,tun_rad_deep)
238 character*50 :: ierrc(im),ierrcm(im)
239 character*50 :: ierrcs(im)
240! ruc variable
241! hfx2 -- sensible heat flux (k m/s), positive upward from sfc
242! qfx2 -- latent heat flux (kg/kg m/s), positive upward from sfc
243! gf needs them in w/m2. define hfx and qfx after simple unit conversion
244 real(kind=kind_phys), dimension (im) :: hfx,qfx
245!$acc declare create(hfx,qfx)
246 real(kind=kind_phys) tem,tem1,tf,tcr,tcrf,psum
247 real(kind=kind_phys) :: cliw_shal,clcw_shal,tem_shal, cliw_both, weight_sum
248 real(kind=kind_phys) :: cliw_deep,clcw_deep,tem_deep, clcw_both
249 integer :: cliw_deep_idx, clcw_deep_idx, cliw_shal_idx, clcw_shal_idx
250
251 real(kind=kind_phys) :: cap_suppress_j(im)
252!$acc declare create(cap_suppress_j)
253 integer :: itime, do_cap_suppress_here
254 logical :: exit_func
255
256 !parameter (tf=243.16, tcr=270.16, tcrf=1.0/(tcr-tf)) ! FV3 original
257 !parameter (tf=263.16, tcr=273.16, tcrf=1.0/(tcr-tf))
258 !parameter (tf=233.16, tcr=263.16, tcrf=1.0/(tcr-tf))
259 parameter(tf=258.16, tcr=273.16, tcrf=1.0/(tcr-tf)) ! as fim, HCB tuning
260 ! initialize ccpp error handling variables
261 errmsg = ''
262 errflg = 0
263
264 ichoice = ichoice_in
265 ichoicem = ichoicem_in
266 ichoice_s = ichoice_s_in
267 if(do_cap_suppress) then
268!$acc serial
269 do itime=1,num_dfi_radar
270 if(ix_dfi_radar(itime)<1) cycle
271 if(fhour<fh_dfi_radar(itime)) cycle
272 if(fhour>=fh_dfi_radar(itime+1)) cycle
273 exit
274 enddo
275!$acc end serial
276 endif
277 if(do_cap_suppress .and. itime<=num_dfi_radar) then
278 do_cap_suppress_here = 1
279!$acc kernels
280 cap_suppress_j(:) = cap_suppress(:,itime)
281!$acc end kernels
282 else
283 do_cap_suppress_here = 0
284!$acc kernels
285 cap_suppress_j(:) = 0
286!$acc end kernels
287 endif
288
289 if(ldiag3d) then
290 if(flag_for_dcnv_generic_tend) then
291 cliw_deep_idx=0
292 clcw_deep_idx=0
293 else
294 cliw_deep_idx=dtidx(100+ntiw,index_of_process_dcnv)
295 clcw_deep_idx=dtidx(100+ntcw,index_of_process_dcnv)
296 endif
297 if(flag_for_scnv_generic_tend) then
298 cliw_shal_idx=0
299 clcw_shal_idx=0
300 else
301 cliw_shal_idx=dtidx(100+ntiw,index_of_process_scnv)
302 clcw_shal_idx=dtidx(100+ntcw,index_of_process_scnv)
303 endif
304 if(cliw_deep_idx>=1 .or. clcw_deep_idx>=1 .or. &
305 cliw_shal_idx>=1 .or. clcw_shal_idx>=1) then
306 allocate(clcw_save(im,km), cliw_save(im,km))
307!$acc enter data create(clcw_save,cliw_save)
308!$acc kernels
309 clcw_save(:,:)=clcw(:,:)
310 cliw_save(:,:)=cliw(:,:)
311!$acc end kernels
312 endif
313 endif
314
315!
316! Scale specific humidity to dry mixing ratio
317!
318!$acc kernels
319 ! state in before physics
320 qv2di = qv2di_spechum/(1.0_kind_phys-qv2di_spechum)
321 ! forcing by dynamics, based on state in
322 forceqv = forceqv_spechum/(1.0_kind_phys-qv2di_spechum)
323 ! current state (updated by preceeding physics)
324 qv = qv_spechum/(1.0_kind_phys-qv_spechum)
325!
326!
327! these should be coming in from outside
328!
329! cactiv(:) = 0
330 if (spp_cu_deep == 0) then
331 rand_mom(:) = 0.
332 rand_vmas(:) = 0.
333 rand_clos(:,:) = 0.
334 else
335 do i=1,im
336 spp_wts_cu_deep_tmp=min(max(-1.0_kind_phys, spp_wts_cu_deep(i,1)),1.0_kind_phys)
337 rand_mom(i) = spp_wts_cu_deep_tmp
338 rand_vmas(i) = spp_wts_cu_deep_tmp
339 rand_clos(i,:) = spp_wts_cu_deep_tmp
340 end do
341 end if
342!$acc end kernels
343!
344 its=1
345 ite=im
346 itf=ite
347 jts=1
348 jte=1
349 jtf=jte
350 kts=1
351 kte=km
352 ktf=kte-1
353!$acc kernels
354!
355 tropics(:)=0
356!
358!
359 tun_rad_shall(:)=.01
360 tun_rad_mid(:)=.3 !.02
361 tun_rad_deep(:)=.3 !.065
362 edt(:)=0.
363 edtm(:)=0.
364 edtd(:)=0.
365 zdd(:,:)=0.
366 flux_tun(:)=5.
367! dx for scale awareness
368!$acc end kernels
369
370 if (imfshalcnv == 3) then
371 ishallow_g3 = 1
372 else
373 ishallow_g3 = 0
374 end if
375 high_resolution=0
376 subcenter=0.
377 iens=1
378!
379! these can be set for debugging
380!
381 ipr=0
382 jpr=0
383 ipr_deep=0
384 jpr_deep= 0 !53322 ! 528196 !0 ! 1136 !0 !421755 !3536
385!
386!
387 ibeg=its
388 iend=ite
389 tcrit=258.
390
391 ztm=0.
392 ztq=0.
393 hfm=0.
394 qfm=0.
395!$acc kernels
396 ud_mf(:,:) =0.
397 dd_mf(:,:) =0.
398 dt_mf(:,:) =0.
399 tau_ecmwf(:)=0.
400!$acc end kernels
401!
402 j=1
403!$acc kernels
404 ht(:)=phil(:,1)/g
405!$acc loop private(zh)
406 do i=its,ite
407 cld1d(i)=0.
408 zo(i,:)=phil(i,:)/g
409 dz8w(i,1)=zo(i,2)-zo(i,1)
410 zh(1)=0.
411 kpbli(i)=2
412 do k=kts+1,ktf
413 dz8w(i,k)=zo(i,k+1)-zo(i,k)
414 enddo
415!$acc loop seq
416 do k=kts+1,ktf
417 zh(k)=zh(k-1)+dz8w(i,k-1)
418 if(zh(k).gt.pbl(i))then
419 kpbli(i)=max(2,k)
420 exit
421 endif
422 enddo
423 enddo
424!$acc end kernels
425
426!$acc kernels
427 do i= its,itf
428 forcing(i,:)=0.
429 forcing2(i,:)=0.
430 ccn_gf(i) = 0.
431 ccn_m(i) = 0.
432
433 ! set aod and ccn
434 if (flag_init .and. .not.flag_restart) then
435 aod_gf(i)=aodc0
436 else
437 if((cactiv(i).eq.0) .and. (cactiv_m(i).eq.0))then
438 if(aodc0>aod_gf(i)) aod_gf(i)=aod_gf(i)+((aodc0-aod_gf(i))*(dt/(aodreturn*60)))
439 if(aod_gf(i)>aodc0) aod_gf(i)=aodc0
440 endif
441 endif
442
443 ccn_gf(i)=max(5., (aod_gf(i)/0.0027)**(1/0.640))
444 ccn_m(i)=ccn_gf(i)
445
446 ccnclean=max(5., (aodc0/0.0027)**(1/0.640))
447
448 hbot(i) =kte
449 htop(i) =kts
450 raincv(i)=0.
451 xlandi(i)=real(xland(i))
452! if(abs(xlandi(i)-1.).le.1.e-3) tun_rad_shall(i)=.15
453! if(abs(xlandi(i)-1.).le.1.e-3) flux_tun(i)=1.5
454 enddo
455 do i= its,itf
456 mconv(i)=0.
457 enddo
458 do k=kts,kte
459 do i= its,itf
460 omeg(i,k)=0.
461 zu(i,k)=0.
462 zum(i,k)=0.
463 zus(i,k)=0.
464 zd(i,k)=0.
465 zdm(i,k)=0.
466 enddo
467 enddo
468
469 psur(:)=0.01*psuri(:)
470 do i=its,itf
471 ter11(i)=max(0.,ht(i))
472 enddo
473 do k=kts,kte
474 do i=its,ite
475 cnvw(i,k)=0.
476 cnvc(i,k)=0.
477 gdc(i,k,1)=0.
478 gdc(i,k,2)=0.
479 gdc(i,k,3)=0.
480 gdc(i,k,4)=0.
481 gdc(i,k,7)=0.
482 gdc(i,k,8)=0.
483 gdc(i,k,9)=0.
484 gdc(i,k,10)=0.
485 gdc2(i,k,1)=0.
486 enddo
487 enddo
488 ierr(:)=0
489 ierrm(:)=0
490 ierrs(:)=0
491 cuten(:)=0.
492 cutenm(:)=0.
493 cutens(:)=0.
494!$acc end kernels
495 ierrc(:)=" "
496!$acc kernels
497
498
499 kbcon(:)=0
500 kbcons(:)=0
501 kbconm(:)=0
502
503 ktop(:)=0
504 ktops(:)=0
505 ktopm(:)=0
506
507 xmb(:)=0.
508 xmb_dumm(:)=0.
509 xmbm(:)=0.
510 xmbs(:)=0.
511 xmbs2(:)=0.
512
513 k22s(:)=0
514 k22m(:)=0
515 k22(:)=0
516
517 jmin(:)=0
518 jminm(:)=0
519
520 pret(:)=0.
521 prets(:)=0.
522 pretm(:)=0.
523
524 umean(:)=0.
525 vmean(:)=0.
526 pmean(:)=0.
527
528 cupclw(:,:)=0.
529 cupclwm(:,:)=0.
530 cupclws(:,:)=0.
531
532 cnvwt(:,:)=0.
533 cnvwts(:,:)=0.
534 cnvwtm(:,:)=0.
535
536 hco(:,:)=0.
537 hcom(:,:)=0.
538 hcdo(:,:)=0.
539 hcdom(:,:)=0.
540
541 outt(:,:)=0.
542 outts(:,:)=0.
543 outtm(:,:)=0.
544
545 outu(:,:)=0.
546 outus(:,:)=0.
547 outum(:,:)=0.
548
549 outv(:,:)=0.
550 outvs(:,:)=0.
551 outvm(:,:)=0.
552
553 outq(:,:)=0.
554 outqs(:,:)=0.
555 outqm(:,:)=0.
556
557 outqc(:,:)=0.
558 outqcs(:,:)=0.
559 outqcm(:,:)=0.
560
561 subm(:,:)=0.
562 dhdt(:,:)=0.
563
564 frhm(:)=0.
565 frhd(:)=0.
566
567 do k=kts,ktf
568 do i=its,itf
569 p2d(i,k)=0.01*p2di(i,k)
570 po(i,k)=p2d(i,k) !*.01
571 rhoi(i,k) = 100.*p2d(i,k)/(287.04*(t2di(i,k)*(1.+0.608*qv2di(i,k))))
572 qcheck(i,k)=qv(i,k)
573 tn(i,k)=t(i,k)!+forcet(i,k)*dt
574 qo(i,k)=max(1.e-16,qv(i,k))!+forceqv(i,k)*dt
575 t2d(i,k)=t2di(i,k)-forcet(i,k)*dt
576 q2d(i,k)=max(1.e-16,qv2di(i,k)-forceqv(i,k)*dt)
577 if(qo(i,k).lt.1.e-16)qo(i,k)=1.e-16
578 tshall(i,k)=t2d(i,k)
579 qshall(i,k)=q2d(i,k)
580 enddo
581 enddo
582!$acc end kernels
583123 format(1x,i2,1x,2(1x,f8.0),1x,2(1x,f8.3),3(1x,e13.5))
584!$acc kernels
585 do i=its,itf
586 do k=kts,kpbli(i)
587 tshall(i,k)=t(i,k)
588 qshall(i,k)=max(1.e-16,qv(i,k))
589 enddo
590 enddo
591!
592! converting hfx2 and qfx2 to w/m2
593! hfx=cp*rho*hfx2
594! qfx=xlv*qfx2
595 do i=its,itf
596 hfx(i)=hfx2(i)*cp*rhoi(i,1)
597 qfx(i)=qfx2(i)*xlv*rhoi(i,1)
598 dx(i) = sqrt(garea(i))
599 enddo
600
601 do i=its,itf
602 do k=kts,kpbli(i)
603 tn(i,k)=t(i,k)
604 qo(i,k)=max(1.e-16,qv(i,k))
605 enddo
606 enddo
607 nbegin=0
608 nend=0
609 do i=its,itf
610 do k=kts,kpbli(i)
611 dhdt(i,k)=cp*(forcet(i,k)+(t(i,k)-t2di(i,k))/dt) + &
612 xlv*(forceqv(i,k)+(qv(i,k)-qv2di(i,k))/dt)
613! tshall(i,k)=t(i,k)
614! qshall(i,k)=qv(i,k)
615 enddo
616 enddo
617!$acc loop collapse(2) independent private(dp)
618 do k= kts+1,ktf-1
619 do i = its,itf
620 if((p2d(i,1)-p2d(i,k)).gt.150.and.p2d(i,k).gt.300)then
621 dp=-.5*(p2d(i,k+1)-p2d(i,k-1))
622!$acc atomic
623 umean(i)=umean(i)+us(i,k)*dp
624!$acc atomic
625 vmean(i)=vmean(i)+vs(i,k)*dp
626!$acc atomic
627 pmean(i)=pmean(i)+dp
628 endif
629 enddo
630 enddo
631 do i = its,itf
632 psum=0.
633 do k=kts,ktf-3
634 if (clcw(i,k) .gt. -999.0 .and. clcw(i,k+1) .gt. -999.0 )then
635 dp=(p2d(i,k)-p2d(i,k+1))
636 psum=psum+dp
637 clwtot = cliw(i,k) + clcw(i,k)
638 if(clwtot.lt.1.e-32)clwtot=0.
639 forcing(i,7)=forcing(i,7)+clwtot*dp
640 endif
641 enddo
642 if(psum.gt.0)forcing(i,7)=forcing(i,7)/psum
643 forcing2(i,7)=forcing(i,7)
644 enddo
645 do k=kts,ktf-1
646 do i = its,itf
647 omeg(i,k)= w(i,k) !-g*rhoi(i,k)*w(i,k)
648 enddo
649 enddo
650 do i = its,itf
651 if(mconv(i).lt.0.)mconv(i)=0.
652 if((dx(i)<6500.).and.do_mynnedmf.and.(maxmf(i).gt.0.))ierr(i)=555
653 enddo
654!$acc end kernels
655 if (dx(its)<6500.) then
656 imid_gf=0
657 endif
658!
659!---- call cumulus parameterization
660!
661 if(ishallow_g3.eq.1)then
662
663!$acc kernels
664 do i=its,ite
665 ierrs(i)=0
666 ierrm(i)=0
667 enddo
668!$acc end kernels
669!
671!
672 call cu_gf_sh_run (us,vs, &
673! input variables, must be supplied
674 zo,t2d,q2d,ter11,tshall,qshall,p2d,psur,dhdt,kpbli, &
675 rhoi,hfx,qfx,xlandi,ichoice_s,tcrit,dt, &
676! input variables. ierr should be initialized to zero or larger than zero for
677! turning off shallow convection for grid points
678 zus,xmbs,kbcons,ktops,k22s,ierrs,ierrcs, &
679! output tendencies
680 outts,outqs,outqcs,outus,outvs,cnvwt,prets,cupclws, &
681! dimesnional variables
682 itf,ktf,its,ite, kts,kte,ipr,tropics)
683
684!$acc kernels
685 do i=its,itf
686 if(xmbs(i).gt.0.)then
687 cutens(i)=1.
688 if (dx(i)<6500.) then
689 ierrm(i)=555
690 ierr(i)=555
691 endif
692 endif
693 enddo
694!$acc end kernels
696 call neg_check('shallow',ipn,dt,qcheck,outqs,outts,outus,outvs, &
697 outqcs,prets,its,ite,kts,kte,itf,ktf,ktops)
698 endif
699
700 ipr=0
701 jpr_deep=0 !340765
703 if(imid_gf == 1)then
704 call cu_gf_deep_run( &
705 itf,ktf,its,ite, kts,kte &
706 ,dicycle_m &
707 ,ichoicem &
708 ,ipr &
709 ,ccn_m &
710 ,ccnclean &
711 ,dt &
712 ,imid_gf &
713 ,kpbli &
714 ,dhdt &
715 ,xlandi &
716 ,zo &
717 ,forcing &
718 ,t2d &
719 ,q2d &
720 ,ter11 &
721 ,tshall &
722 ,qshall &
723 ,p2d &
724 ,psur &
725 ,us &
726 ,vs &
727 ,rhoi &
728 ,hfx &
729 ,qfx &
730 ,dx &
731 ,mconv &
732 ,omeg &
733 ,cactiv_m &
734 ,cnvwtm &
735 ,zum &
736 ,zdm & ! hli
737 ,zdd &
738 ,edtm &
739 ,edtd & ! hli
740 ,xmbm &
741 ,xmb_dumm &
742 ,xmbs &
743 ,pretm &
744 ,outum &
745 ,outvm &
746 ,outtm &
747 ,outqm &
748 ,outqcm &
749 ,kbconm &
750 ,ktopm &
751 ,cupclwm &
752 ,frhm &
753 ,ierrm &
754 ,ierrcm &
755 ,nchem &
756 ,fscav &
757 ,chem3d &
758 ,wetdpc_mid &
759 ,do_smoke_transport &
760! the following should be set to zero if not available
761 ,rand_mom & ! for stochastics mom, if temporal and spatial patterns exist
762 ,rand_vmas & ! for stochastics vertmass, if temporal and spatial patterns exist
763 ,rand_clos & ! for stochastics closures, if temporal and spatial patterns exist
764 ,spp_cu_deep & ! flag to what you want perturbed
765 ! 1 = momentum transport
766 ! 2 = normalized vertical mass flux profile
767 ! 3 = closures
768 ! more is possible, talk to developer or
769 ! implement yourself. pattern is expected to be
770 ! betwee -1 and +1
771 ,do_cap_suppress_here,cap_suppress_j &
772 ,k22m &
773 ,jminm,kdt,tropics)
774!$acc kernels
775 do i=its,itf
776 do k=kts,ktf
777 qcheck(i,k)=qv(i,k) +outqs(i,k)*dt
778 enddo
779 enddo
780!$acc end kernels
782 call neg_check('mid',ipn,dt,qcheck,outqm,outtm,outum,outvm, &
783 outqcm,pretm,its,ite,kts,kte,itf,ktf,ktopm)
784 endif
786 if(ideep.eq.1)then
787 call cu_gf_deep_run( &
788 itf,ktf,its,ite, kts,kte &
789
790 ,dicycle &
791 ,ichoice &
792 ,ipr &
793 ,ccn_gf &
794 ,ccnclean &
795 ,dt &
796 ,0 &
797
798 ,kpbli &
799 ,dhdt &
800 ,xlandi &
801
802 ,zo &
803 ,forcing2 &
804 ,t2d &
805 ,q2d &
806 ,ter11 &
807 ,tn &
808 ,qo &
809 ,p2d &
810 ,psur &
811 ,us &
812 ,vs &
813 ,rhoi &
814 ,hfx &
815 ,qfx &
816 ,dx &
817 ,mconv &
818 ,omeg &
819 ,cactiv &
820 ,cnvwt &
821 ,zu &
822 ,zd &
823 ,zdm & ! hli
824 ,edt &
825 ,edtm & ! hli
826 ,xmb &
827 ,xmbm &
828 ,xmbs &
829 ,pret &
830 ,outu &
831 ,outv &
832 ,outt &
833 ,outq &
834 ,outqc &
835 ,kbcon &
836 ,ktop &
837 ,cupclw &
838 ,frhd &
839 ,ierr &
840 ,ierrc &
841 ,nchem &
842 ,fscav &
843 ,chem3d &
844 ,wetdpc_deep &
845 ,do_smoke_transport &
846! the following should be set to zero if not available
847 ,rand_mom & ! for stochastics mom, if temporal and spatial patterns exist
848 ,rand_vmas & ! for stochastics vertmass, if temporal and spatial patterns exist
849 ,rand_clos & ! for stochastics closures, if temporal and spatial patterns exist
850 ,spp_cu_deep & ! flag to what you want perturbed
851 ! 1 = momentum transport
852 ! 2 = normalized vertical mass flux profile
853 ! 3 = closures
854 ! more is possible, talk to developer or
855 ! implement yourself. pattern is expected to be
856 ! betwee -1 and +1
857 ,do_cap_suppress_here,cap_suppress_j &
858 ,k22 &
859 ,jmin,kdt,tropics)
860 jpr=0
861 ipr=0
862!$acc kernels
863 do i=its,itf
864 do k=kts,ktf
865 qcheck(i,k)=qv(i,k) +(outqs(i,k)+outqm(i,k))*dt
866 enddo
867 enddo
868!$acc end kernels
870 call neg_check('deep',ipn,dt,qcheck,outq,outt,outu,outv, &
871 outqc,pret,its,ite,kts,kte,itf,ktf,ktop)
872!
873 endif
874!$acc kernels
875 do i=its,itf
876 kcnv(i)=0
877 if(pretm(i).gt.0.)then
878 kcnv(i)= 1 !jmin(i)
879 cutenm(i)=1.
880 else
881 kbconm(i)=0
882 ktopm(i)=0
883 cutenm(i)=0.
884 endif ! pret > 0
885
886 if(pret(i).gt.0.)then
887 cuten(i)=1.
888 cutenm(i)=0.
889 pretm(i)=0.
890 kcnv(i)= 1 !jmin(i)
891 ktopm(i)=0
892 kbconm(i)=0
893 else
894 kbcon(i)=0
895 ktop(i)=0
896 cuten(i)=0.
897 endif ! pret > 0
898 enddo
899!$acc end kernels
900!
901!$acc parallel loop private(kstop,dtime_max,massflx,trcflx_in1,clw_in1,po_cup)
902 do i=its,itf
903 massflx(:)=0.
904 trcflx_in1(:)=0.
905 clw_in1(:)=0.
906 do k=kts,ktf
907 clw_ten(i, k)=0.
908 enddo
909 po_cup(:)=0.
910 kstop=kts
911 if(ktopm(i).gt.kts .or. ktop(i).gt.kts)kstop=max(ktopm(i),ktop(i))
912 if(ktops(i).gt.kts)kstop=max(kstop,ktops(i))
913 if(kstop.gt.2)then
914 htop(i)=kstop
915 if(kbcon(i).gt.2 .or. kbconm(i).gt.2)then
916 hbot(i)=max(kbconm(i),kbcon(i)) !jmin(i)
917 endif
918
919 dtime_max=dt
920 forcing2(i,3)=0.
921 do k=kts,kstop
922 cnvc(i,k) = 0.04 * log(1. + 675. * zu(i,k) * xmb(i)) + &
923 0.04 * log(1. + 675. * zum(i,k) * xmbm(i)) + &
924 0.04 * log(1. + 675. * zus(i,k) * xmbs(i))
925 cnvc(i,k) = min(cnvc(i,k), 0.6)
926 cnvc(i,k) = max(cnvc(i,k), 0.0)
927 cnvw(i,k)=cnvwt(i,k)*xmb(i)*dt+cnvwts(i,k)*xmbs(i)*dt+cnvwtm(i,k)*xmbm(i)*dt
928 ud_mf(i,k)=cuten(i)*zu(i,k)*xmb(i)*dt
929 dd_mf(i,k)=cuten(i)*zd(i,k)*edt(i)*xmb(i)*dt
930 t(i,k)=t(i,k)+dt*(cutens(i)*outts(i,k)+cutenm(i)*outtm(i,k)+outt(i,k)*cuten(i))
931 qv(i,k)=max(1.e-16,qv(i,k)+dt*(cutens(i)*outqs(i,k)+cutenm(i)*outqm(i,k)+outq(i,k)*cuten(i)))
932 gdc(i,k,7)=sqrt(us(i,k)**2 +vs(i,k)**2)
933 us(i,k)=us(i,k)+outu(i,k)*cuten(i)*dt +outum(i,k)*cutenm(i)*dt +outus(i,k)*cutens(i)*dt
934 vs(i,k)=vs(i,k)+outv(i,k)*cuten(i)*dt +outvm(i,k)*cutenm(i)*dt +outvs(i,k)*cutens(i)*dt
935
936 gdc(i,k,1)= max(0.,tun_rad_shall(i)*cupclws(i,k)*cutens(i)) ! my mod
937 !gdc2(i,k,1)=max(0.,tun_rad_deep(i)*(cupclwm(i,k)*cutenm(i)+cupclw(i,k)*cuten(i)))
938 gdc2(i,k,1)=max(0.,tun_rad_mid(i)*cupclwm(i,k)*cutenm(i)+frhd(i)*cupclw(i,k)*cuten(i)+tun_rad_shall(i)*cupclws(i,k)*cutens(i))
939 !gdc2(i,k,1) = min(0.1, max(0.01, tun_rad_mid(i)*frhm(i)))*cupclwm(i,k)*cutenm(i) + min(0.1, max(0.01, tun_rad_deep(i)*(frhd(i))))*cupclw(i,k)*cuten(i) + tun_rad_shall(i)*cupclws(i,k)*cutens(i)
940 qci_conv(i,k)=gdc2(i,k,1)
941 gdc(i,k,2)=(outt(i,k))*86400.
942 gdc(i,k,3)=(outtm(i,k))*86400.
943 gdc(i,k,4)=(outts(i,k))*86400.
944 gdc(i,k,7)=-(gdc(i,k,7)-sqrt(us(i,k)**2 +vs(i,k)**2))/dt
945 !gdc(i,k,8)=(outq(i,k))*86400.*xlv/cp
946 gdc(i,k,8)=(outqm(i,k)+outqs(i,k)+outq(i,k))*86400.*xlv/cp
947 gdc(i,k,9)=gdc(i,k,2)+gdc(i,k,3)+gdc(i,k,4)
948
950 dp=100.*(p2d(i,k)-p2d(i,k+1))
951 dtime_max=min(dtime_max,.5*dp)
952 po_cup(k)=.5*(p2d(i,k)+p2d(i,k+1))
953 if (clcw(i,k) .gt. -999.0 .and. clcw(i,k+1) .gt. -999.0 )then
954 clwtot = cliw(i,k) + clcw(i,k)
955 if(clwtot.lt.1.e-32)clwtot=0.
956 clwtot1= cliw(i,k+1) + clcw(i,k+1)
957 if(clwtot1.lt.1.e-32)clwtot1=0.
958 clw_in1(k)=clwtot
959 massflx(k)=-(xmb(i) *( zu(i,k)- edt(i)* zd(i,k))) &
960 -(xmbm(i)*(zdm(i,k)-edtm(i)*zdm(i,k))) &
961 -(xmbs(i)*zus(i,k))
962 trcflx_in1(k)=massflx(k)*.5*(clwtot+clwtot1)
963 forcing2(i,3)=forcing2(i,3)+clwtot
964 endif
965 enddo
966
967 massflx(1)=0.
968 trcflx_in1(1)=0.
969 call fct1d3 (kstop,kte,dtime_max,po_cup, &
970 clw_in1,massflx,trcflx_in1,clw_ten(i,:),g)
971
972 do k=1,kstop
973 tem = dt*(outqcs(i,k)*cutens(i)+outqc(i,k)*cuten(i) &
974 +outqcm(i,k)*cutenm(i) &
975 +clw_ten(i,k) &
976 )
977 tem1 = max(0.0, min(1.0, (tcr-t(i,k))*tcrf))
978 if (clcw(i,k) .gt. -999.0) then
979 cliw(i,k) = max(0.,cliw(i,k) + tem * tem1) ! ice
980 clcw(i,k) = max(0.,clcw(i,k) + tem *(1.0-tem1)) ! water
981 else
982 cliw(i,k) = max(0.,cliw(i,k) + tem)
983 endif
984
985 enddo
986
987 gdc(i,1,10)=forcing(i,1)
988 gdc(i,2,10)=forcing(i,2)
989 gdc(i,3,10)=forcing(i,3)
990 gdc(i,4,10)=forcing(i,4)
991 gdc(i,5,10)=forcing(i,5)
992 gdc(i,6,10)=forcing(i,6)
993 gdc(i,7,10)=forcing(i,7)
994 gdc(i,8,10)=forcing(i,8)
995 gdc(i,10,10)=xmb(i)
996 gdc(i,11,10)=xmbm(i)
997 gdc(i,12,10)=xmbs(i)
998 gdc(i,13,10)=hfx(i)
999 gdc(i,15,10)=qfx(i)
1000 gdc(i,16,10)=pret(i)*3600.
1001
1002 maxupmf(i)=0.
1003 if(forcing2(i,6).gt.0.)then
1004 maxupmf(i)=maxval(xmb(i)*zu(i,kts:ktf)/forcing2(i,6))
1005 endif
1006
1007 if(ktop(i).gt.2 .and.pret(i).gt.0.)dt_mf(i,ktop(i)-1)=ud_mf(i,ktop(i))
1008 endif
1009 enddo
1010!$acc end parallel
1011!$acc kernels
1012 do i=its,itf
1013 if(pret(i).gt.0.)then
1014 cactiv(i)=1
1015 raincv(i)=.001*(cutenm(i)*pretm(i)+cutens(i)*prets(i)+cuten(i)*pret(i))*dt
1016 else
1017 cactiv(i)=0
1018 if(pretm(i).gt.0)raincv(i)=.001*cutenm(i)*pretm(i)*dt
1019 endif ! pret > 0
1020
1021 if(pretm(i).gt.0)then
1022 cactiv_m(i)=1
1023 else
1024 cactiv_m(i)=0
1025 endif
1026
1027 ! Unify ccn
1028 if(ccn_m(i).lt.ccn_gf(i))then
1029 ccn_gf(i)=ccn_m(i)
1030 endif
1031
1032 if(ccn_gf(i)<0) ccn_gf(i)=0
1033
1034 ! Convert ccn back to aod
1035 aod_gf(i)=0.0027*(ccn_gf(i)**0.64)
1036 if(aod_gf(i)<0.007)then
1037 aod_gf(i)=0.007
1038 ccn_gf(i)=(aod_gf(i)/0.0027)**(1/0.640)
1039 elseif(aod_gf(i)>aodc0)then
1040 aod_gf(i)=aodc0
1041 ccn_gf(i)=(aod_gf(i)/0.0027)**(1/0.640)
1042 endif
1043 enddo
1044!$acc end kernels
1045 100 continue
1046!
1047! Scale dry mixing ratios for water wapor and cloud water to specific humidy / moist mixing ratios
1048!
1049!$acc kernels
1050 qv_spechum = qv/(1.0_kind_phys+qv)
1051 cnvw_moist = cnvw/(1.0_kind_phys+qv)
1052!$acc end kernels
1053!
1054! Diagnostic tendency updates
1055!
1056 if(ldiag3d) then
1057 if(ishallow_g3.eq.1 .and. .not.flag_for_scnv_generic_tend) then
1058 uidx=dtidx(index_of_x_wind,index_of_process_scnv)
1059 vidx=dtidx(index_of_y_wind,index_of_process_scnv)
1060 tidx=dtidx(index_of_temperature,index_of_process_scnv)
1061 qidx=dtidx(100+ntqv,index_of_process_scnv)
1062 if(uidx>=1) then
1063!$acc kernels
1064 do k=kts,ktf
1065 dtend(:,k,uidx) = dtend(:,k,uidx) + cutens(:)*outus(:,k) * dt
1066 enddo
1067!$acc end kernels
1068 endif
1069 if(vidx>=1) then
1070!$acc kernels
1071 do k=kts,ktf
1072 dtend(:,k,vidx) = dtend(:,k,vidx) + cutens(:)*outvs(:,k) * dt
1073 enddo
1074!$acc end kernels
1075 endif
1076 if(tidx>=1) then
1077!$acc kernels
1078 do k=kts,ktf
1079 dtend(:,k,tidx) = dtend(:,k,tidx) + cutens(:)*outts(:,k) * dt
1080 enddo
1081!$acc end kernels
1082 endif
1083 if(qidx>=1) then
1084!$acc kernels
1085 do k=kts,ktf
1086 do i=its,itf
1087 tem = cutens(i)*outqs(i,k)* dt
1088 tem = tem/(1.0_kind_phys+tem)
1089 dtend(i,k,qidx) = dtend(i,k,qidx) + tem
1090 enddo
1091 enddo
1092!$acc end kernels
1093 endif
1094 endif
1095 if((ideep.eq.1. .or. imid_gf.eq.1) .and. .not.flag_for_dcnv_generic_tend) then
1096 uidx=dtidx(index_of_x_wind,index_of_process_dcnv)
1097 vidx=dtidx(index_of_y_wind,index_of_process_dcnv)
1098 tidx=dtidx(index_of_temperature,index_of_process_dcnv)
1099 if(uidx>=1) then
1100!$acc kernels
1101 do k=kts,ktf
1102 dtend(:,k,uidx) = dtend(:,k,uidx) + (cuten*outu(:,k)+cutenm*outum(:,k)) * dt
1103 enddo
1104!$acc end kernels
1105 endif
1106 if(vidx>=1) then
1107!$acc kernels
1108 do k=kts,ktf
1109 dtend(:,k,vidx) = dtend(:,k,vidx) + (cuten*outv(:,k)+cutenm*outvm(:,k)) * dt
1110 enddo
1111!$acc end kernels
1112 endif
1113 if(tidx>=1) then
1114!$acc kernels
1115 do k=kts,ktf
1116 dtend(:,k,tidx) = dtend(:,k,tidx) + (cuten*outt(:,k)+cutenm*outtm(:,k)) * dt
1117 enddo
1118!$acc end kernels
1119 endif
1120
1121 qidx=dtidx(100+ntqv,index_of_process_dcnv)
1122 if(qidx>=1) then
1123!$acc kernels
1124 do k=kts,ktf
1125 do i=its,itf
1126 tem = (cuten(i)*outq(i,k) + cutenm(i)*outqm(i,k))* dt
1127 tem = tem/(1.0_kind_phys+tem)
1128 dtend(i,k,qidx) = dtend(i,k,qidx) + tem
1129 enddo
1130 enddo
1131!$acc end kernels
1132 endif
1133 endif
1134 if(allocated(clcw_save)) then
1135!$acc parallel loop collapse(2) private(tem_shal,tem_deep,tem,tem1,weight_sum,cliw_both,clcw_both)
1136 do k=kts,ktf
1137 do i=its,itf
1138 tem_shal = dt*(outqcs(i,k)*cutens(i)+outqcm(i,k)*cutenm(i))
1139 tem_deep = dt*(outqc(i,k)*cuten(i)+clw_ten(i,k))
1140 tem = tem_shal+tem_deep
1141 tem1 = max(0.0, min(1.0, (tcr-t(i,k))*tcrf))
1142 weight_sum = abs(tem_shal)+abs(tem_deep)
1143 if(weight_sum<1e-12) then
1144 cycle
1145 endif
1146
1147 if (clcw_save(i,k) .gt. -999.0) then
1148 cliw_both = max(0.,cliw_save(i,k) + tem * tem1) - cliw_save(i,k)
1149 clcw_both = max(0.,clcw_save(i,k) + tem) - clcw_save(i,k)
1150 else if(cliw_idx>=1) then
1151 cliw_both = max(0.,cliw_save(i,k) + tem) - cliw_save(i,k)
1152 clcw_both = 0
1153 endif
1154 if(cliw_deep_idx>=1) then
1155 dtend(i,k,cliw_deep_idx) = dtend(i,k,cliw_deep_idx) + abs(tem_deep)/weight_sum*cliw_both
1156 endif
1157 if(clcw_deep_idx>=1) then
1158 dtend(i,k,clcw_deep_idx) = dtend(i,k,clcw_deep_idx) + abs(tem_deep)/weight_sum*clcw_both
1159 endif
1160 if(cliw_shal_idx>=1) then
1161 dtend(i,k,cliw_shal_idx) = dtend(i,k,cliw_shal_idx) + abs(tem_shal)/weight_sum*cliw_both
1162 endif
1163 if(clcw_shal_idx>=1) then
1164 dtend(i,k,clcw_shal_idx) = dtend(i,k,clcw_shal_idx) + abs(tem_shal)/weight_sum*clcw_both
1165 endif
1166 enddo
1167 enddo
1168!$acc end parallel
1169 endif
1170 endif
1171 end subroutine cu_gf_driver_run
1173end module cu_gf_driver
subroutine fct1d3(ktop, n, dt, z, tracr, massflx, trflx_in, dellac, g)
Calculates tracer fluxes due to subsidence, only up-stream differencing is currently used but flux co...
subroutine cu_gf_deep_run(itf, ktf, its, ite, kts, kte, dicycle, ichoice, ipr, ccn, ccnclean, dtime, imid, kpbl, dhdt, xland, zo, forcing, t, q, z1, tn, qo, po, psur, us, vs, rho, hfx, qfx, dx, mconv, omeg, csum, cnvwt, zuo, zdo, zdm, edto, edtm, xmb_out, xmbm_in, xmbs_in, pre, outu, outv, outt, outq, outqc, kbcon, ktop, cupclw, frh_out, ierr, ierrc, nchem, fscav, chem3d, wetdpc_deep, do_smoke_transport, rand_mom, rand_vmas, rand_clos, nranflag, do_capsuppress, cap_suppress_j, k22, jmin, kdt, tropics)
Driver for the deep or congestus GF routine.
subroutine neg_check(name, j, dt, q, outq, outt, outu, outv, outqc, pret, its, ite, kts, kte, itf, ktf, ktop)
Checks for negative or excessive tendencies and corrects in a mass conversing way by adjusting the cl...
subroutine, public cu_gf_driver_init(imfshalcnv, imfshalcnv_gf, imfdeepcnv, imfdeepcnv_gf, mpirank, mpiroot, errmsg, errflg)
subroutine, public cu_gf_driver_run(ntracer, garea, im, km, dt, flag_init, flag_restart, cactiv, cactiv_m, g, cp, xlv, r_v, forcet, forceqv_spechum, phil, raincv, qv_spechum, t, cld1d, us, vs, t2di, w, qv2di_spechum, p2di, psuri, hbot, htop, kcnv, xland, hfx2, qfx2, aod_gf, cliw, clcw, pbl, ud_mf, dd_mf, dt_mf, cnvw_moist, cnvc, imfshalcnv, flag_for_scnv_generic_tend, flag_for_dcnv_generic_tend, dtend, dtidx, ntqv, ntiw, ntcw, index_of_temperature, index_of_x_wind, index_of_y_wind, index_of_process_scnv, index_of_process_dcnv, fhour, fh_dfi_radar, ix_dfi_radar, num_dfi_radar, cap_suppress, dfi_radar_max_intervals, ldiag3d, qci_conv, do_cap_suppress, maxupmf, maxmf, do_mynnedmf, ichoice_in, ichoicem_in, ichoice_s_in, spp_cu_deep, spp_wts_cu_deep, nchem, chem3d, fscav, wetdpc_deep, do_smoke_transport, kdt, errmsg, errflg)
This is the Grell-Freitas convection scheme driver module.
subroutine cu_gf_sh_run(us, vs, zo, t, q, z1, tn, qo, po, psur, dhdt, kpbl, rho, hfx, qfx, xland, ichoice, tcrit, dtime, zuo, xmb_out, kbcon, ktop, k22, ierr, ierrc, outt, outq, outqc, outu, outv, cnvwt, pre, cupclw, itf, ktf, its, ite, kts, kte, ipr, tropics)
Definition cu_gf_sh.F90:71
This module contains the Grell_Freitas deep convection scheme.
Definition cu_gf_deep.F90:5
This module contains the scale-aware Grell-Freitas cumulus scheme driver.
This module contains the Grell-Freitas shallow convection scheme.
Definition cu_gf_sh.F90:5