1 !WRF:MEDIATION_LAYER:PHYSICS
4 MODULE module_shallowcu_driver
8 real, private, parameter :: QCLOUD_MIN = 0.0
9 real, private, parameter :: QVAPOR_MIN = 0.0
12 SUBROUTINE shallowcu_driver( &
13 ! Order dependent args for domain, mem, and tile dims
14 ids,ide, jds,jde, kds,kde &
15 ,ims,ime, jms,jme, kms,kme &
16 ,ips,ipe, jps,jpe, kps,kpe &
17 ,i_start,i_end,j_start,j_end,kts,kte,num_tiles &
18 ! Order independent args (use VAR= in call)
24 ,itimestep,dt,dx,dx2d,area2d &
25 ,cudt,curr_secs,adapt_step_flag &
26 ,rainsh,pratesh,nca,rainshv &
27 ,z,z_at_w,dz8w,mavail,pblh,p8w &
29 ,cldfra,cldfra_old,cldfra_old_mp,cldfra_conv &
32 ! Package selection variables
34 ! Optional moisture tracers
35 ,qv_curr, qc_curr, qr_curr &
36 ,qi_curr, qs_curr, qg_curr &
41 ! Optional output arguments for CAMZM scheme
42 ,dlf, rliq, rliq2,dlf2 &
44 ! Optional output arguments for CAMUW scheme
45 ,cush, snowsh, icwmrsh, rprdsh, cbmf, cmfsl &
47 ! Optional moisture and other tendencies
48 ,rqvshten,rqcshten,rqrshten &
49 ,rqishten,rqsshten,rqgshten &
50 ,rqcnshten,rqinshten &
53 ,rthshten,rthraten,rthblten,rthften &
54 ! Optional moisture tracer flags
57 ,ht,shfrc3d,is_CAMMGMP_used &
58 ! for grims shallow convection with ysupbl
59 ,wstar,delta,kpbl,znu,raincv &
60 ! for nscv shallow convection
61 ,w,xland,hfx,qfx,mp_physics,pgcon &
62 ,RDCASHTEN, RQCDCSHTEN, W0AVG &
64 ,cldareaa, cldareab, cldliqa, cldliqb &
65 ,cldfra_sh, ca_rad, cw_rad &
67 ,rainshvb, capesave, radsave &
70 ,el_pbl, rthratenlw, rthratensw, exch_h &
71 ,dnw, xtime, xtime1, gmt &
72 ,qke,PBLHAVG, TKEAVG &
80 ,pert_deng_qv, pert_deng_qc, pert_deng_t &
83 !----------------------------------------------------------------------
84 USE module_model_constants
85 USE module_state_description, ONLY: CAMUWSHCUSCHEME &
93 ! *** add new modules of schemes here
95 USE module_shcu_camuwshcu_driver, ONLY : camuwshcu_driver
96 USE module_shcu_grims , ONLY : grims
97 USE module_shcu_nscv , ONLY : shcu_nscv
98 USE module_shcu_deng, ONLY: deng_shcu_driver
101 USE module_domain, ONLY: domain
102 #if ( WRF_CHEM == 1 )
103 USE module_state_description, ONLY: num_chem
106 ! This driver calls subroutines for the shallow cumulus
109 ! 1. G3 shallow cumulus
110 ! 2. UW shallow cumulus from CAM
111 ! 3. GRIMs shallow cumulus from GRIMs (available only with ysupbl)
112 ! 4. NCEP shallow cumulus scheme (separated from NSAS cumulus scheme)
113 ! (Han and Pan 2011, YSU)
115 !----------------------------------------------------------------------
117 !======================================================================
118 ! Grid structure in physics part of WRF
119 !----------------------------------------------------------------------
120 ! The horizontal velocities used in the physics are unstaggered
121 ! relative to temperature/moisture variables. All predicted
122 ! variables are carried at half levels except w, which is at full
123 ! levels. Some arrays with names (*8w) are at w (full) levels.
125 !----------------------------------------------------------------------
126 ! In WRF, kms (smallest number) is the bottom level and kme (largest
127 ! number) is the top level. In your scheme, if 1 is at the top level,
128 ! then you have to reverse the order in the k direction.
130 ! kme - half level (no data at this level)
131 ! kme ----- full level
133 ! kme-1 ----- full level
138 ! kms+2 ----- full level
140 ! kms+1 ----- full level
142 ! kms ----- full level
144 !======================================================================
147 ! Rho_d dry density (kg/m^3)
148 ! Theta_m moist potential temperature (K)
149 ! Qv water vapor mixing ratio (kg/kg)
150 ! Qc cloud water mixing ratio (kg/kg)
151 ! Qr rain water mixing ratio (kg/kg)
152 ! Qi cloud ice mixing ratio (kg/kg)
153 ! Qs snow mixing ratio (kg/kg)
154 !-----------------------------------------------------------------
155 !-- DT time step (second)
156 !-- CUDT cumulus time step (minute)
157 !-- curr_secs current forecast time (seconds)
158 !-- itimestep number of time step (integer)
159 !-- DX horizontal space interval (m)
160 !-- rr dry air density (kg/m^3)
162 !-- RUSHTEN Zonal wind tendency due to shallow
163 ! cumulus scheme precipitation (m/s/s)
164 !-- RVSHTEN Meridional wind tendency due to
165 ! cumulus scheme precipitation (m/s/s)
166 !-- RTHSHTEN Theta tendency due to shallow
167 ! cumulus scheme precipitation (K/s)
168 !-- RQVSHTEN Qv tendency due to shallow
169 ! cumulus scheme precipitation (kg/kg/s)
170 !-- RQRSHTEN Qr tendency due to shallow
171 ! cumulus scheme precipitation (kg/kg/s)
172 !-- RQCSHTEN Qc tendency due to shallow
173 ! cumulus scheme precipitation (kg/kg/s)
174 !-- RQSSHTEN Qs tendency due to shallow
175 ! cumulus scheme precipitation (kg/kg/s)
176 !-- RQISHTEN Qi tendency due to shallow
177 ! cumulus scheme precipitation (kg/kg/s)
178 !-- RQGSHTEN Qg tendency due to shallow
179 ! cumulus scheme precipitation (kg/kg/s)
181 !-- RAINSH accumulated total shallow cumulus scheme precipitation (mm)
182 !-- RAINSHV time-step shallow cumulus scheme precipitation (mm)
183 !-- PRATESH precipitiation rate from shallow cumulus scheme (mm/s)
184 !-- NCA counter of the cloud relaxation
185 ! time in KF cumulus scheme (integer)
186 !-- u_phy u-velocity interpolated to theta points (m/s)
187 !-- v_phy v-velocity interpolated to theta points (m/s)
188 !-- th_phy potential temperature (K)
189 !-- t_phy temperature (K)
190 !-- tke_pbl turbulent kinetic energy from PBL scheme (m2/s2)
191 !-- w vertical velocity (m/s)
192 !-- moist moisture array (4D - last index is species) (kg/kg)
193 !-- z height above sea level at middle of layers (m)
194 !-- z_at_w height above sea level at layer interfaces (m)
195 !-- dz8w dz between full levels (m)
196 !-- pblh planetary boundary layer height (m)
197 !-- mavail soil moisture availability
198 !-- p8w pressure at full levels (Pa)
199 !-- p_phy pressure (Pa)
200 !-- pi_phy the exner function, (p/p0)**(R/Cp) (dimensionless)
201 ! points (dimensionless)
202 !-- hfx upward heat flux at surface (W/m2)
203 !-- RTHRATEN radiative temp forcing for Grell-Devenyi scheme
204 !-- RTHBLTEN PBL temp forcing for Grell-Devenyi scheme
205 !-- RQVBLTEN PBL moisture forcing for Grell-Devenyi scheme
208 !-- cldfra cloud fraction
209 !-- cldfra_old cloud fraction from previous time step
210 !-- cldfrash cloud fraction from shallow Cu
211 !-- cldfra_old_mp cloud fraction from previous time step if CAMMGMP microphysics is used
212 !-- cldfra_conv convective cloud fraction
213 !-- rho density (kg/m^3)
214 !-- XLV0 latent heat of vaporization constant
215 ! used in temperature dependent formula (J/kg)
216 !-- XLV1 latent heat of vaporization constant
217 ! used in temperature dependent formula (J/kg/K)
218 !-- XLS0 latent heat of sublimation constant
219 ! used in temperature dependent formula (J/kg)
220 !-- XLS1 latent heat of sublimation constant
221 ! used in temperature dependent formula (J/kg/K)
222 !-- R_d gas constant for dry air ( 287. J/kg/K)
223 !-- R_v gas constant for water vapor (461 J/k/kg)
224 !-- Cp specific heat at constant pressure (1004 J/k/kg)
225 !-- rvovrd R_v divided by R_d (dimensionless)
226 !-- G acceleration due to gravity (m/s^2)
227 !-- EP_1 constant for virtual temperature
228 ! (R_v/R_d - 1) (dimensionless)
229 !--shfrc3d Shallow cloud fraction
230 !-- cmfmc Deep + Shallow Convective m
231 !-- ids start index for i in domain
232 !-- ide end index for i in domain
233 !-- jds start index for j in domain
234 !-- jde end index for j in domain
235 !-- kds start index for k in domain
236 !-- kde end index for k in domain
237 !-- ims start index for i in memory
238 !-- ime end index for i in memory
239 !-- jms start index for j in memory
240 !-- jme end index for j in memory
241 !-- kms start index for k in memory
242 !-- kme end index for k in memory
243 !-- i_start start indices for i in tile
244 !-- i_end end indices for i in tile
245 !-- j_start start indices for j in tile
246 !-- j_end end indices for j in tile
247 !-- kts start index for k in tile
248 !-- kte end index for k in tile
249 !-- num_tiles number of tiles
250 !-- HBOT index of lowest model layer with convection
251 !-- HTOP index of highest model layer with convection
252 !-- LBOT index of lowest model layer with convection
253 !-- LTOP index of highest model layer with convection
254 !-- periodic_x T/F this is using periodic lateral boundaries in the X direction
255 !-- periodic_y T/F this is using periodic lateral boundaries in the Y-direction
257 !======================================================================
258 LOGICAL, INTENT(IN ) :: is_CAMMGMP_used !BSINGH:01/31/2013: Added for CAMUWSHCU
259 INTEGER, INTENT(IN ) :: &
260 ids,ide, jds,jde, kds,kde, &
261 ims,ime, jms,jme, kms,kme, &
264 #if ( WRF_CHEM == 1 )
265 INTEGER, INTENT(IN ) :: chem_opt
268 INTEGER, DIMENSION(num_tiles), INTENT(IN) :: &
269 & i_start,i_end,j_start,j_end
271 INTEGER, INTENT(IN ) :: &
274 INTEGER, INTENT(IN ) :: shcu_physics
276 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_moist ), &
279 #if ( WRF_CHEM == 1 )
280 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
286 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
304 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
308 REAL, DIMENSION( ims:ime , jms:jme ), INTENT(IN) :: &
312 REAL, DIMENSION( ims:ime , jms:jme ), &
313 INTENT(INOUT) :: RAINSH &
319 REAL, DIMENSION( ims:ime , jms:jme ),INTENT(INOUT),OPTIONAL :: &
321 REAL, DIMENSION( ims:ime , jms:jme ) :: tmppratesh
323 REAL, INTENT(IN ) :: DT, DX
324 REAL, DIMENSION( ims:ime , jms:jme ), INTENT(IN), OPTIONAL :: &
326 INTEGER, INTENT(IN ),OPTIONAL :: &
327 ips,ipe, jps,jpe, kps,kpe
328 REAL, INTENT(IN ),OPTIONAL :: CUDT
329 REAL, INTENT(IN ),OPTIONAL :: CURR_SECS
330 LOGICAL,INTENT(IN ),OPTIONAL :: adapt_step_flag
331 REAL :: cudt_pass, curr_secs_pass
332 LOGICAL :: adapt_step_flag_pass
337 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
338 OPTIONAL, INTENT(INOUT) :: &
339 ! optional moisture tracers
340 qv_curr, qc_curr, qr_curr &
341 ,qi_curr, qs_curr, qg_curr &
342 ! optional scalar tracers !BSINGH
344 ! optional moisture and other tendencies
345 ,rqvshten,rqcshten,rqrshten &
346 ,rqishten,rqsshten,rqgshten &
347 ,rqcnshten,rqinshten &
350 ,rthften,rushten,rvshten,rthshten
352 REAL, DIMENSION( ims:ime , jms:jme ), &
353 OPTIONAL, INTENT(INOUT) :: &
356 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
357 OPTIONAL, INTENT(INOUT) :: &
358 cldfrash, cmfsl, cmflq, icwmrsh, &
360 cmfmc, cmfmc2, rprdsh
361 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
363 dlf2 ! Required by CAMMGMP Microphysics
364 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
366 shfrc3d ! Required by wet scavenging code in WRF_CHEM
367 ! for grims shallow convection with ysupbl
369 REAL, DIMENSION( ims:ime, jms:jme ) , &
370 OPTIONAL, INTENT(IN ) :: wstar
371 REAL, DIMENSION( ims:ime, jms:jme ) , &
372 OPTIONAL, INTENT(IN ) :: delta
373 REAL, DIMENSION( ims:ime, jms:jme ) , &
374 OPTIONAL, INTENT(IN ) :: raincv
375 REAL, DIMENSION( kms:kme ) , &
376 OPTIONAL, INTENT(IN ) :: znu
377 INTEGER, DIMENSION( ims:ime , jms:jme ) , &
378 OPTIONAL, INTENT(IN) :: kpbl
380 ! for nscv shallow convection
382 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ) , &
384 REAL, DIMENSION( ims:ime, jms:jme ) , &
386 REAL, DIMENSION( ims:ime, jms:jme ) , &
387 OPTIONAL, INTENT(IN ) :: hfx,qfx
388 INTEGER, INTENT(IN), OPTIONAL :: mp_physics
389 REAL, INTENT(IN), OPTIONAL :: pgcon
392 ! PSU-DENG Shallow cu variables
394 INTEGER, INTENT(IN ) :: bl_pbl_physics
395 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
398 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
399 INTENT(INOUT) :: RDCASHTEN, RQCDCSHTEN
401 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
402 INTENT(INOUT) :: cldareaa, cldareab, &
403 cldliqa, cldliqb, wub, &
405 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
406 INTENT( OUT) :: cldfra_sh, ca_rad, cw_rad
408 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
409 INTENT(INOUT) :: el_pbl, rthratenlw, rthratensw, exch_h
411 REAL, DIMENSION( kms:kme ), INTENT(IN ) :: &
414 REAL, INTENT(IN ) :: xtime, gmt
416 REAL, DIMENSION( ims:ime, jms:jme ), &
418 cldtopb, pblmax, rainshvb, capesave, radsave, &
419 clddpthb, xtime1, PBLHAVG, MAVAIL
421 REAL, DIMENSION( ims:ime, jms:jme ), &
424 INTEGER, DIMENSION( ims:ime, jms:jme ), &
426 ltopb, kdcldtop, kdcldbas
428 REAL, DIMENSION( ims:ime, 1:100, jms:jme ), &
432 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ) :: tke_scr, kth_scr, bbls_scr
435 integer, intent (in) :: multi_perturb
436 logical, intent (in) :: pert_deng
437 real, intent(in) :: pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w
438 real, dimension(ims:ime, kms:kme, jms:jme), optional, intent(inout) :: perts_qvapor, &
439 perts_qcloud, perts_th, perts_w
442 ! End PSU-DENG Shallow cu variables
446 ! Flags relating to the optional tendency arrays declared above
447 ! Models that carry the optional tendencies will provdide the
448 ! optional arguments at compile time; these flags all the model
449 ! to determine at run-time whether a particular tracer is in
452 LOGICAL, INTENT(IN), OPTIONAL :: &
463 INTEGER :: i,j,k,its,ite,jts,jte,ij
464 CHARACTER(len=200) :: message
467 !-----------------------------------------------------------------
469 if (.not. PRESENT(CURR_SECS)) then
472 curr_secs_pass = curr_secs
475 if (.not. PRESENT(CUDT)) then
481 if (.not. PRESENT(adapt_step_flag)) then
482 adapt_step_flag_pass = .false.
484 adapt_step_flag_pass = adapt_step_flag
487 ! Initialize tmppratesh to pratesh
489 if ( PRESENT ( pratesh ) ) then
490 tmppratesh(:,:) = pratesh(:,:)
496 IF (shcu_physics .eq. 0) return
499 ! DON'T JUDGE TIME STEP HERE, SINCE KF NEEDS ACCUMULATED W FIELD.
500 ! DO IT INSIDE THE INDIVIDUAL CUMULUS SCHEME
502 ! SET START AND END POINTS FOR TILES
504 !$OMP PRIVATE ( ij ,its,ite,jts,jte, i,j,k)
505 DO ij = 1 , num_tiles
512 scps_select: SELECT CASE(shcu_physics)
515 ! This setting takes the place of ishallow in v3.1.1+
517 CASE (CAMUWSHCUSCHEME)
518 CALL wrf_debug(100,'in camuw_scps')
520 WRITE( message , * ) 'This shallow cumulus option requires ice microphysics option: f_qi = ', f_qi
521 CALL wrf_error_fatal ( message )
523 CALL camuwshcu_driver( &
524 IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde &
525 ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme &
526 ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte &
527 ,NUM_MOIST=num_moist, DT=dt &
528 ,P=p, P8W=p8w, PI_PHY=pi &
529 ,Z=z, Z_AT_W=z_at_w, DZ8W=dz8w &
530 ,T_PHY=t, U_PHY=u, V_PHY=v &
531 ,MOIST=moist, QV=qv_curr, QC=qc_curr, QI=qi_curr &
532 ,QNC=qnc_curr, QNI=qni_curr &
533 #if ( WRF_CHEM == 1 )
534 ,CHEM=chem, CHEM_OPT=chem_opt &
536 ,PBLH_IN=pblh, TKE_PBL=tke_pbl &
537 ,CLDFRA=cldfra, CLDFRA_OLD=cldfra_old &
538 ,CLDFRA_OLD_MP=cldfra_old_mp &
539 ,CLDFRA_CONV=cldfra_conv,IS_CAMMGMP_USED=is_CAMMGMP_used &
541 ,CUSH_INOUT=cush, PRATESH=tmppratesh &
543 ,ICWMRSH=icwmrsh, CMFMC=cmfmc, CMFMC2_INOUT=cmfmc2 &
544 ,RPRDSH_INOUT=rprdsh, CBMF_INOUT=cbmf &
545 ,CMFSL=cmfsl, CMFLQ=cmflq, DLF=dlf,DLF2=dlf2 & !DLF2 is required by CAMMGMP microphysics
546 ,EVAPCSH_INOUT=evapcsh &
547 ,RLIQ=rliq, RLIQ2_INOUT=rliq2, CUBOT=hbot, CUTOP=htop &
548 ,RUSHTEN=rushten, RVSHTEN=rvshten, RTHSHTEN=rthshten &
549 ,RQVSHTEN=rqvshten, RQCSHTEN=rqcshten, RQRSHTEN=rqrshten &
550 ,RQISHTEN=rqishten, RQSSHTEN=rqsshten, RQGSHTEN=rqgshten &
551 ,RQCNSHTEN=rqcnshten,RQINSHTEN=rqinshten &
552 ,HT=ht,SHFRC3D=shfrc3d,ITIMESTEP=itimestep &
555 CASE (GRIMSSHCUSCHEME)
556 CALL wrf_debug(100,'in grims_scps')
557 IF ( PRESENT( wstar ) ) THEN
560 ,P3DI=p8w,P3D=p,PI3D=pi,Z3DI=Z_AT_W &
561 ,WSTAR=wstar,HPBL=pblh,DELTA=delta &
562 ,RTHSHTEN=rthshten,RQVSHTEN=rqvshten &
563 ,DT=dt,G=g,XLV=xlv,RD=r_d,RV=r_v &
564 ,RCP=rcp,P1000MB=p1000mb &
565 ,KPBL2D=kpbl,ZNU=znu,RAINCV=raincv &
566 ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde &
567 ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme &
568 ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte &
572 CASE (NSCVSHCUSCHEME)
573 CALL wrf_debug(100,'in nscv_scps')
574 IF ( PRESENT ( QFX ) .AND. PRESENT( HFX ) ) THEN
576 DT=dt,P3DI=p8w,P3D=p,PI3D=pi &
577 ,QC3D=QC_CURR,QI3D=QI_CURR,RHO3D=rho &
578 ,QV3D=QV_CURR,T3D=t &
580 ,XLAND=XLAND,DZ8W=dz8w,W=w,U3D=u,V3D=v &
581 ,HPBL=pblh,HFX=hfx,QFX=qfx &
582 ,MP_PHYSICS=mp_physics &
584 ,CP=cp,CLIQ=cliq,CPV=cpv,G=g,XLV=xlv,R_D=r_d &
585 ,R_V=r_v,EP_1=ep_1,EP_2=EP_2 &
586 ,CICE=cice,XLS=xls,PSAT=psat &
587 ,F_QI=f_qi,F_QC=f_qc &
588 ,RTHSHTEN=RTHSHTEN,RQVSHTEN=RQVSHTEN &
589 ,RQCSHTEN=RQCSHTEN,RQISHTEN=RQISHTEN &
590 ,RUSHTEN=RUSHTEN,RVSHTEN=RVSHTEN &
591 ,PRATESH=tmppratesh &
592 ,HBOT=HBOT,HTOP=HTOP &
593 ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde &
594 ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme &
595 ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte &
598 CALL wrf_error_fatal('Lacking arguments for SHCU_NSCV in shallow cumulus driver')
601 CASE (DENGSHCUSCHEME)
602 CALL wrf_debug(100,'in PSU-DENG Shallow Cu CPS')
605 IF ( bl_pbl_physics == 2 ) THEN ! MYJ
609 tke_scr(I,K,J) = tke_pbl(I,K,J)
610 kth_scr(I,K,J) = exch_h(I,K,J)
611 bbls_scr(I,K,J) = el_pbl(I,K,J)
615 ELSE IF ( bl_pbl_physics == 5 ) THEN ! MYNN
619 tke_scr(I,K,J) = 0.5 * qke(I,K,J)
620 kth_scr(I,K,J) = exch_h(I,K,J)
621 bbls_scr(I,K,J) = el_pbl(I,K,J)
626 WRITE(message,*) 'PSU DENG ShCu currently does not support PBL option: ', &
628 CALL wrf_error_fatal ( message )
632 if (pert_deng .and. multi_perturb == 1) &
633 call Add_multi_perturb_shcu_perturbations (perts_qvapor, perts_qcloud, &
634 perts_th, perts_w, pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w, &
635 th, qv_curr, qc_curr, w, its, ite, jts, jte, ims, ime, jms, jme, kms, kme, &
638 CALL deng_shcu_driver( &
639 IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde &
640 ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme &
641 ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte &
642 ,DT=dt ,KTAU=itimestep ,DX=dx, xtime=xtime, gmt=gmt &
644 ,XLV0=XLV0 ,XLV1=XLV1 ,XLS0=XLS0 ,XLS1=XLS1 &
646 ,SVP1=SVP1 ,SVP2=SVP2 ,SVP3=SVP3 ,SVPT0=SVPT0 &
647 ,ADAPT_STEP_FLAG=ADAPT_STEP_FLAG, DSIGMA=dnw &
648 ,XLONG=xlong, HT=ht, PBLH=pblh &
649 ,U=u,V=v, w=w, TH=th ,T=t &
650 ,QV=qv_curr, QC=qc_curr, QR=qr_curr, DZ8W=dz8w &
651 ,PCPS=p, RHO=rho, Z_AT_W=z_at_w, PI=pi &
652 ,TKE=tke_scr, kth=kth_scr, bbls=bbls_scr &
653 ,ten_radl=rthratenlw, ten_rads=rthratensw &
654 ,RAINSH=RAINSH, RAINSHV=RAINSHV,rainshvb=rainshvb &
655 ,pblmax=pblmax, capesave=capesave, xtime1=xtime1 &
656 ,radsave=radsave, clddpthb=clddpthb, cldtopb=cldtopb &
657 ,MAVAIL=MAVAIL, PBLHAVG=PBLHAVG &
659 ,ltopb=ltopb, kdcldtop=kdcldtop, kdcldbas=kdcldbas &
660 ,W0AVG=W0AVG, TKEAVG=TKEAVG &
661 ,cldareaa=cldareaa, cldareab=cldareab &
662 ,cldliqa=cldliqa, cldliqb=cldliqb &
663 ,cldfra_sh=cldfra_sh,ca_rad=ca_rad, cw_rad=cw_rad, wub=wub &
664 ,RUSHTEN=rushten, RVSHTEN=rvshten, RTHSHTEN=rthshten &
665 ,RQVSHTEN=rqvshten, RQCSHTEN=rqcshten, RQRSHTEN=rqrshten &
666 ,RDCASHTEN=RDCASHTEN, RQCDCSHTEN=RQCDCSHTEN)
668 if (pert_deng .and. multi_perturb == 1) &
669 call Remove_multi_perturb_shcu_perturbations (perts_qvapor, perts_qcloud, &
670 perts_th, perts_w, pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w, &
671 th, qv_curr, qc_curr, w, its, ite, jts, jte, ims, ime, jms, jme, kms, kme, &
675 WRITE( message , * ) 'The shallow cumulus option does not exist: shcu_physics = ', shcu_physics
676 CALL wrf_error_fatal ( message )
678 END SELECT scps_select
681 !$OMP END PARALLEL DO
684 ! Copy pratesh back to output array, if necessary.
686 if (PRESENT(PRATESH)) then
687 pratesh(:,:) = tmppratesh(:,:)
688 if (PRESENT(RAINSHV)) then
689 rainshv(:,:) = pratesh(:,:)*dt
693 END SUBROUTINE shallowcu_driver
695 subroutine Add_multi_perturb_shcu_perturbations (perts_qvapor, perts_qcloud, &
696 perts_th, perts_w, pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w, &
697 th, qv_curr, qc_curr, w, its, ite, jts, jte, ims, ime, jms, jme, kms, kme, &
702 integer, intent(in) :: its, ite, jts, jte, ims, ime, jms, jme, kms, kme, kts, kte
703 real, intent(in) :: pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w
704 real, dimension(ims:ime, kms:kme, jms:jme), optional, intent (in) :: perts_qvapor, &
705 perts_qcloud, perts_th, perts_w
706 real, dimension(ims:ime, kms:kme, jms:jme), intent (inout) :: th, w, qv_curr, qc_curr
714 qc_curr(i, k, j) = max ((1.0 + perts_qcloud(i, k, j) * pert_deng_qc) * qc_curr(i, k, j), QCLOUD_MIN)
715 qv_curr(i, k, j) = max ((1.0 + perts_qvapor(i, k, j) * pert_deng_qv) * qv_curr(i, k, j), QVAPOR_MIN)
716 th(i, k, j) = (1.0 + perts_th(i, k, j) * pert_deng_t) * th(i, k, j)
724 w(i, k, j) = (1.0 + perts_w(i, k, j) * pert_deng_t) * w(i, k, j)
729 end subroutine Add_multi_perturb_shcu_perturbations
731 subroutine Remove_multi_perturb_shcu_perturbations (perts_qvapor, perts_qcloud, &
732 perts_th, perts_w, pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w, &
733 th, qv_curr, qc_curr, w, its, ite, jts, jte, ims, ime, jms, jme, kms, kme, &
738 integer, intent(in) :: its, ite, jts, jte, ims, ime, jms, jme, kms, kme, kts, kte
739 real, intent(in) :: pert_deng_qv, pert_deng_qc, pert_deng_t, pert_deng_w
740 real, dimension(ims:ime, kms:kme, jms:jme), optional, intent (in) :: perts_qvapor, &
741 perts_qcloud, perts_th, perts_w
742 real, dimension(ims:ime, kms:kme, jms:jme), intent (inout) :: th, w, qv_curr, qc_curr
750 qc_curr(i, k, j) = max (qc_curr(i, k, j) / (1.0 + perts_qcloud(i, k, j) * pert_deng_qc), QCLOUD_MIN)
751 qv_curr(i, k, j) = max (qv_curr(i, k, j) / (1.0 + perts_qvapor(i, k, j) * pert_deng_qv), QVAPOR_MIN)
752 th(i, k, j) = th(i, k, j) / (1.0 + perts_th(i, k, j) * pert_deng_t)
760 w(i, k, j) = w(i, k, j) / (1.0 + perts_w(i, k, j) * pert_deng_w)
765 end subroutine Remove_multi_perturb_shcu_perturbations
767 END MODULE module_shallowcu_driver