[2759] | 1 | !WRF:MEDIATION_LAYER:PHYSICS |
---|
| 2 | ! |
---|
| 3 | |
---|
| 4 | MODULE module_cumulus_driver |
---|
| 5 | CONTAINS |
---|
| 6 | SUBROUTINE cumulus_driver(grid & |
---|
| 7 | ! Order dependent args for domain, mem, and tile dims |
---|
| 8 | ,ids,ide, jds,jde, kds,kde & |
---|
| 9 | ,ims,ime, jms,jme, kms,kme & |
---|
| 10 | ,ips,ipe, jps,jpe, kps,kpe & |
---|
| 11 | ,i_start,i_end,j_start,j_end,kts,kte,num_tiles & |
---|
| 12 | ! Order independent args (use VAR= in call) |
---|
| 13 | ! --Prognostic |
---|
| 14 | ,u,v,th,t,w & |
---|
| 15 | ,p,pi,rho & |
---|
| 16 | ! --Other arguments |
---|
| 17 | ,itimestep,dt,dx,cudt,curr_secs,adapt_step_flag & |
---|
| 18 | ,rainc,raincv,pratec,nca & |
---|
| 19 | ,dz8w,p8w,forcet,forceq & |
---|
| 20 | ,w0avg,stepcu,gsw & |
---|
| 21 | ,cldefi,lowlyr,xland,cu_act_flag,warm_rain & |
---|
| 22 | ,htop,hbot,kpbl,ht & |
---|
| 23 | ,ensdim,maxiens,maxens,maxens2,maxens3 & |
---|
| 24 | ,periodic_x,periodic_y & |
---|
| 25 | ! Package selection variable |
---|
| 26 | ,cu_physics & |
---|
| 27 | ! Optional moisture tracers |
---|
| 28 | ,qv_curr, qc_curr, qr_curr & |
---|
| 29 | ,qi_curr, qs_curr, qg_curr & |
---|
| 30 | ,qv_prev, qc_prev, qr_prev & |
---|
| 31 | ,qi_prev, qs_prev, qg_prev & |
---|
| 32 | ! Optional arguments for GD scheme |
---|
| 33 | ,apr_gr,apr_w,apr_mc,apr_st,apr_as,apr_capma & |
---|
| 34 | ,apr_capme,apr_capmi,edt_out,clos_choice & |
---|
| 35 | ,mass_flux,xf_ens,pr_ens,cugd_avedx,imomentum & |
---|
| 36 | ,cugd_tten,cugd_qvten,cugd_qcten & |
---|
| 37 | ,cugd_ttens,cugd_qvtens & |
---|
| 38 | ,gd_cloud,gd_cloud2 & |
---|
| 39 | ! Optional moisture and other tendencies |
---|
| 40 | ,rqvcuten,rqccuten,rqrcuten & |
---|
| 41 | ,rqicuten,rqscuten,rqgcuten & |
---|
| 42 | ,rqvblten,rqvften & |
---|
| 43 | ,rthcuten,rthraten,rthblten,rthften & |
---|
| 44 | ! Optional moisture tracer flags |
---|
| 45 | ,f_qv,f_qc,f_qr & |
---|
| 46 | ,f_qi,f_qs,f_qg & |
---|
| 47 | ) |
---|
| 48 | !---------------------------------------------------------------------- |
---|
| 49 | USE module_model_constants |
---|
| 50 | USE module_state_description, ONLY: KFSCHEME,BMJSCHEME & |
---|
| 51 | ,KFETASCHEME,GDSCHEME & |
---|
| 52 | ,G3SCHEME & |
---|
| 53 | ,SASSCHEME |
---|
| 54 | |
---|
| 55 | ! *** add new modules of schemes here |
---|
| 56 | |
---|
| 57 | USE module_cu_kf |
---|
| 58 | USE module_cu_bmj |
---|
| 59 | USE module_dm |
---|
| 60 | USE module_domain, ONLY: domain |
---|
| 61 | USE module_cu_kfeta |
---|
| 62 | USE module_cu_gd, ONLY : GRELLDRV |
---|
| 63 | USE module_cu_g3, ONLY : G3DRV,CONV_GRELL_SPREAD3D |
---|
| 64 | USE module_cu_sas |
---|
| 65 | |
---|
| 66 | ! This driver calls subroutines for the cumulus parameterizations. |
---|
| 67 | ! |
---|
| 68 | ! 1. Kain & Fritsch (1993) |
---|
| 69 | ! 2. Betts-Miller-Janjic (Janjic, 1994) |
---|
| 70 | ! |
---|
| 71 | !---------------------------------------------------------------------- |
---|
| 72 | IMPLICIT NONE |
---|
| 73 | !====================================================================== |
---|
| 74 | ! Grid structure in physics part of WRF |
---|
| 75 | !---------------------------------------------------------------------- |
---|
| 76 | ! The horizontal velocities used in the physics are unstaggered |
---|
| 77 | ! relative to temperature/moisture variables. All predicted |
---|
| 78 | ! variables are carried at half levels except w, which is at full |
---|
| 79 | ! levels. Some arrays with names (*8w) are at w (full) levels. |
---|
| 80 | ! |
---|
| 81 | !---------------------------------------------------------------------- |
---|
| 82 | ! In WRF, kms (smallest number) is the bottom level and kme (largest |
---|
| 83 | ! number) is the top level. In your scheme, if 1 is at the top level, |
---|
| 84 | ! then you have to reverse the order in the k direction. |
---|
| 85 | ! |
---|
| 86 | ! kme - half level (no data at this level) |
---|
| 87 | ! kme ----- full level |
---|
| 88 | ! kme-1 - half level |
---|
| 89 | ! kme-1 ----- full level |
---|
| 90 | ! . |
---|
| 91 | ! . |
---|
| 92 | ! . |
---|
| 93 | ! kms+2 - half level |
---|
| 94 | ! kms+2 ----- full level |
---|
| 95 | ! kms+1 - half level |
---|
| 96 | ! kms+1 ----- full level |
---|
| 97 | ! kms - half level |
---|
| 98 | ! kms ----- full level |
---|
| 99 | ! |
---|
| 100 | !====================================================================== |
---|
| 101 | ! Definitions |
---|
| 102 | !----------- |
---|
| 103 | ! Rho_d dry density (kg/m^3) |
---|
| 104 | ! Theta_m moist potential temperature (K) |
---|
| 105 | ! Qv water vapor mixing ratio (kg/kg) |
---|
| 106 | ! Qc cloud water mixing ratio (kg/kg) |
---|
| 107 | ! Qr rain water mixing ratio (kg/kg) |
---|
| 108 | ! Qi cloud ice mixing ratio (kg/kg) |
---|
| 109 | ! Qs snow mixing ratio (kg/kg) |
---|
| 110 | !----------------------------------------------------------------- |
---|
| 111 | !-- DT time step (second) |
---|
| 112 | !-- CUDT cumulus time step (minute) |
---|
| 113 | !-- curr_secs current forecast time (seconds) |
---|
| 114 | !-- itimestep number of time step (integer) |
---|
| 115 | !-- DX horizontal space interval (m) |
---|
| 116 | !-- rr dry air density (kg/m^3) |
---|
| 117 | ! |
---|
| 118 | !-- RTHCUTEN Theta tendency due to |
---|
| 119 | ! cumulus scheme precipitation (K/s) |
---|
| 120 | !-- RQVCUTEN Qv tendency due to |
---|
| 121 | ! cumulus scheme precipitation (kg/kg/s) |
---|
| 122 | !-- RQRCUTEN Qr tendency due to |
---|
| 123 | ! cumulus scheme precipitation (kg/kg/s) |
---|
| 124 | !-- RQCCUTEN Qc tendency due to |
---|
| 125 | ! cumulus scheme precipitation (kg/kg/s) |
---|
| 126 | !-- RQSCUTEN Qs tendency due to |
---|
| 127 | ! cumulus scheme precipitation (kg/kg/s) |
---|
| 128 | !-- RQICUTEN Qi tendency due to |
---|
| 129 | ! cumulus scheme precipitation (kg/kg/s) |
---|
| 130 | ! |
---|
| 131 | !-- RAINC accumulated total cumulus scheme precipitation (mm) |
---|
| 132 | !-- RAINCV cumulus scheme precipitation (mm) |
---|
| 133 | !-- PRATEC precipitiation rate from cumulus scheme (mm/s) |
---|
| 134 | !-- NCA counter of the cloud relaxation |
---|
| 135 | ! time in KF cumulus scheme (integer) |
---|
| 136 | !-- u_phy u-velocity interpolated to theta points (m/s) |
---|
| 137 | !-- v_phy v-velocity interpolated to theta points (m/s) |
---|
| 138 | !-- th_phy potential temperature (K) |
---|
| 139 | !-- t_phy temperature (K) |
---|
| 140 | !-- w vertical velocity (m/s) |
---|
| 141 | !-- moist moisture array (4D - last index is species) (kg/kg) |
---|
| 142 | !-- dz8w dz between full levels (m) |
---|
| 143 | !-- p8w pressure at full levels (Pa) |
---|
| 144 | !-- p_phy pressure (Pa) |
---|
| 145 | !-- pi_phy exner function (dimensionless) |
---|
| 146 | ! points (dimensionless) |
---|
| 147 | !-- RTHRATEN radiative temp forcing for Grell-Devenyi scheme |
---|
| 148 | !-- RTHBLTEN PBL temp forcing for Grell-Devenyi scheme |
---|
| 149 | !-- RQVBLTEN PBL moisture forcing for Grell-Devenyi scheme |
---|
| 150 | !-- RTHFTEN |
---|
| 151 | !-- RQVFTEN |
---|
| 152 | !-- MASS_FLUX |
---|
| 153 | !-- XF_ENS |
---|
| 154 | !-- PR_ENS |
---|
| 155 | !-- warm_rain |
---|
| 156 | !-- CU_ACT_FLAG |
---|
| 157 | !-- W0AVG average vertical velocity, (for KF scheme) (m/s) |
---|
| 158 | !-- rho density (kg/m^3) |
---|
| 159 | !-- CLDEFI precipitation efficiency (for BMJ scheme) (dimensionless) |
---|
| 160 | !-- STEPCU # of fundamental timesteps between convection calls |
---|
| 161 | !-- XLAND land-sea mask (1.0 for land; 2.0 for water) |
---|
| 162 | !-- LOWLYR index of lowest model layer above the ground |
---|
| 163 | !-- XLV0 latent heat of vaporization constant |
---|
| 164 | ! used in temperature dependent formula (J/kg) |
---|
| 165 | !-- XLV1 latent heat of vaporization constant |
---|
| 166 | ! used in temperature dependent formula (J/kg/K) |
---|
| 167 | !-- XLS0 latent heat of sublimation constant |
---|
| 168 | ! used in temperature dependent formula (J/kg) |
---|
| 169 | !-- XLS1 latent heat of sublimation constant |
---|
| 170 | ! used in temperature dependent formula (J/kg/K) |
---|
| 171 | !-- R_d gas constant for dry air ( 287. J/kg/K) |
---|
| 172 | !-- R_v gas constant for water vapor (461 J/k/kg) |
---|
| 173 | !-- Cp specific heat at constant pressure (1004 J/k/kg) |
---|
| 174 | !-- rvovrd R_v divided by R_d (dimensionless) |
---|
| 175 | !-- G acceleration due to gravity (m/s^2) |
---|
| 176 | !-- EP_1 constant for virtual temperature |
---|
| 177 | ! (R_v/R_d - 1) (dimensionless) |
---|
| 178 | !-- pi_phy the exner function, (p/p0)**(R/Cp) (none unit) |
---|
| 179 | !-- ids start index for i in domain |
---|
| 180 | !-- ide end index for i in domain |
---|
| 181 | !-- jds start index for j in domain |
---|
| 182 | !-- jde end index for j in domain |
---|
| 183 | !-- kds start index for k in domain |
---|
| 184 | !-- kde end index for k in domain |
---|
| 185 | !-- ims start index for i in memory |
---|
| 186 | !-- ime end index for i in memory |
---|
| 187 | !-- jms start index for j in memory |
---|
| 188 | !-- jme end index for j in memory |
---|
| 189 | !-- kms start index for k in memory |
---|
| 190 | !-- kme end index for k in memory |
---|
| 191 | !-- i_start start indices for i in tile |
---|
| 192 | !-- i_end end indices for i in tile |
---|
| 193 | !-- j_start start indices for j in tile |
---|
| 194 | !-- j_end end indices for j in tile |
---|
| 195 | !-- kts start index for k in tile |
---|
| 196 | !-- kte end index for k in tile |
---|
| 197 | !-- num_tiles number of tiles |
---|
| 198 | !-- HBOT index of lowest model layer with convection |
---|
| 199 | !-- HTOP index of highest model layer with convection |
---|
| 200 | !-- LBOT index of lowest model layer with convection |
---|
| 201 | !-- LTOP index of highest model layer with convection |
---|
| 202 | !-- KPBL layer index of the PBL |
---|
| 203 | !-- periodic_x T/F this is using periodic lateral boundaries in the X direction |
---|
| 204 | !-- periodic_y T/F this is using periodic lateral boundaries in the Y-direction |
---|
| 205 | ! |
---|
| 206 | !====================================================================== |
---|
| 207 | |
---|
| 208 | INTEGER, INTENT(IN ) :: & |
---|
| 209 | ids,ide, jds,jde, kds,kde, & |
---|
| 210 | ims,ime, jms,jme, kms,kme, & |
---|
| 211 | kts,kte, & |
---|
| 212 | itimestep, num_tiles |
---|
| 213 | LOGICAL periodic_x, periodic_y |
---|
| 214 | TYPE(domain) , INTENT(INOUT) :: grid |
---|
| 215 | INTEGER, DIMENSION(num_tiles), INTENT(IN) :: & |
---|
| 216 | & i_start,i_end,j_start,j_end |
---|
| 217 | |
---|
| 218 | INTEGER, INTENT(IN ) :: & |
---|
| 219 | ensdim,maxiens,maxens,maxens2,maxens3 |
---|
| 220 | |
---|
| 221 | INTEGER, OPTIONAL, INTENT(IN ) :: & |
---|
| 222 | cugd_avedx,clos_choice |
---|
| 223 | |
---|
| 224 | INTEGER, INTENT(IN ) :: cu_physics |
---|
| 225 | INTEGER, INTENT(IN ) :: STEPCU |
---|
| 226 | LOGICAL, INTENT(IN ) :: warm_rain |
---|
| 227 | |
---|
| 228 | INTEGER,DIMENSION( ims:ime, jms:jme ), & |
---|
| 229 | INTENT(IN ) :: LOWLYR |
---|
| 230 | |
---|
| 231 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), & |
---|
| 232 | INTENT(IN ) :: & |
---|
| 233 | dz8w & |
---|
| 234 | , p8w & |
---|
| 235 | , p & |
---|
| 236 | , pi & |
---|
| 237 | , u & |
---|
| 238 | , v & |
---|
| 239 | , th & |
---|
| 240 | , t & |
---|
| 241 | , rho & |
---|
| 242 | , w |
---|
| 243 | |
---|
| 244 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), & |
---|
| 245 | INTENT(INOUT) :: & |
---|
| 246 | W0AVG |
---|
| 247 | |
---|
| 248 | REAL, DIMENSION( ims:ime , jms:jme ), INTENT(IN) :: & |
---|
| 249 | GSW,HT,XLAND |
---|
| 250 | |
---|
| 251 | REAL, DIMENSION( ims:ime , jms:jme ), & |
---|
| 252 | INTENT(INOUT) :: RAINC & |
---|
| 253 | , RAINCV & |
---|
| 254 | , NCA & |
---|
| 255 | , HTOP & |
---|
| 256 | , HBOT & |
---|
| 257 | , CLDEFI |
---|
| 258 | |
---|
| 259 | |
---|
| 260 | REAL, DIMENSION( ims:ime , jms:jme ),INTENT(INOUT),OPTIONAL :: & |
---|
| 261 | PRATEC |
---|
| 262 | REAL, DIMENSION( ims:ime , jms:jme ) :: tmppratec |
---|
| 263 | |
---|
| 264 | INTEGER, DIMENSION( ims:ime , jms:jme ), & |
---|
| 265 | INTENT(IN) :: KPBL |
---|
| 266 | |
---|
| 267 | |
---|
| 268 | LOGICAL, DIMENSION( ims:ime , jms:jme ), & |
---|
| 269 | INTENT(INOUT) :: CU_ACT_FLAG |
---|
| 270 | |
---|
| 271 | REAL, INTENT(IN ) :: DT, DX |
---|
| 272 | INTEGER, INTENT(IN ),OPTIONAL :: & |
---|
| 273 | ips,ipe, jps,jpe, kps,kpe,imomentum |
---|
| 274 | REAL, INTENT(IN ),OPTIONAL :: CUDT |
---|
| 275 | REAL, INTENT(IN ),OPTIONAL :: CURR_SECS |
---|
| 276 | LOGICAL,INTENT(IN ),OPTIONAL :: adapt_step_flag |
---|
| 277 | REAL :: cudt_pass, curr_secs_pass |
---|
| 278 | LOGICAL :: adapt_step_flag_pass |
---|
| 279 | |
---|
| 280 | ! |
---|
| 281 | ! optional arguments |
---|
| 282 | ! |
---|
| 283 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), & |
---|
| 284 | OPTIONAL, INTENT(INOUT) :: & |
---|
| 285 | ! optional moisture tracers |
---|
| 286 | ! 2 time levels; if only one then use CURR |
---|
| 287 | qv_curr, qc_curr, qr_curr & |
---|
| 288 | ,qi_curr, qs_curr, qg_curr & |
---|
| 289 | ,qv_prev, qc_prev, qr_prev & |
---|
| 290 | ,qi_prev, qs_prev, qg_prev & |
---|
| 291 | ! optional moisture and other tendencies |
---|
| 292 | ,rqvcuten,rqccuten,rqrcuten & |
---|
| 293 | ,rqicuten,rqscuten,rqgcuten & |
---|
| 294 | ,rqvblten,rqvften & |
---|
| 295 | ,rthraten,rthblten & |
---|
| 296 | ,cugd_tten,cugd_qvten,cugd_qcten & |
---|
| 297 | ,cugd_ttens,cugd_qvtens & |
---|
| 298 | , forcet & |
---|
| 299 | , forceq & |
---|
| 300 | ,rthften,rthcuten |
---|
| 301 | |
---|
| 302 | REAL, DIMENSION( ims:ime , jms:jme ), & |
---|
| 303 | OPTIONAL, & |
---|
| 304 | INTENT(INOUT) :: & |
---|
| 305 | apr_gr,apr_w,apr_mc,apr_st,apr_as,apr_capma & |
---|
| 306 | ,apr_capme,apr_capmi,edt_out & |
---|
| 307 | , MASS_FLUX |
---|
| 308 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), & |
---|
| 309 | OPTIONAL, INTENT(INOUT) :: & |
---|
| 310 | GD_CLOUD,GD_CLOUD2 |
---|
| 311 | REAL, DIMENSION( ims:ime , jms:jme , 1:ensdim ), & |
---|
| 312 | OPTIONAL, & |
---|
| 313 | INTENT(INOUT) :: XF_ENS, PR_ENS |
---|
| 314 | |
---|
| 315 | ! |
---|
| 316 | ! Flags relating to the optional tendency arrays declared above |
---|
| 317 | ! Models that carry the optional tendencies will provdide the |
---|
| 318 | ! optional arguments at compile time; these flags all the model |
---|
| 319 | ! to determine at run-time whether a particular tracer is in |
---|
| 320 | ! use or not. |
---|
| 321 | ! |
---|
| 322 | LOGICAL, INTENT(IN), OPTIONAL :: & |
---|
| 323 | f_qv & |
---|
| 324 | ,f_qc & |
---|
| 325 | ,f_qr & |
---|
| 326 | ,f_qi & |
---|
| 327 | ,f_qs & |
---|
| 328 | ,f_qg |
---|
| 329 | |
---|
| 330 | |
---|
| 331 | ! LOCAL VAR |
---|
| 332 | |
---|
| 333 | INTEGER :: i,j,k,its,ite,jts,jte,ij |
---|
| 334 | |
---|
| 335 | !----------------------------------------------------------------- |
---|
| 336 | |
---|
| 337 | if (.not. PRESENT(CURR_SECS)) then |
---|
| 338 | curr_secs_pass = -1 |
---|
| 339 | else |
---|
| 340 | curr_secs_pass = curr_secs |
---|
| 341 | endif |
---|
| 342 | |
---|
| 343 | if (.not. PRESENT(CUDT)) then |
---|
| 344 | cudt_pass = -1 |
---|
| 345 | else |
---|
| 346 | cudt_pass = cudt |
---|
| 347 | endif |
---|
| 348 | |
---|
| 349 | if (.not. PRESENT(adapt_step_flag)) then |
---|
| 350 | adapt_step_flag_pass = .false. |
---|
| 351 | else |
---|
| 352 | adapt_step_flag_pass = adapt_step_flag |
---|
| 353 | endif |
---|
| 354 | |
---|
| 355 | ! Initialize tmppratec to pratec |
---|
| 356 | |
---|
| 357 | if ( PRESENT ( pratec ) ) then |
---|
| 358 | tmppratec(:,:) = pratec(:,:) |
---|
| 359 | else |
---|
| 360 | tmppratec(:,:) = 0. |
---|
| 361 | end if |
---|
| 362 | |
---|
| 363 | |
---|
| 364 | IF (cu_physics .eq. 0) return |
---|
| 365 | #if ( EM_CORE == 1 ) |
---|
| 366 | if(cu_physics .eq. 5 ) then |
---|
| 367 | !$OMP PARALLEL DO & |
---|
| 368 | !$OMP PRIVATE ( ij,i,j,k,its,ite,jts,jte ) |
---|
| 369 | |
---|
| 370 | DO ij = 1 , num_tiles |
---|
| 371 | its = i_start(ij) |
---|
| 372 | ite = i_end(ij) |
---|
| 373 | jts = j_start(ij) |
---|
| 374 | jte = j_end(ij) |
---|
| 375 | do j=jts,min(jte,jde-1) |
---|
| 376 | do k=kts,kte |
---|
| 377 | do i=its,min(ite,ide-1) |
---|
| 378 | RTHFTEN(i,k,j)=(RTHFTEN(i,k,j)+RTHRATEN(i,k,j) & |
---|
| 379 | +RTHBLTEN(i,k,j))*pi(i,k,j) |
---|
| 380 | RQVFTEN(i,k,j)=RQVFTEN(i,k,j)+RQVBLTEN(i,k,j) |
---|
| 381 | enddo |
---|
| 382 | enddo |
---|
| 383 | enddo |
---|
| 384 | ENDDO |
---|
| 385 | !$OMP END PARALLEL DO |
---|
| 386 | #ifdef DM_PARALLEL |
---|
| 387 | #include "HALO_CUP_G3_IN.inc" |
---|
| 388 | #endif |
---|
| 389 | endif |
---|
| 390 | #endif |
---|
| 391 | |
---|
| 392 | ! DON'T JUDGE TIME STEP HERE, SINCE KF NEEDS ACCUMULATED W FIELD. |
---|
| 393 | ! DO IT INSIDE THE INDIVIDUAL CUMULUS SCHEME |
---|
| 394 | |
---|
| 395 | ! SET START AND END POINTS FOR TILES |
---|
| 396 | !$OMP PARALLEL DO & |
---|
| 397 | !$OMP PRIVATE ( ij ,its,ite,jts,jte, i,j,k) |
---|
| 398 | |
---|
| 399 | DO ij = 1 , num_tiles |
---|
| 400 | its = i_start(ij) |
---|
| 401 | ite = i_end(ij) |
---|
| 402 | jts = j_start(ij) |
---|
| 403 | jte = j_end(ij) |
---|
| 404 | |
---|
| 405 | |
---|
| 406 | cps_select: SELECT CASE(cu_physics) |
---|
| 407 | |
---|
| 408 | CASE (KFSCHEME) |
---|
| 409 | CALL wrf_debug(100,'in kfcps') |
---|
| 410 | |
---|
| 411 | CALL KFCPS( & |
---|
| 412 | ! order independent arguments |
---|
| 413 | DT=dt ,KTAU=itimestep ,DX=dx , CUDT=cudt_pass & |
---|
| 414 | ,CURR_SECS=curr_secs_pass & |
---|
| 415 | ,ADAPT_STEP_FLAG=adapt_step_flag_pass & |
---|
| 416 | ,RHO=rho & |
---|
| 417 | ,U=u ,V=v ,TH=th ,T=t ,W=w & |
---|
| 418 | ,PCPS=p ,PI=pi & |
---|
| 419 | ,XLV0=xlv0 ,XLV1=xlv1 ,XLS0=xls0 ,XLS1=xls1 & |
---|
| 420 | ,RAINCV=raincv, PRATEC=tmppratec, NCA=nca & |
---|
| 421 | ,DZ8W=dz8w & |
---|
| 422 | ,W0AVG=w0avg & |
---|
| 423 | ,CP=cp ,R=r_d ,G=g ,EP1=ep_1 ,EP2=ep_2 & |
---|
| 424 | ,SVP1=svp1 ,SVP2=svp2 ,SVP3=svp3 ,SVPT0=svpt0 & |
---|
| 425 | ,STEPCU=stepcu & |
---|
| 426 | ,CU_ACT_FLAG=cu_act_flag & |
---|
| 427 | ,WARM_RAIN=warm_rain & |
---|
| 428 | ,QV=qv_curr & |
---|
| 429 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 430 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 431 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 432 | ! optionals |
---|
| 433 | ,RTHCUTEN=rthcuten ,RQVCUTEN=rqvcuten & |
---|
| 434 | ,RQCCUTEN=rqccuten ,RQRCUTEN=rqrcuten & |
---|
| 435 | ,RQICUTEN=rqicuten ,RQSCUTEN=rqscuten & |
---|
| 436 | ,F_QV=f_qv,F_QC=f_qc,F_QR=f_qr & |
---|
| 437 | ,F_QI=f_qi,F_QS=f_qs & |
---|
| 438 | ) |
---|
| 439 | |
---|
| 440 | CASE (BMJSCHEME) |
---|
| 441 | CALL wrf_debug(100,'in bmj_cps') |
---|
| 442 | CALL BMJDRV( & |
---|
| 443 | TH=th,T=T ,RAINCV=raincv, PRATEC=tmppratec & |
---|
| 444 | ,RHO=rho & |
---|
| 445 | ,DT=dt ,ITIMESTEP=itimestep ,STEPCU=stepcu & |
---|
| 446 | ,CUDT=cudt_pass & |
---|
| 447 | ,CURR_SECS=curr_secs_pass & |
---|
| 448 | ,ADAPT_STEP_FLAG=adapt_step_flag_pass & |
---|
| 449 | ,CUTOP=htop, CUBOT=hbot, KPBL=kpbl & |
---|
| 450 | ,DZ8W=dz8w ,PINT=p8w, PMID=p, PI=pi & |
---|
| 451 | ,CP=cp ,R=r_d ,ELWV=xlv ,ELIV=xls ,G=g & |
---|
| 452 | ,TFRZ=svpt0 ,D608=ep_1 ,CLDEFI=cldefi & |
---|
| 453 | ,LOWLYR=lowlyr ,XLAND=xland & |
---|
| 454 | ,CU_ACT_FLAG=cu_act_flag & |
---|
| 455 | ,QV=qv_curr & |
---|
| 456 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 457 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 458 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 459 | ! optionals |
---|
| 460 | ,RTHCUTEN=rthcuten ,RQVCUTEN=rqvcuten & |
---|
| 461 | ) |
---|
| 462 | |
---|
| 463 | CASE (KFETASCHEME) |
---|
| 464 | CALL wrf_debug(100,'in kf_eta_cps') |
---|
| 465 | CALL KF_ETA_CPS( & |
---|
| 466 | U=u ,V=v ,TH=th ,T=t ,W=w ,RHO=rho & |
---|
| 467 | ,CUDT=cudt_pass & |
---|
| 468 | ,CURR_SECS=curr_secs_pass & |
---|
| 469 | ,ADAPT_STEP_FLAG=adapt_step_flag_pass & |
---|
| 470 | ,RAINCV=raincv, PRATEC=tmppratec, NCA=nca & |
---|
| 471 | ,DZ8W=dz8w & |
---|
| 472 | ,PCPS=p, PI=pi ,W0AVG=W0AVG & |
---|
| 473 | ,CUTOP=HTOP,CUBOT=HBOT & |
---|
| 474 | ,XLV0=XLV0 ,XLV1=XLV1 ,XLS0=XLS0 ,XLS1=XLS1 & |
---|
| 475 | ,CP=CP ,R=R_d ,G=G ,EP1=EP_1 ,EP2=EP_2 & |
---|
| 476 | ,SVP1=SVP1 ,SVP2=SVP2 ,SVP3=SVP3 ,SVPT0=SVPT0 & |
---|
| 477 | ,DT=dt ,KTAU=itimestep ,DX=dx & |
---|
| 478 | ,STEPCU=stepcu & |
---|
| 479 | ,CU_ACT_FLAG=cu_act_flag ,WARM_RAIN=warm_rain & |
---|
| 480 | ,QV=qv_curr & |
---|
| 481 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 482 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 483 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 484 | ! optionals |
---|
| 485 | ,RTHCUTEN=rthcuten & |
---|
| 486 | ,RQVCUTEN=rqvcuten ,RQCCUTEN=rqccuten & |
---|
| 487 | ,RQRCUTEN=rqrcuten ,RQICUTEN=rqicuten & |
---|
| 488 | ,RQSCUTEN=rqscuten & |
---|
| 489 | ,F_QV=f_qv,F_QC=f_qc,F_QR=f_qr & |
---|
| 490 | ,F_QI=f_qi,F_QS=f_qs & |
---|
| 491 | ) |
---|
| 492 | |
---|
| 493 | CASE (GDSCHEME) |
---|
| 494 | CALL wrf_debug(100,'in grelldrv') |
---|
| 495 | CALL GRELLDRV( & |
---|
| 496 | DT=dt, ITIMESTEP=itimestep, DX=dx & |
---|
| 497 | ,U=u,V=v,T=t,W=w ,RHO=rho & |
---|
| 498 | ,P=p,PI=pi ,Q=qv_curr ,RAINCV=raincv & |
---|
| 499 | ,DZ8W=dz8w,P8W=p8w,XLV=xlv,CP=cp,G=g,R_V=r_v & |
---|
| 500 | ,PRATEC=tmppratec & |
---|
| 501 | ,APR_GR=apr_gr,APR_W=apr_w,APR_MC=apr_mc & |
---|
| 502 | ,APR_ST=apr_st,APR_AS=apr_as & |
---|
| 503 | ,APR_CAPMA=apr_capma,APR_CAPME=apr_capme & |
---|
| 504 | ,APR_CAPMI=apr_capmi,MASS_FLUX=mass_flux & |
---|
| 505 | ,XF_ENS=xf_ens,PR_ENS=pr_ens,HT=ht & |
---|
| 506 | ,xland=xland,gsw=gsw & |
---|
| 507 | ,GDC=gd_cloud,GDC2=gd_cloud2 & |
---|
| 508 | ,ENSDIM=ensdim,MAXIENS=maxiens,MAXENS=maxens & |
---|
| 509 | ,MAXENS2=maxens2,MAXENS3=maxens3 & |
---|
| 510 | ,STEPCU=STEPCU,htop=htop,hbot=hbot & |
---|
| 511 | ,CU_ACT_FLAG=CU_ACT_FLAG,warm_rain=warm_rain & |
---|
| 512 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 513 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 514 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 515 | ,PERIODIC_X=periodic_x,PERIODIC_Y=periodic_y & |
---|
| 516 | ! optionals |
---|
| 517 | #if (NMM_CORE == 1 ) |
---|
| 518 | ,RTHCUTEN=RTHCUTEN ,RTHFTEN=forcet & |
---|
| 519 | ,RQICUTEN=RQICUTEN ,RQVFTEN=forceq & |
---|
| 520 | #else |
---|
| 521 | ,RTHCUTEN=RTHCUTEN ,RTHFTEN=RTHFTEN & |
---|
| 522 | ,RQICUTEN=RQICUTEN ,RQVFTEN=RQVFTEN & |
---|
| 523 | #endif |
---|
| 524 | ,RTHRATEN=RTHRATEN,RTHBLTEN=RTHBLTEN & |
---|
| 525 | ,RQVCUTEN=RQVCUTEN,RQCCUTEN=RQCCUTEN & |
---|
| 526 | ,RQVBLTEN=RQVBLTEN & |
---|
| 527 | ,F_QV=f_qv,F_QC=f_qc,F_QR=f_qr & |
---|
| 528 | ,F_QI=f_qi,F_QS=f_qs & |
---|
| 529 | ) |
---|
| 530 | CASE (SASSCHEME) |
---|
| 531 | |
---|
| 532 | IF ( adapt_step_flag_pass ) THEN |
---|
| 533 | WRITE( wrf_err_message , * ) 'The SAS cumulus option will not work properly with an adaptive time step' |
---|
| 534 | CALL wrf_error_fatal ( wrf_err_message ) |
---|
| 535 | END IF |
---|
| 536 | CALL wrf_debug(100,'in cu_sas') |
---|
| 537 | CALL CU_SAS( & |
---|
| 538 | DT=dt,ITIMESTEP=itimestep,STEPCU=STEPCU & |
---|
| 539 | ,RAINCV=RAINCV,PRATEC=tmpPRATEC,HTOP=HTOP,HBOT=HBOT & |
---|
| 540 | ,U3D=u,V3D=v,W=w,T3D=t,PI3D=pi,RHO3D=rho & |
---|
| 541 | ,QV3D=QV_CURR,QC3D=QC_CURR,QI3D=QI_CURR & |
---|
| 542 | ,DZ8W=dz8w,PCPS=p,P8W=p8w,XLAND=XLAND & |
---|
| 543 | ,CU_ACT_FLAG=CU_ACT_FLAG & |
---|
| 544 | ,CUDT=cudt_pass & |
---|
| 545 | ,CURR_SECS=curr_secs_pass & |
---|
| 546 | ,ADAPT_STEP_FLAG=adapt_step_flag_pass & |
---|
| 547 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 548 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 549 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 550 | ! optionals |
---|
| 551 | ,RTHCUTEN=RTHCUTEN,RQVCUTEN=RQVCUTEN & |
---|
| 552 | ,RQCCUTEN=RQCCUTEN,RQICUTEN=RQICUTEN & |
---|
| 553 | ,F_QV=f_qv,F_QC=f_qc,F_QR=f_qr & |
---|
| 554 | ,F_QI=f_qi,F_QS=f_qs & |
---|
| 555 | ) |
---|
| 556 | CASE (G3SCHEME) |
---|
| 557 | CALL wrf_debug(100,'in grelldrv') |
---|
| 558 | CALL G3DRV( & |
---|
| 559 | DT=dt, ITIMESTEP=itimestep, DX=dx & |
---|
| 560 | ,U=u,V=v,T=t,W=w ,RHO=rho & |
---|
| 561 | ,P=p,PI=pi,Q=qv_curr,RAINCV=raincv & |
---|
| 562 | ,DZ8W=dz8w ,P8W=p8w,XLV=xlv,CP=cp,G=g,R_V=r_v & |
---|
| 563 | ,APR_GR=apr_gr,APR_W=apr_w,APR_MC=apr_mc & |
---|
| 564 | ,APR_ST=apr_st,APR_AS=apr_as,PRATEC=tmppratec & |
---|
| 565 | ,APR_CAPMA=apr_capma,APR_CAPME=apr_capme & |
---|
| 566 | ,APR_CAPMI=apr_capmi,MASS_FLUX=mass_flux & |
---|
| 567 | ,XF_ENS=xf_ens,PR_ENS=pr_ens,HT=ht & |
---|
| 568 | ,xland=xland,gsw=gsw,edt_out=edt_out & |
---|
| 569 | ,GDC=gd_cloud,GDC2=gd_cloud2 & |
---|
| 570 | ,cugd_tten=cugd_tten,cugd_qvten=cugd_qvten & |
---|
| 571 | ,cugd_ttens=cugd_ttens,cugd_qvtens=cugd_qvtens & |
---|
| 572 | ,cugd_qcten=cugd_qcten,cugd_avedx=cugd_avedx & |
---|
| 573 | ,imomentum=imomentum & |
---|
| 574 | ,ENSDIM=ensdim,MAXIENS=maxiens,MAXENS=maxens & |
---|
| 575 | ,MAXENS2=maxens2,MAXENS3=maxens3,ichoice=clos_choice & |
---|
| 576 | ,STEPCU=STEPCU,htop=htop,hbot=hbot & |
---|
| 577 | ,CU_ACT_FLAG=CU_ACT_FLAG,warm_rain=warm_rain & |
---|
| 578 | ,IDS=ids,IDE=ide,JDS=jds,JDE=jde,KDS=kds,KDE=kde & |
---|
| 579 | ,IMS=ims,IME=ime,JMS=jms,JME=jme,KMS=kms,KME=kme & |
---|
| 580 | ,IPS=ips,IPE=ipe,JPS=jps,JPE=jpe,KPS=kps,KPE=kpe & |
---|
| 581 | ,ITS=its,ITE=ite,JTS=jts,JTE=jte,KTS=kts,KTE=kte & |
---|
| 582 | ,PERIODIC_X=periodic_x,PERIODIC_Y=periodic_y & |
---|
| 583 | ! optionals |
---|
| 584 | #if (NMM_CORE == 1 ) |
---|
| 585 | ,RTHCUTEN=RTHCUTEN ,RTHFTEN=forcet & |
---|
| 586 | ,RQICUTEN=RQICUTEN ,RQVFTEN=forceq & |
---|
| 587 | #else |
---|
| 588 | ,RTHCUTEN=RTHCUTEN ,RTHFTEN=RTHFTEN & |
---|
| 589 | ,RQICUTEN=RQICUTEN ,RQVFTEN=RQVFTEN & |
---|
| 590 | #endif |
---|
| 591 | ,RQVCUTEN=RQVCUTEN,RQCCUTEN=RQCCUTEN & |
---|
| 592 | ,F_QV=f_qv,F_QC=f_qc,F_QR=f_qr & |
---|
| 593 | ,F_QI=f_qi,F_QS=f_qs & |
---|
| 594 | ) |
---|
| 595 | |
---|
| 596 | CASE DEFAULT |
---|
| 597 | |
---|
| 598 | WRITE( wrf_err_message , * ) 'The cumulus option does not exist: cu_physics = ', cu_physics |
---|
| 599 | CALL wrf_error_fatal ( wrf_err_message ) |
---|
| 600 | |
---|
| 601 | END SELECT cps_select |
---|
| 602 | |
---|
| 603 | ENDDO |
---|
| 604 | !$OMP END PARALLEL DO |
---|
| 605 | #if ( EM_CORE == 1 ) |
---|
| 606 | IF(cu_physics .eq. 5 )then |
---|
| 607 | #ifdef DM_PARALLEL |
---|
| 608 | # include "HALO_CUP_G3_OUT.inc" |
---|
| 609 | #endif |
---|
| 610 | call conv_grell_spread3d(rthcuten=rthcuten,rqvcuten=rqvcuten & |
---|
| 611 | & ,rqccuten=rqccuten,raincv=raincv,cugd_avedx=cugd_avedx & |
---|
| 612 | & ,cugd_tten=cugd_tten,cugd_qvten=cugd_qvten,rqicuten=rqicuten & |
---|
| 613 | & ,cugd_ttens=cugd_ttens,cugd_qvtens=cugd_qvtens & |
---|
| 614 | & ,cugd_qcten=cugd_qcten,pi_phy=pi,moist_qv=qv_curr & |
---|
| 615 | & ,PRATEC=tmppratec,dt=dt,num_tiles=num_tiles & |
---|
| 616 | & ,imomentum=imomentum & |
---|
| 617 | & ,F_QV=F_QV,F_QC=F_QC,F_QR=F_QR,F_QI=F_QI,F_QS=F_QS & |
---|
| 618 | & ,ids=IDS,ide=IDE, jds=JDS,jde=JDE, kds=KDS,kde=KDE & |
---|
| 619 | & ,ips=IPS,ipe=IPE, jps=JPS,jpe=JPE, kps=KPS,kpe=KPE & |
---|
| 620 | & ,ims=IMS,ime=IME, jms=JMS,jme=JME, kms=KMS,kme=KME & |
---|
| 621 | & ,i_start=i_start,i_end=i_end & |
---|
| 622 | & ,j_start=j_start,j_end=j_end & |
---|
| 623 | & ,kts=kts, kte=kte) |
---|
| 624 | endif |
---|
| 625 | #endif |
---|
| 626 | |
---|
| 627 | ! |
---|
| 628 | ! Copy pratec back to output array, if necessary. |
---|
| 629 | ! |
---|
| 630 | if (PRESENT(PRATEC)) then |
---|
| 631 | pratec(:,:) = tmppratec(:,:) |
---|
| 632 | endif |
---|
| 633 | END SUBROUTINE cumulus_driver |
---|
| 634 | |
---|
| 635 | END MODULE module_cumulus_driver |
---|