Adjusting include paths for removal of redundant code
[WRF.git] / chem / chem_driver.F
blob4c8268df1bdafecd6cddbe0a7a5c92083a1edd98
2 !   WRF-chem V3.1 : Original version of chem_driver written by Georg Grell (ESRL/GSD)
3 !                   Further developments, bugfixes and improvements  by
4 !                   William Gustafson (PNNL), Marc Salzmann (GFDL), and Georg Grell
5 ! 10/12/2011 - Ravan Ahmadov (NOAA) updated to include the RACM_SOA_VBS option
6 ! 10/08/2014 - Kai Wang and Yang Zhang (NCSU) updated to include the CB05_MADE/SORGAM and CB05_MADE/VBS options
8 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
9     subroutine chem_driver ( grid , config_flags   &
11 #include "dummy_new_args.inc"
13                  )
14 !----------------------------------------------------------------------
15   USE module_domain , only : domain
16   USE module_configure
17   USE module_driver_constants
18   USE module_machine
19   USE module_tiles
20   USE module_dm
21   USE module_model_constants
22   USE module_state_description
23   USE module_data_radm2
24   USE module_data_sorgam
25   USE module_radm
26   USE module_dep_simple
27   USE module_bioemi_simple
28   USE module_phot_mad
29   USE module_phot_tuv,    only : tuv_timestep_init
30   USE module_ftuv_driver, only : ftuv_timestep_init
31   USE module_aerosols_sorgam
32   USE module_chem_utilities
33   USE module_gocart_so2so4
34   USE module_aer_opt_out,only: aer_opt_out
35   USE module_ctrans_grell
36   USE module_data_soa_vbs, only: ldrog_vbs
37   USE module_dust_load
38   USE module_chem_cup, only: chem_cup_driver !BSINGH - For WRFCuP scheme
39   USE module_dry_dep_driver
40   USE module_emissions_driver
41   USE module_input_tracer, only: set_tracer
42   USE module_wetscav_driver, only: wetscav_driver
43   USE module_wetdep_ls, only:wetdep_ls
44   USE module_uoc_dustwd ! Claudia, 3 April 2014  [mklose 03082015]
45   USE module_input_chem_data, only: last_chem_time, &
46 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
47                                      chem_dbg,      &
48 #endif
49                                      mozcart_lbc_set, &
50                                      bdy_chem_value_top_pv,PVS
51   USE module_chem_share,      only: get_last_gas
52   USE module_upper_bc_driver, only: upper_bc_driver
53   USE module_tropopause,      only: tropopause_driver
54   USE modal_aero_data, only: ntot_amode_cam_mam => ntot_amode
55   USE module_cam_support, only: gas_pcnst => gas_pcnst_modal_aero,gas_pcnst_pos => gas_pcnst_modal_aero_pos, &
56        pcnst =>pcnst_runtime, numgas_mam, cam_mam_aerosols
57   USE module_cu_camzm_driver, only: zm_conv_tend_2
58   USE module_cam_mam_gas_wetdep_driver, only: cam_mam_gas_wetdep_driver
59   USE module_trajectory, only: trajectory_dchm_tstep_init, trajectory_dchm_tstep_set
61   IMPLICIT NONE
64   !BSINGH(PNNL)- Lahey compiler forces to declare the following interface
65   interface
66      SUBROUTINE sum_pm_driver ( config_flags,                       &
67           alt, chem, h2oaj, h2oai,                                  &
68           pm2_5_dry, pm2_5_water, pm2_5_dry_ec, pm10,               &
69           tsoa,asoa,bsoa,                                             &
70           hoa_a01,hoa_a02,hoa_a03,hoa_a04,                          &
71           hoa_a05,hoa_a06,hoa_a07,hoa_a08,                          & !BSINGH Added 4 more bins for each species for SAPRC 8 bin version
72           bboa_a01,bboa_a02,bboa_a03,bboa_a04,                      &
73           bboa_a05,bboa_a06,bboa_a07,bboa_a08,                      &
74           soa_a01,soa_a02,soa_a03,soa_a04,                          &
75           soa_a05,soa_a06,soa_a07,soa_a08,                          &
76           bbsoa_a01,bbsoa_a02,bbsoa_a03,bbsoa_a04,                  &
77           bbsoa_a05,bbsoa_a06,bbsoa_a07,bbsoa_a08,                  &
78           hsoa_a01,hsoa_a02,hsoa_a03,hsoa_a04,                      &
79           hsoa_a05,hsoa_a06,hsoa_a07,hsoa_a08,                      &
80           biog_a01,biog_a02,biog_a03,biog_a04,                      &
81           biog_a05,biog_a06,biog_a07,biog_a08,                      &
82           asmpsoa_a01,asmpsoa_a02,asmpsoa_a03,asmpsoa_a04,          &
83           arosoa_a01,arosoa_a02,arosoa_a03,arosoa_a04,              &
84           arosoa_a05,arosoa_a06,arosoa_a07,arosoa_a08,              &
85           totoa_a01,totoa_a02,totoa_a03,totoa_a04,                  &
86           totoa_a05,totoa_a06,totoa_a07,totoa_a08,                  & 
87           hsoa_c,hsoa_o,bbsoa_c,bbsoa_o,                            &
88           biog_v1,biog_v2,biog_v3,biog_v4,                          &
89           ant_v1,ant_v2,ant_v3,ant_v4,                              &
90           smpa_v1,smpbb_v1,                                         &
91           !BSINGH - Added cw aerosols(for VBS)
92           hoa_cw01,    hoa_cw02,    hoa_cw03,    hoa_cw04,          &
93           hoa_cw05,    hoa_cw06,    hoa_cw07,    hoa_cw08,          &
94           bboa_cw01,   bboa_cw02,   bboa_cw03,   bboa_cw04,         &
95           bboa_cw05,   bboa_cw06,   bboa_cw07,   bboa_cw08,         &
96           soa_cw01,    soa_cw02,    soa_cw03,    soa_cw04,          &
97           soa_cw05,    soa_cw06,    soa_cw07,    soa_cw08,          &
98           bbsoa_cw01,  bbsoa_cw02,  bbsoa_cw03,  bbsoa_cw04,        &
99           bbsoa_cw05,  bbsoa_cw06,  bbsoa_cw07,  bbsoa_cw08,        &
100           biog_cw01,   biog_cw02,   biog_cw03,   biog_cw04,         &
101           biog_cw05,   biog_cw06,   biog_cw07,   biog_cw08,         &
102           hsoa_cw01,   hsoa_cw02,   hsoa_cw03,   hsoa_cw04,         &
103           hsoa_cw05,   hsoa_cw06,   hsoa_cw07,   hsoa_cw08,         &
104           arosoa_cw01, arosoa_cw02, arosoa_cw03, arosoa_cw04,       &
105           arosoa_cw05, arosoa_cw06, arosoa_cw07, arosoa_cw08,       &
106           totoa_cw01,  totoa_cw02,  totoa_cw03,  totoa_cw04,        &
107           totoa_cw05,  totoa_cw06,  totoa_cw07,  totoa_cw08,        &        
108           hsoa_cw_c,   hsoa_cw_o,   bbsoa_cw_c,  bbsoa_cw_o,        &
109           biog_cw_v1,                                               &
110           ant_cw_v1,                                                &
111           !BSINGH -ENDS
112           ids,ide, jds,jde, kds,kde,                                &
113           ims,ime, jms,jme, kms,kme,                                &
114           its,ite, jts,jte, kts,kte                                 )
115        
116        
117        USE module_configure
118        USE module_aerosols_sorgam, only: sum_pm_sorgam
119        USE module_mosaic_driver, only: sum_pm_mosaic,sum_pm_mosaic_vbs2,sum_pm_mosaic_vbs0,sum_vbs9,sum_vbs2,sum_vbs0
120        USE module_gocart_aerosols, only: sum_pm_gocart
121        USE module_aerosols_soa_vbs, only: sum_pm_soa_vbs
122        USE module_aerosols_sorgam_vbs, only: sum_pm_sorgam_vbs
123        
124        IMPLICIT NONE
125        
126        INTEGER,      INTENT(IN   )    ::                                   &
127             ids,ide, jds,jde, kds,kde,       &
128             ims,ime, jms,jme, kms,kme,       &
129             its,ite, jts,jte, kts,kte
130        
131        REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ),             &
132             INTENT(IN ) :: chem
133        
134        REAL, DIMENSION( ims:ime, kms:kme, jms:jme ),                       &
135             INTENT(IN ) :: alt
136        REAL, DIMENSION( ims:ime, kms:kme, jms:jme ),                       &
137             OPTIONAL,                                                     &
138             INTENT(IN ) :: h2oaj,h2oai
139        
140        REAL, DIMENSION( ims:ime, kms:kme, jms:jme ),                  &
141             OPTIONAL,                                                 &
142             INTENT(OUT) :: pm2_5_dry,pm2_5_water,pm2_5_dry_ec,pm10,   &
143             tsoa,asoa,bsoa,                                           &
144             hoa_a01,hoa_a02,hoa_a03,hoa_a04,                          &
145             hoa_a05,hoa_a06,hoa_a07,hoa_a08,                          &!BSINGH Added 4 more bins for each species for SAPRC 8 bin version             
146             bboa_a01,bboa_a02,bboa_a03,bboa_a04,                      &
147             bboa_a05,bboa_a06,bboa_a07,bboa_a08,                      &     
148             soa_a01,soa_a02,soa_a03,soa_a04,                          &
149             soa_a05,soa_a06,soa_a07,soa_a08,                          &        
150             bbsoa_a01,bbsoa_a02,bbsoa_a03,bbsoa_a04,                  &
151             bbsoa_a05,bbsoa_a06,bbsoa_a07,bbsoa_a08,                  &     
152             hsoa_a01,hsoa_a02,hsoa_a03,hsoa_a04,                      &
153             hsoa_a05,hsoa_a06,hsoa_a07,hsoa_a08,                      &       
154             biog_a01,biog_a02,biog_a03,biog_a04,                      &
155             biog_a05,biog_a06,biog_a07,biog_a08,                      &       
156             arosoa_a01,arosoa_a02,arosoa_a03,arosoa_a04,              &
157             arosoa_a05,arosoa_a06,arosoa_a07,arosoa_a08,              &             
158             totoa_a01,totoa_a02,totoa_a03,totoa_a04,                  &
159             totoa_a05,totoa_a06,totoa_a07,totoa_a08,                  &    
160             hsoa_c,hsoa_o,bbsoa_c,bbsoa_o,                            &
161             biog_v1,biog_v2,biog_v3,biog_v4,                          &
162             ant_v1,ant_v2,ant_v3,ant_v4,                              &
163             smpa_v1,                                                  &
164             smpbb_v1,                                                 &
165             asmpsoa_a01,asmpsoa_a02,asmpsoa_a03,asmpsoa_a04,          &!BSINGH - Not adding 5-8 bins for asmpsoa, as it is not req. for 8 bin SAPRC
166             !BSINGH - Added cw aerosols(for VBS)
167             hoa_cw01,    hoa_cw02,    hoa_cw03,    hoa_cw04,          &
168             hoa_cw05,    hoa_cw06,    hoa_cw07,    hoa_cw08,          &
169             bboa_cw01,   bboa_cw02,   bboa_cw03,   bboa_cw04,         &
170             bboa_cw05,   bboa_cw06,   bboa_cw07,   bboa_cw08,         &
171             soa_cw01,    soa_cw02,    soa_cw03,    soa_cw04,          &
172             soa_cw05,    soa_cw06,    soa_cw07,    soa_cw08,          &
173             bbsoa_cw01,  bbsoa_cw02,  bbsoa_cw03,  bbsoa_cw04,        &
174             bbsoa_cw05,  bbsoa_cw06,  bbsoa_cw07,  bbsoa_cw08,        &
175             biog_cw01,   biog_cw02,   biog_cw03,   biog_cw04,         &
176             biog_cw05,   biog_cw06,   biog_cw07,   biog_cw08,         &
177             hsoa_cw01,   hsoa_cw02,   hsoa_cw03,   hsoa_cw04,         &
178             hsoa_cw05,   hsoa_cw06,   hsoa_cw07,   hsoa_cw08,         &
179             arosoa_cw01, arosoa_cw02, arosoa_cw03, arosoa_cw04,       &
180             arosoa_cw05, arosoa_cw06, arosoa_cw07, arosoa_cw08,       &
181             totoa_cw01,  totoa_cw02,  totoa_cw03,  totoa_cw04,        &
182             totoa_cw05,  totoa_cw06,  totoa_cw07,  totoa_cw08,        &        
183             hsoa_cw_c,   hsoa_cw_o,   bbsoa_cw_c,  bbsoa_cw_o,        &
184             biog_cw_v1,                                               &
185             ant_cw_v1                                                 
186              !BSINGH -ENDS
187        
188        
189        TYPE(grid_config_rec_type),  INTENT(IN   )    :: config_flags
190        
191      end SUBROUTINE sum_pm_driver
192   end interface
193   
194   
195    !  Input data.
197    TYPE(domain) , TARGET          :: grid
198    !
199    !  Definitions of dummy arguments to solve
200 # include "dummy_new_decl.inc"
201 # define NO_I1_OLD
203    TYPE(grid_config_rec_type),  INTENT(IN   )    :: config_flags
205    INTEGER                     :: ids,ide, jds,jde, kds,kde,    &
206                                   ims,ime, jms,jme, kms,kme,    &
207                                   ips,ipe, jps,jpe, kps,kpe,    &
208                                   its,ite, jts,jte, kts,kte
209 ! ..
210 ! .. Local Scalars ..
211 #if( WRF_CHEM == 1 && WRF_KPP == 1 )
212    REAL, PARAMETER  ::  navgdro = 6.022e23   ! molecules/mol
213    REAL, PARAMETER  ::  mw_air = 28.97       !molecular weights
214 !dens2con air
215    REAL, PARAMETER ::          dens2con_a = 1.e-3    &! kg/m3 -> g/cm3
216                                * (1./mw_air)         &! -> mole/cm3
217                                * navgdro              ! -> molec/cm3
218 #endif
219    INTEGER :: stepave,i,j,k,l,numgas,nv,n, nr,ktau,k_start,k_end,idf,jdf,kdf
220    INTEGER :: ijulian
221 ! UoC dust scheme option
222    INTEGER :: imod   
223 !  REAL    :: convtrans_avglen_m
225 ! ................................................................
226       real, dimension(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33) ::vcsulf_old,vcso2_old,vch2o2_old
227       real, dimension(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,ldrog) ::vdrog3
229       real, dimension(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,ldrog_vbs) ::vdrog3_vbs
231       real, dimension(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33) ::n2o5_het 
232       REAL,DIMENSION(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33) ::              &
233                                                               p_phy,u_phy,v_phy                   &
234                                                              ,t_phy,dz8w,t8w,p8w                  &
235                                                              ,rho,rri,z_at_w,vvel,zmid,rh
236       REAL,DIMENSION(grid%sm31:grid%em31,grid%sm33:grid%em33) :: pbl_h
237       REAL,DIMENSION(grid%sm31:grid%em31,grid%sm33:grid%em33, 5) :: seasin,dustin
238       REAL,DIMENSION(grid%sm32:grid%em32-1) :: QL,TL
239       REAL,DIMENSION(grid%sm31:grid%em31,grid%sm33:grid%em33) :: REXNSFC,FACTRS                &
240                                         ,TOT,TSFC
242       ! temporary arrays for old chemistry values
243       REAL,DIMENSION(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,num_chem_ct) :: chem_old
244       INTEGER,DIMENSION(num_chem_ct) :: chem_ct_indices
246 ! Variables for adaptive time steps...
247       TYPE(WRFU_TimeInterval) :: tmpTimeInterval
248       REAL(KIND=8) :: curr_secs
249       REAL(KIND=8) :: real_time_r8                                 !ext. function in adapt_timestep_em.F
250       LOGICAL      :: adapt_step_flag, do_chemstep, do_photstep
252       LOGICAL      :: chm_is_mozart
254       REAL :: DAYI,DPL,FICE,FRAIN,HOUR,PLYR          &
255      &       ,QI,QR,QW,RADT,TIMES,WC,TDUM,WMSK,RWMSK
258       INTEGER                         :: ij 
259       INTEGER                         :: im , num_3d_m , ic , num_3d_c, num_3d_s
260       INTEGER                         :: ijds, ijde
261       INTEGER                         :: astat
262       INTEGER                         :: num_irr_diag
263       INTEGER                         :: ksubt
265       REAL :: chem_minval, dtstepc
267       REAL, DIMENSION(:,:,:,:), ALLOCATABLE :: irr_rates
269       REAL, DIMENSION(grid%sm31:grid%em31, grid%sm32:grid%em32, grid%sm33:grid%em33, numgas_mam) :: gas_aqfrac ! fraction of gas that is in cloud water
270       
271       REAL, DIMENSION( grid%sm31:grid%em31, grid%sm32:grid%em32, grid%sm33:grid%em33, ntot_amode_cam_mam ) :: &
272          wetdens_ap  ! modal diameters and wet densities for cam_mam packages
274       REAL, DIMENSION( grid%sm31:grid%em31, grid%sm32:grid%em32, grid%sm33:grid%em33 ) :: &
275          del_h2so4_gasprod  ! change to h2so4 during gas-phase chemistry, cam_mam packages
277       !Balwinder.Singh@pnnl.gov:dvmrdt_sv13d,dvmrcwdt_sv13d are the tendencies which are passsed on from the CAM-MAM cloud chemistry
278       !to gasaerexch subroutine in cam_mam_aerchem_driver
279       REAL, DIMENSION( grid%sm31:grid%em31, grid%sm32:grid%em32, grid%sm33:grid%em33,gas_pcnst_pos) :: dvmrdt_sv13d,dvmrcwdt_sv13d 
281       LOGICAL :: haveaer
282       CHARACTER (LEN=1000) :: msg 
283       CHARACTER (LEN=256) :: current_date_char 
284       integer :: current_month
285 ! ..
286 ! .. Intrinsic Functions ..
287       INTRINSIC max, min
289       REAL, ALLOCATABLE, DIMENSION(:,:,:) :: si_zsigf, si_zsig
290 !<DESCRIPTION>
291 !<pre>
292 ! chem_driver is the main driver for handling chemistry related tasks
293 ! for a particular timestep. chem_driver is a mediation-layer routine ->
294 ! DM and SM calls are made where needed for parallel processing.
296 ! The main sections of chem_driver are:
298 ! (1) Initialization of meteorology variables as needed for chemistry
300 ! (2) Calls to drivers for the various chemistry tasks:
301 !        emissions_driver
302 !        photolysis_driver
303 !        dry_dep_driver
304 !        uoc_dustwd_driver
305 !        grelldrvct (convective tracer transport)
306 !        mechanism_driver (gases)
307 !        cloud_chem_driver
308 !        aerosols_driver
309 !        wetscav_driver
310 !        sum_pm_driver
312 ! Handling of tile indices in chem_driver is as close as possible to
313 ! what is done in solve_em. For subroutines called from chem_driver,
314 ! the its:ite, jts:jte, and kts:kte variables represent the extent of
315 ! the domain that each processor should loop over. For example, a do
316 ! loop in the vertical for the chem array should go from kts to kte.
317 ! For the EM core, kte=kde-1. For the NMM core, kte=kde-2.
319 ! Note that the tile indices for the chemistry initialization differ
320 ! from the integration loop indices in that the initializataion routines
321 ! use kte=kde. Go figure, this is how the met. folks set things up.
323 !</pre>
324 !</DESCRIPTION>
326 ! ..
328 ! Number of levels to exclude from the chem calculations counting from
329 ! the model top.
330 !  ksubt=0
333 ! Setup the adaptive timestep for the chem routines. Most of this follows
334 ! what is in solve_em, except for the call to adjust time_step.
336   !The necesssary variables exist for the EM core and using the adaptive
337   !techniques will work even with a constant time step. In fact, they
338   !prevent issues with restarts and changed time steps. So, we will
339   !always use them with the EM core.
340   adapt_step_flag = .TRUE.
341   ktau = grid%itimestep
342   tmpTimeInterval = domain_get_time_since_sim_start(grid)
343   curr_secs = real_time_r8(tmpTimeInterval)
344   ijulian=ifix(grid%julian)
346   do_photstep = .false.
347   IF ( ktau==1 ) then
348      do_photstep = .true.
349   ELSE IF ( adapt_step_flag ) THEN
350      IF ( (grid%photdt<=0) .or. &
351           ( curr_secs+real(grid%dt,8)+0.01 >= &
352           ( INT( curr_secs/real(grid%photdt*60.,8)+1,8 )*real(grid%photdt*60.,8) ) ) &
353           ) then
354           !NOTE: The 0.01 added to the LHS of the conditional is to compensate
355           !      for occasional round off errors that prevented the >= from
356           !      ever testing as true. This adjustment has been added to the
357           !      other checks within the chem directory as well. wig, 30-Sep-2008
358         do_photstep = .true.
359      ENDIF
360   ELSE IF ( (MOD(ktau,grid%stepphot)==0) .or. (grid%stepphot==1) ) THEN
361      do_photstep = .true.
362   ENDIF
364   if( ktau==1 ) then
365      dtstepc = grid%dt
366   else
367      tmpTimeInterval = domain_get_current_time(grid) - last_chem_time(grid%id)
368      dtstepc = real(real_time_r8(tmpTimeInterval),4)
369   end if
370      
371   ! initializing diagnostics and macros
373   if( ktau==1 ) then
374      grid%conv_ct(:,:,:,:)   = 0.
375      grid%chem_ct(:,:,:,:)   = 0.
376      grid%vmix_ct(:,:,:,:)   = 0.  
377   endif
378   if(config_flags%chemdiag == USECHEMDIAG)then
379   ! modify tendency list here
380       chem_ct_indices(p_chem_co)   = p_co
381       chem_ct_indices(p_chem_o3)   = p_o3
382       chem_ct_indices(p_chem_no)   = p_no
383       chem_ct_indices(p_chem_no2)  = p_no2
384       chem_ct_indices(p_chem_hno3) = p_hno3
385       chem_ct_indices(p_chem_iso)  = p_iso
386       chem_ct_indices(p_chem_ho)   = p_ho
387       chem_ct_indices(p_chem_ho2)  = p_ho2
388   endif
389   
390   do_chemstep = .false.
391   IF ( ktau==1 ) then
392      do_chemstep = .true.
393      grid%ktauc = 1
394   ELSE IF ( adapt_step_flag ) THEN
395      IF ( (grid%chemdt<=0) .or. &
396           ( curr_secs+real(grid%dt,8)+0.01 >= &
397           ( INT( curr_secs/real(grid%chemdt*60.,8)+1,8 )*real(grid%chemdt*60.,8) ) ) &
398           ) then
399         do_chemstep = .true.
400         grid%ktauc = grid%ktauc+1
401         last_chem_time(grid%id) = domain_get_current_time( grid )
402         call WRFU_TimeGet( last_chem_time(grid%id),         &
403                            YY = grid%last_chem_time_year,   &
404                            MM = grid%last_chem_time_month,  &
405                            DD = grid%last_chem_time_day,    &
406                            H  = grid%last_chem_time_hour,   &
407                            M  = grid%last_chem_time_minute, &
408                            S  = grid%last_chem_time_second  )
409      ENDIF
410   ELSE IF ( (MOD(ktau,grid%stepchem)==0) .or. (grid%stepchem==1) ) THEN
411      do_chemstep = .true.
412      grid%ktauc=max(ktau/grid%stepchem,1)
413   ENDIF
415 ! ..
417 ! Some stuff for convective transport...
419 !  convtrans_avglen_m = 30.  !Averaging time for convective transport in min.
420 ! stepave=convtrans_avglen_m*60./grid%dt
422   CALL get_ijk_from_grid (  grid ,                   &
423                             ids, ide, jds, jde, kds, kde,    &
424                             ims, ime, jms, jme, kms, kme,    &
425                             ips, ipe, jps, jpe, kps, kpe    )
427 ! following two lines needed for MEGAN
428   CALL domain_clock_get( grid, current_timestr=current_date_char )
429   read(current_date_char(6:7),FMT='(I2)') current_month
431 !initialize
433   seasin(:,:,:)=0.
434   dustin(:,:,:)=0.
436   if(config_flags%cu_diag == 0 ) grid%raincv_b(:,:) = grid%raincv(:,:)
438   num_3d_m        = num_moist
439   num_3d_c        = num_chem
440   num_3d_s        = num_scalar
441   numgas          = get_last_gas(config_flags%chem_opt)
444    !  Compute these starting and stopping locations for each tile and number of tiles.
445   CALL set_tiles ( grid , ids , ide , jds , jde , ips , ipe , jps , jpe )
446   k_start         = kps
447   k_end           = kpe
449   ijds = min(ids, jds)
450   ijde = max(ide, jde)
451    chem_minval = epsilc !chem_minval can be case dependant and set below...
452    chem_select: SELECT CASE(config_flags%chem_opt)
453      CASE (RADM2)
454        CALL wrf_debug(15,'calling radm2 from chem_driver')
455        haveaer = .false.
456      CASE (RADM2_KPP)
457        CALL wrf_debug(15,'calling radm2_kpp from chem_driver')
458        haveaer = .false.
459      CASE (CRIMECH_KPP)
460        CALL wrf_debug(15,'calling crimech_kpp from chem_driver')
461        haveaer = .false.
462      CASE (RADM2SORG)
463        CALL wrf_debug(15,'calling radm2sorg aerosols driver from chem_driver')
464        haveaer = .true.
465      CASE (RADM2SORG_KPP)
466        CALL wrf_debug(15,'calling radm2sorg aerosols driver from chem_driver')
467        haveaer = .false.
468      CASE (RADM2SORG_AQ)
469        CALL wrf_debug(15,'calling radm2sorg_aq aerosols driver from chem_driver')
470        haveaer = .true.
471      CASE (RACMSORG_AQ)
472        CALL wrf_debug(15,'calling racmsorg_aq aerosols driver from chem_driver')
473        haveaer = .true.
474      CASE (RADM2SORG_AQCHEM)
475        CALL wrf_debug(15,'calling radm2sorg_aqchem aerosols driver from chem_driver')
476        haveaer = .true.
477      CASE (RACMSORG_AQCHEM_KPP)
478        CALL wrf_debug(15,'calling racmsorg_aqchem_kpp aerosols driver from chem_driver')
479        haveaer = .true.
480      CASE (RACM_ESRLSORG_AQCHEM_KPP)
481        CALL wrf_debug(15,'calling racm_esrlsorg_aqchem_kpp aerosols driver from chem_driver')
482        haveaer = .true.
483      CASE (RACM_KPP)
484        CALL wrf_debug(15,'calling racm_kpp from chem_driver')
485      CASE (RACMPM_KPP)
486        CALL wrf_debug(15,'calling racmpm_kpp from chem_driver')
487        haveaer = .false.
488      CASE (RACM_MIM_KPP)
489        CALL wrf_debug(15,'calling racm_mim_kpp from chem_driver')
490        haveaer = .false.
491      CASE (RACM_ESRLSORG_KPP)
492        CALL wrf_debug(15,'calling racmsorgesrl_kpp aerosols driver from chem_driver')
493        haveaer = .false.
494      CASE (RACMSORG_KPP)
495        CALL wrf_debug(15,'calling racmsorg_kpp aerosols driver from chem_driver')
496        haveaer = .false.
497      CASE (RACM_SOA_VBS_KPP)
498        CALL wrf_debug(15,'calling racm_soa_vbs_kpp aerosols driver from chem_driver')
499        haveaer = .false.
500 !!! TUCCELLA
501      CASE (RACM_SOA_VBS_AQCHEM_KPP)
502        CALL wrf_debug(15,'calling racm_soa_vbs_aqchem_kpp aerosols driver from chem_driver')
503        haveaer = .false.
504      CASE (RACM_SOA_VBS_HET_KPP)
505        CALL wrf_debug(15,'calling racm_soa_vbs_het_kpp aerosols driver from chem_driver')
506        haveaer = .false.
507      CASE (GOCART_SIMPLE)
508        CALL wrf_debug(15,'calling only gocart aerosols driver from chem_driver')
509        haveaer = .false.
510      CASE (GOCARTRACM_KPP)
511        CALL wrf_debug(15,'calling gocart and racm driver from chem_driver')
512        haveaer = .false.
513      CASE (GOCARTRADM2)
514        CALL wrf_debug(15,'calling gocart and radm driver from chem_driver')
515        haveaer = .false.
516      CASE (SAPRC99_KPP)
517        CALL wrf_debug(15,'calling saprc99_kpp from chem_driver')
518        haveaer = .false.
519      CASE (SAPRC99_MOSAIC_4BIN_VBS2_KPP)
520        CALL wrf_debug(15,'calling saprc99_mosaic_4bin_vbs2_kpp from chem_driver')
521        haveaer = .false.
522      CASE (MOZART_MOSAIC_4BIN_KPP) 
523        CALL wrf_debug(15,'calling mozart_mosaic_4bin_kpp from chem_driver')
524 #if( WRF_CHEM == 1 && WRF_KPP == 1 )
525        IF( grid%irr_opt == 1 ) then
526          ALLOCATE( irr_rates(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,num_irr_diag_mozart_mosaic_4bin),stat=astat )
527          IF( astat /= 0 ) THEN
528            write(msg,'(''chem_driver: Failed to allocate irr_rates; error = '',i8)') astat
529            CALL wrf_error_fatal( trim(msg) )
530          ENDIF
531          num_irr_diag = num_irr_diag_mozart_mosaic_4bin
532        ENDIF
533 #endif
534        haveaer = .true.
535      CASE (MOZART_MOSAIC_4BIN_AQ_KPP)
536        CALL wrf_debug(15,'calling mozart_mosaic_4bin_aq_kpp from chem_driver')
537 #if( WRF_CHEM == 1 && WRF_KPP == 1 )
538        IF( grid%irr_opt == 1 ) then
539          ALLOCATE( irr_rates(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,num_irr_diag_mozart_mosaic_4bin_aq),stat=astat )
540          IF( astat /= 0 ) THEN
541            write(msg,'(''chem_driver: Failed to allocate irr_rates; error = '',i8)') astat
542            CALL wrf_error_fatal( trim(msg) )
543          ENDIF
544          num_irr_diag = num_irr_diag_mozart_mosaic_4bin_aq
545        ENDIF
546 #endif
547        haveaer = .true.
548        !BSINGH -  Added case statement for SAPRC 8 bin
549      CASE (SAPRC99_MOSAIC_8BIN_VBS2_AQ_KPP)
550        CALL wrf_debug(15,'calling saprc99_mosaic_8bin_vbs2_aq_kpp from chem_driver')
551        haveaer = .false.
552      CASE (SAPRC99_MOSAIC_8BIN_VBS2_KPP) !BSINGH(04/04/2014)
553        CALL wrf_debug(15,'calling saprc99_mosaic_8bin_vbs2_kpp from chem_driver')
554        haveaer = .false.
555        !BSINGH - ENDS
556      CASE (CBMZSORG)
557        CALL wrf_debug(15,'calling cbmzsorg aerosols from chem_driver')
558        haveaer = .true.
559      CASE (CBMZSORG_AQ)
560        CALL wrf_debug(15,'calling cbmzsorg_aq aerosols from chem_driver')
561        haveaer = .true.
562      CASE (CBMZ)
563        CALL wrf_debug(15,'calling cbmz from chem_driver')
564        haveaer = .false.
565      CASE (CBMZ_BB)
566        CALL wrf_debug(15,'calling cbmz_bb from chem_driver')
567        haveaer = .false.
568      CASE (CBMZ_BB_KPP)
569        CALL wrf_debug(15,'calling cbmz_bb_kpp from chem_driver')
570        haveaer = .false.
571      CASE (CBMZ_MOSAIC_KPP)
572        CALL wrf_debug(15,'calling cbmz_mosaic_kpp from chem_driver')
573        haveaer = .false.
574      CASE (CBMZ_MOSAIC_4BIN)
575        CALL wrf_debug(15,'calling cbmz_mosaic_4bin aerosols driver from chem_driver')
576        haveaer = .true.
577      CASE (CBMZ_MOSAIC_8BIN)
578        CALL wrf_debug(15,'calling cbmz_mosaic_8bin aerosols driver from chem_driver')
579        haveaer = .true.
580      CASE (CBMZ_MOSAIC_4BIN_AQ)
581        CALL wrf_debug(15,'calling cbmz_mosaic_4bin_aq aerosols driver from chem_driver')
582        haveaer = .true.
583      CASE (CBMZ_MOSAIC_8BIN_AQ)
584        CALL wrf_debug(15,'calling cbmz_mosaic_8bin_aq aerosols driver from chem_driver')
585        haveaer = .true.
586      CASE (CBMZ_MOSAIC_DMS_4BIN)
587        CALL wrf_debug(15,'calling cbmz_mosaic_dms_4bin aerosols driver from chem_driver')
588        haveaer = .true.
589      CASE (CBMZ_MOSAIC_DMS_8BIN)
590        CALL wrf_debug(15,'calling cbmz_mosaic_dms_8bin aerosols driver from chem_driver')
591        haveaer = .true.
592      CASE (CBMZ_MOSAIC_DMS_4BIN_AQ)
593        CALL wrf_debug(15,'calling cbmz_mosaic_dms_4bin_aq aerosols driver from chem_driver')
594        haveaer = .true.
595      CASE (CBMZ_MOSAIC_DMS_8BIN_AQ)
596        CALL wrf_debug(15,'calling cbmz_mosaic_dms_8bin_aq aerosols driver from chem_driver')
597        haveaer = .true.
598      CASE (CRI_MOSAIC_8BIN_AQ_KPP)
599        CALL wrf_debug(15,'calling cri_mosaic_8bin_aq_kpp aerosols driver from chem_driver')
600        haveaer = .true.
601      CASE (CRI_MOSAIC_4BIN_AQ_KPP)
602        CALL wrf_debug(15,'calling cri_mosaic_4bin_aq_kpp aerosols driver from chem_driver')
603        haveaer = .true.
604      CASE (MOZART_KPP)
605        CALL wrf_debug(15,'calling mozart driver from chem_driver')
606      CASE (MOZCART_KPP)
607        CALL wrf_debug(15,'calling mozcart driver from chem_driver')
608 #if( WRF_CHEM == 1 && WRF_KPP == 1 )
609        IF( grid%irr_opt == 1 ) then
610          ALLOCATE( irr_rates(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,num_irr_diag_mozcart),stat=astat )
611          IF( astat /= 0 ) THEN
612            write(msg,'(''chem_driver: Failed to allocate irr_rates; error = '',i8)') astat
613            CALL wrf_error_fatal( trim(msg) )
614          ENDIF
615          num_irr_diag = num_irr_diag_mozcart
616        ENDIF
617 #endif
618      CASE (T1_MOZCART_KPP)
619        CALL wrf_debug(15,'calling t1_mozcart driver from chem_driver')
620 #if( WRF_CHEM == 1 && WRF_KPP == 1 )
621        IF( grid%irr_opt == 1 ) then
622          ALLOCATE( irr_rates(grid%sm31:grid%em31,grid%sm32:grid%em32,grid%sm33:grid%em33,num_irr_diag_t1_mozcart),stat=astat )
623          IF( astat /= 0 ) THEN
624            write(msg,'(''chem_driver: Failed to allocate irr_rates; error = '',i8)') astat
625            CALL wrf_error_fatal( trim(msg) )
626          ENDIF
627          num_irr_diag = num_irr_diag_t1_mozcart
628        ENDIF
629 #endif
630      CASE (CB05_SORG_AQ_KPP)
631        CALL wrf_debug(15,'calling cb05_sorg_aq_kpp from chem_driver')
632        haveaer = .true.
633      CASE (CB05_SORG_VBS_AQ_KPP)
634        CALL wrf_debug(15,'calling cb05_sorg_vbs_aq_kpp from chem_driver')
635        haveaer = .true.
636     CASE (CHEM_TRACER,CHEM_TRACE2)
637        CALL wrf_debug(15,'tracer mode: only doing emissions and dry dep in chem_driver')
638     CASE (CHEM_VOLC)
639        CALL wrf_debug(15,'Full Volcanic Ash mode: doing emissions (SO2 + ASH), settling, and subgrid transport in chem_driver')
640     CASE (CHEM_VOLC_4BIN)
641        CALL wrf_debug(15,'4bin Volcanic Ash mode: doing emissions (ASH), settling, and subgrid transport in chem_driver')
642     CASE (CHEM_VASH)
643        CALL wrf_debug(15,'Volcanic Ash mode: only doing emissions, settling, and subgrid transport in chem_driver')
644     CASE (DUST)
645        CALL wrf_debug(15,'Dust only mode: only doing emissions, settling, and subgrid transport chem_driver')
646     CASE (CO2_TRACER,GHG_TRACER)
647       CALL wrf_debug(15,'Greenhouse gas mode: fluxes and transport of GHG')
648     CASE DEFAULT
649        if(config_flags%tracer_opt > 0 )then
650        CALL wrf_debug(15,'only doing tracer transport in chem_driver')
651        else
652        CALL wrf_debug(15,'calling chem_opt=? from chem_driver')
653        endif
654    END SELECT chem_select                              
655    tracer_select: SELECT CASE(config_flags%tracer_opt)
656     CASE (TRACER_SMOKE)
657        CALL wrf_debug(15,'tracer mode: 1 tracer for fires')
658     CASE (TRACER_TEST1)
659        CALL wrf_debug(15,'tracer mode: 8 tracers')
660     CASE (TRACER_TEST2)
661        CALL wrf_debug(15,'tracer mode: 8 tracers')
662     CASE (TRACER_TEST3)
663        CALL wrf_debug(15,'tracer mode: 10 tracers')
664      CASE DEFAULT
665        CALL wrf_debug(15,'calling chem_opt=? from chem_driver')
666    END SELECT tracer_select
668 ! initialize cam_mam local arrays
669    if ((config_flags%chem_opt == CBMZ_CAM_MAM3_NOAQ) .or. &
670        (config_flags%chem_opt == CBMZ_CAM_MAM3_AQ  ) .or. &
671        (config_flags%chem_opt == CBMZ_CAM_MAM7_NOAQ) .or. &
672        (config_flags%chem_opt == CBMZ_CAM_MAM7_AQ  )) then
673       grid%dgnum4d(:,:,:,:) = 0.0 !Balwinder.Singh@pnnl: dgnum is now defined in Registry as dgnum4d
674       grid%dgnumwet4d(:,:,:,:) = 0.0 !Balwinder.Singh@pnnl: dgnumwet is now defined in Registry as dgnumwet4d
675       wetdens_ap(:,:,:,:) = 0.0
676       !BSINGH: 01/31/2013: Sanity check for numgas_mam
677       if(numgas_mam < numgas) then
678          write(msg,*)'CHEM_DRIVER - NUMGAS_MAM is should be equal to numgas (check chemics_init.F), numgas_mam=',numgas_mam,' and numgas=',numgas
679          call wrf_error_fatal( msg )
680       endif
681       !BSINGH:02/01/2013 - Sanity check for cam_mam_aerols, it should be true (updated in module_cam_mam_init.F)
682       if(.NOT.cam_mam_aerosols) then
683          write(msg,*)'CHEM_DRIVER - cam_mam_aerosol should be TRUE (check module_physics_init.F), module_cam_mam_aerosol=',cam_mam_aerosols
684          call wrf_error_fatal(msg)
685       endif
686    end if
691       do nv=1,num_chem
692          do j=jps,jpe
693             do k=kps,kpe
694                do i=ips,ipe
695                   chem(i,k,j,nv)=max(chem(i,k,j,nv),chem_minval)
696                enddo
697             enddo
698          enddo
699       enddo
700       select case (config_flags%chem_opt)
701 !!! TUCCELLA
702       case (RADM2SORG, RADM2SORG_AQ, RADM2SORG_AQCHEM, RADM2SORG_KPP, &
703             RACM_ESRLSORG_KPP,RACMSORG_AQ,RACMSORG_KPP, RACMSORG_AQCHEM_KPP, RACM_ESRLSORG_AQCHEM_KPP, &
704             RACM_SOA_VBS_KPP,RACM_SOA_VBS_AQCHEM_KPP,RACM_SOA_VBS_HET_KPP)
705          do j=jps,jpe
706             do k=kps,kpe
707                do i=ips,ipe
708                   if(chem(i,k,j,p_nu0).lt.1.e07) then
709                      chem(i,k,j,p_nu0)=1.e7
710                   endif
711                enddo
712             enddo
713          enddo
715 ! Special treatment of CH4 in SAPRC99
716       case (SAPRC99_KPP,SAPRC99_MOSAIC_4BIN_VBS2_KPP, &
717            SAPRC99_MOSAIC_8BIN_VBS2_AQ_KPP,SAPRC99_MOSAIC_8BIN_VBS2_KPP) !BSINGH -  Added for SAPRC 8 bin and non-aq on (04/04/2014)
718          CALL wrf_debug ( 15 , ' fixing ch4 conc using co conc' )
719          do j=jps,jpe
720             do k=kps,kpe
721                do i=ips,ipe
722                   chem(i,k,j,p_ch4)=1.74
723                enddo
724             enddo
725          enddo
726       end select
728       vdrog3=0.
729       do j=jps,min(jde-1,jpe)
730          do k=kps,kpe
731             do i=ips,min(ide-1,ipe)
732               vvel(i,k,j)=grid%w_2(i,k,j)
733               zmid(i,k,j)=grid%z(i,k,j)
734             enddo
735          enddo
736       enddo
737       do j=jps,min(jde-1,jpe)
738          do k=kps,min(kde-1,kpe)
739             do i=ips,min(ide-1,ipe)
740               rri(i,k,j)=grid%alt(i,k,j)
741             enddo
742          enddo
743       enddo
744       do j=jps,min(jde-1,jpe)
745          do i=ips,min(ide-1,ipe)
746             pbl_h(i,j)=grid%pblh(i,j)
747          enddo
748       enddo
750       chm_is_mozart = config_flags%chem_opt == MOZART_KPP .or. &                               
751                       config_flags%chem_opt == MOZCART_KPP .or. &
752                       config_flags%chem_opt == T1_MOZCART_KPP .or. &
753                       config_flags%chem_opt == MOZART_MOSAIC_4BIN_KPP .or. &
754                       config_flags%chem_opt == MOZART_MOSAIC_4BIN_AQ_KPP
755 !------------------------------------------------------------------------
756 ! setup ftuv column density timing
757 !------------------------------------------------------------------------
758       if( chm_is_mozart .and. config_flags%phot_opt == FTUV ) then
759         CALL ftuv_timestep_init( grid%id, grid%julday )
760       endif
762 !------------------------------------------------------------------------
763 ! Main chemistry tile loop
764 !------------------------------------------------------------------------
765      
766 !$OMP PARALLEL DO   &
767 !$OMP PRIVATE ( ij, its, ite, jts, jte )
768    chem_tile_loop_1: DO ij = 1 , grid%num_tiles
769        its = grid%i_start(ij) 
770        ite = min(grid%i_end(ij),ide-1)
771        jts = grid%j_start(ij)
772        jte = min(grid%j_end(ij),jde-1) 
774        kts=k_start
775        kte=min(k_end,kde-1)
780          CALL wrf_debug ( 15 , ' call chem_prep' )
781          CALL chem_prep ( config_flags,                                               &
782                          grid%u_2, grid%v_2, grid%p, grid%pb,             &
783                          grid%alt,grid%ph_2, grid%phb, grid%t_2,          &
784                          moist, num_3d_m, rho,                                        &
785                          p_phy,  u_phy, v_phy,                                        &
786                          p8w, t_phy, t8w, grid%z, z_at_w,                          &
787                          dz8w, rh, grid%fnm, grid%fnp,                              &
788                          ids, ide, jds, jde, kds, kde,                                &
789                          ims, ime, jms, jme, kms, kme,                                &
790                          its,ite,jts,jte,                                             &
791                          k_start, k_end                                               )
793 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
794        if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and.                          &
795            (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and.                          &
796            (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
797           call wrf_debug(15,"calling chem_dbg at top of chem_driver")
798           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
799                dz8w,t_phy,p_phy,rho,chem,emis_ant,                       &
800                ids,ide, jds,jde, kds,kde,                                             &
801                ims,ime, jms,jme, kms,kme,                                             &
802                its,ite, jts,jte,kts,kte,                                              &
803                config_flags%kemit                                                     )
804        end if
805 #endif
808 !--- emissions
809 !  MEGAN2_CLM uses driver below even if other flags are zero.  BJG 5/29/2019
811     if(config_flags%emiss_inpt_opt > 0 .or. config_flags%dust_opt > 0  &
812        .or. config_flags%tracer_opt > 0 .or. config_flags%bio_emiss_opt == MEGAN2_CLM ) then
813       call wrf_debug(15,'calling emissions driver')
815       call emissions_driver(grid%id,ktau,grid%dt,grid%DX,                                  &
816               adapt_step_flag, curr_secs,                                                  &
817               grid%plumerisefire_frq,grid%stepfirepl,                                      &
818               grid%bioemdt,grid%stepbioe,                                                  &
819               config_flags,                                                                &
820               grid%gmt,ijulian,rri,t_phy,moist,p8w,t8w,u_phy,v_phy,vvel,               &
821               grid%e_bio,p_phy,chem,rho,dz8w,grid%ne_area,emis_ant,emis_vol,grid%tsk,      &
822               grid%erod,grid%erod_dri,grid%lai_vegmask,                                    &
823               g,emis_seas,emis_dust,tracer,                                                &
824               emis_seas2,      &
825               ebu, ebu_in,grid%mean_fct_agtf,grid%mean_fct_agef,grid%mean_fct_agsv,       &
826               grid%mean_fct_aggr,grid%firesize_agtf, &
827               grid%firesize_agef,grid%firesize_agsv,grid%firesize_aggr,                    &
828               grid%u10,grid%v10,grid%ivgtyp,grid%isltyp,grid%gsw,grid%vegfra,grid%rmol,    &
829               grid%ust,grid%znt,grid%dms_0,grid%erup_beg,grid%erup_end,                    &
830               grid%xland,grid%xlat,grid%xlong,                                             &
831               z_at_w,zmid,grid%smois,dustin,seasin,                                        &
832               grid%sebio_iso,grid%sebio_oli,grid%sebio_api,grid%sebio_lim,                 &
833               grid%sebio_xyl,grid%sebio_hc3,grid%sebio_ete,grid%sebio_olt,                 &
834               grid%sebio_ket,grid%sebio_ald,grid%sebio_hcho,grid%sebio_eth,                &
835               grid%sebio_ora2,grid%sebio_co,grid%sebio_nr,                                 &
836               grid%sebio_sesq,grid%sebio_mbo,                                              &
837               grid%noag_grow,grid%noag_nongrow,grid%nononag,grid%slai,                     &
838               grid%ebio_iso,grid%ebio_oli,grid%ebio_api,grid%ebio_lim,grid%ebio_xyl,       &
839               grid%ebio_hc3,grid%ebio_ete,grid%ebio_olt,grid%ebio_ket,grid%ebio_ald,       &
840               grid%ebio_hcho,grid%ebio_eth,grid%ebio_ora2,grid%ebio_co,grid%ebio_nr,       &
841               grid%ebio_no,grid%ebio_sesq,grid%ebio_mbo,grid%ebio_bpi,grid%ebio_myrc,      &
842               grid%ebio_c10h16,grid%ebio_tol,grid%ebio_bigalk,                             &
843               grid%ebio_ch3oh,grid%ebio_acet,grid%ebio_nh3,grid%ebio_no2,                  &
844               grid%ebio_c2h5oh,grid%ebio_ch3cooh,grid%ebio_mek,grid%ebio_bigene,           &
845               grid%ebio_c2h6,grid%ebio_c2h4,grid%ebio_c3h6,grid%ebio_c3h8,grid%ebio_so2,   &
846               grid%ebio_dms,grid%ebio_hcn,                                                 &
847               grid%ebio_alk3, grid%ebio_alk4, grid%ebio_alk5, grid%ebio_ole1, grid%ebio_ole2,            &    
848               grid%ebio_aro1, grid%ebio_aro2, grid%ebio_ccho, grid%ebio_meoh,             &    
849               grid%ebio_ethene, grid%ebio_hcooh, grid%ebio_terp, grid%ebio_bald,          &    
850               grid%ebio_cco_oh, grid%ebio_rco_oh,                                         &    
851               grid%clayfrac,grid%sandfrac,grid%dust_alpha,grid%dust_gamma,grid%dust_smtune, grid%dust_ustune, &
852               grid%clayfrac_nga,grid%sandfrac_nga,                                        &
853               grid%snowh,grid%zs,grid%afwa_dustloft,                                       &
854               grid%tot_dust,grid%tot_edust,grid%vis_dust,                                  &
855               grid%soilctop, grid%ust_t, grid%rough_cor, grid%smois_cor,                  &              
856               grid%ebio_c5h8,grid%ebio_apinene,grid%ebio_bpinene,grid%ebio_toluene,        &
857               grid%ebio_ch3cho,grid%ebio_ch3co2h,grid%ebio_tbut2ene,                      &
858               grid%ebio_c2h5cho,grid%ebio_nc4h10,                                                             &
859               grid%T2,grid%swdown,                                                         &
860               grid%nmegan,grid%EFmegan,                                                    &
861               grid%msebio_isop,                                                            &
862               grid%mlai,                                                                   &
863               grid%pftp_bt, grid%pftp_nt, grid%pftp_sb, grid%pftp_hb,                      &
864               grid%mtsa,                                                                   &
865               grid%mswdown,                                                                &
866               grid%mebio_isop,grid%mebio_apin,grid%mebio_bpin, grid%mebio_bcar,            &
867               grid%mebio_acet,grid%mebio_mbo,grid%mebio_no,                                &
868               current_month,                                                               &
869          ! stuff for LNOx emissions
870              grid%ht, grid%refl_10cm, grid%ic_flashrate, grid%cg_flashrate,                 &
871          ! stuff for aircraft emissions
872               emis_aircraft,                                                               &
873          ! stuff for the GHG fluxes
874               vprm_in,grid%rad_vprm,grid%lambda_vprm,                                      &
875               grid%alpha_vprm,grid%resp_vprm,grid%xtime,                                   &
876               grid%TSLB, wet_in,grid%RAINC,grid%RAINNC,                                    &
877               grid%potevp,grid%SFCEVP,grid%LU_INDEX,                                       &
878               grid%biomt_par,grid%emit_par,grid%ebio_co2oce,                               &
879               eghg_bio,                                                                    &
880               grid%seas_flux,                                                              &
881               ids,ide, jds,jde, kds,kde,                                                   &
882               ims,ime, jms,jme, kms,kme,                                                   &
883               its,ite,jts,jte,kts,kte)
884               if( chm_is_mozart ) then
885                  call mozcart_lbc_set( chem, num_chem, grid%id, &
886                                        ims, ime, jms, jme, kms, kme,    &                  
887                                        its, ite, jts, jte, kts )                           
888               end if
890 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
891        if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and.                          &
892            (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and.                          &
893            (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
894           call wrf_debug(15,'calling chem_dbg after emissions_driver')
895           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
896                dz8w,t_phy,p_phy,rho,chem,emis_ant,                                    &
897                ids,ide, jds,jde, kds,kde,                                             &
898                ims,ime, jms,jme, kms,kme,                                             &
899                its,ite, jts,jte, kts, kte,                                            &
900                config_flags%kemit                                                     )
901        end if
902 #endif 
903     endif
906 ! calculate aerosol optical properities
908       if( do_photstep .and. &
909            config_flags%chem_opt /= CHEM_TRACER .and. &
910            config_flags%chem_opt /= CHEM_VASH .and. &
911            config_flags%chem_opt /= CHEM_VOLC .and. &
912            config_flags%chem_opt /= CHEM_VOLC_4BIN .and. &
913            config_flags%chem_opt /= DUST .and. &
914            config_flags%chem_opt /= CHEM_TRACE2 .and. &
915            config_flags%chem_opt /= CO2_TRACER  .and. &
916            config_flags%chem_opt /= GHG_TRACER ) then
917          call wrf_debug(15,'calling optical driver')
918          call optical_driver (grid%id,curr_secs,grid%dt,config_flags,haveaer,         &
919               chem,dz8w,rri,rh,                                                       &
920               grid%h2oai,grid%h2oaj,                                                  &
921               grid%tauaer1,grid%tauaer2,grid%tauaer3,grid%tauaer4,                    &
922              !czhao
923               grid%extaer1,grid%extaer2,grid%extaer3,grid%extaer4,                    &
924               grid%gaer1,grid%gaer2,grid%gaer3,grid%gaer4,                            &
925               grid%waer1,grid%waer2,grid%waer3,grid%waer4,                            &
926               grid%bscoef1,grid%bscoef2,grid%bscoef3,grid%bscoef4,                    &
927               grid%l2aer,grid%l3aer,grid%l4aer,grid%l5aer,grid%l6aer,grid%l7aer,      &
928               grid%totoa_a01,grid%totoa_a02,grid%totoa_a03,grid%totoa_a04,            &
929               grid%totoa_a05,grid%totoa_a06,grid%totoa_a07,grid%totoa_a08,            &
930               grid%extaerlw1,grid%extaerlw2,grid%extaerlw3,grid%extaerlw4,grid%extaerlw5, &
931               grid%extaerlw6,grid%extaerlw7,grid%extaerlw8,grid%extaerlw9,grid%extaerlw10, &
932               grid%extaerlw11,grid%extaerlw12,grid%extaerlw13,grid%extaerlw14,grid%extaerlw15, &
933               grid%extaerlw16,    &
934               grid%tauaerlw1,grid%tauaerlw2,grid%tauaerlw3,grid%tauaerlw4,grid%tauaerlw5,  &
935               grid%tauaerlw6,grid%tauaerlw7,grid%tauaerlw8,grid%tauaerlw9,grid%tauaerlw10,  &
936               grid%tauaerlw11,grid%tauaerlw12,grid%tauaerlw13,grid%tauaerlw14,grid%tauaerlw15,  &
937               grid%tauaerlw16,    &
938               ids,ide, jds,jde, kds,kde,                                              &
939               ims,ime, jms,jme, kms,kme,                                              &
940               its,ite, jts,jte, kts,kte)
941       endif
944 ! calculate photolysis rates
946       if( do_photstep .and. &
947            config_flags%chem_opt /= CHEM_TRACER .and. &
948            config_flags%chem_opt /= CHEM_VASH .and. &
949            config_flags%chem_opt /= CHEM_VOLC .and. &
950            config_flags%chem_opt /= CHEM_VOLC_4BIN .and. &
951            config_flags%chem_opt /= DUST .and. &
952            config_flags%chem_opt /= CHEM_TRACE2 .and. &
953            config_flags%chem_opt /= CO2_TRACER  .and. &
954            config_flags%chem_opt /= GHG_TRACER  ) then
955          call wrf_debug(15,'calling photolysis driver')
956          call photolysis_driver (grid%id,curr_secs,ktau,grid%dt,                      &
957               config_flags,haveaer,                                                   &
958               grid%dt_cld,grid%af_dir,grid%af_dn,grid%af_up,grid%ph_par,grid%ph_erythema,   &
959               grid%gmt,ijulian,t_phy,moist,grid%aerwrf,p8w,t8w,p_phy,                 &
960               chem,rho,dz8w,grid%xlat,grid%xlong,                                     &
961               zmid,z_at_w,                                                            &
962               grid%qc_cu,grid%qi_cu,                                                  &
963               grid%qc_bl,grid%qi_bl,                                                  &
964               grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                     &
965               grid%ph_clno2,                                                          &
966               grid%ph_no3o2,                                                          &
967               grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,       &
968               grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,            &
969               grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,                      &
970               grid%ph_hcochest,grid%ph_ch3o2h,grid%ph_ch3coo2h,grid%ph_ch3ono2,       &
971               grid%ph_hcochob,grid%ph_n2o5,grid%ph_o2,grid%ph_n2o,                    &
972               grid%ph_pan,grid%ph_mpan,grid%ph_acetol,grid%ph_gly,                    &
973               grid%ph_open,grid%ph_mek,grid%ph_etooh,grid%ph_prooh,grid%ph_pooh,      &
974               grid%ph_acetp,grid%ph_xooh,grid%ph_isooh,grid%ph_alkooh,                &
975               grid%ph_mekooh,grid%ph_tolooh,grid%ph_terpooh,grid%ph_mvk,              &
976               grid%ph_glyald,grid%ph_hyac,                                            &
977               grid%ph_cl2,grid%ph_hocl,grid%ph_fmcl,                                  &
978               config_flags%track_tuv_lev,                                             &
979               config_flags%track_rad_num,                                             &
980               config_flags%track_tuv_num,                                             &
981               grid%radfld,grid%adjcoe,grid%phrate,                                    &
982               grid%track_wc,grid%track_zref,                                          &
983               grid%tauaer1,grid%tauaer2,grid%tauaer3,grid%tauaer4,                    &
984               grid%gaer1,grid%gaer2,grid%gaer3,grid%gaer4,                            &
985               grid%waer1,grid%waer2,grid%waer3,grid%waer4,                            &
986               grid%bscoef1,grid%bscoef2,grid%bscoef3,grid%bscoef4,                    &
987               grid%l2aer,grid%l3aer,grid%l4aer,grid%l5aer,grid%l6aer,grid%l7aer,      &
988               grid%pm2_5_dry,grid%pm2_5_water,grid%uvrad,grid%ivgtyp,                 &
989               grid%o3_gfs_du,                                                         &
990               ids,ide, jds,jde, kds,kde,                                              &
991               ims,ime, jms,jme, kms,kme,                                              &
992               its,ite,jts,jte,kts,kte)
994 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
995          if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and. &
996              (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and. &
997              (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
998           call wrf_debug(15,'calling chem_dbg after photolysis_driver')
999           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
1000                dz8w,t_phy,p_phy,rho,chem,emis_ant,                                    &
1001                ids,ide, jds,jde, kds,kde,                                             &
1002                ims,ime, jms,jme, kms,kme,                                             &
1003                its,ite, jts,jte, kts,kte,                                             &
1004                config_flags%kemit,                                                    &
1005                grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                    &
1006                grid%ph_cl2,grid%ph_hocl,grid%ph_clno2,grid%ph_fmcl,                   &
1007                grid%ph_no3o2,                                                         &
1008                grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,      &
1009                grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,           &
1010                grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,                     &
1011                grid%ph_hcochest,grid%ph_ch3o2h,grid%ph_ch3coo2h,                      &
1012                grid%ph_ch3ono2,grid%ph_hcochob,grid%ph_n2o5,grid%ph_o2                &
1013                                                                                       )
1014        end if
1015 #endif
1016       endif
1019 ! do vertical mixing with dry deposition
1021 ! save old concentrations for vertical mixing tendencies
1022    DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1023      chem_old(its:ite,kts:kte,jts:jte,nv) = chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv))
1024    ENDDO
1026 ! UoC dust scheme option
1027    scheme_select: SELECT CASE(config_flags%dust_schme)
1028    CASE (SHAO_2001)
1029      imod = 1
1030    CASE (SHAO_2004)
1031      imod = 2
1032    CASE (SHAO_2011)
1033      imod = 3
1034    CASE DEFAULT
1035      imod = 2
1036    END SELECT scheme_select
1038       if (config_flags%vertmix_onoff>0) then
1039          if (ktau.gt.2) then
1040             call wrf_debug(15,'calling dry_deposition_driver')
1041             call dry_dep_driver(grid%id,curr_secs,ktau,grid%dt,config_flags,          &
1042                  grid%gmt,ijulian,t_phy,moist,scalar,p8w,t8w,vvel,                &
1043                  rri,p_phy,chem,tracer,rho,dz8w,rh,grid%exch_h,grid%hfx,grid%dx,      & 
1044                  grid%cldfra, grid%cldfra_old,grid%raincv_b,seasin,dustin,            &
1045                  grid%ccn1, grid%ccn2, grid%ccn3, grid%ccn4, grid%ccn5, grid%ccn6,    &
1046                  grid%qndropsource,grid%ivgtyp,grid%tsk,grid%gsw,grid%vegfra,pbl_h,   &
1047                  grid%rmol,grid%ust,grid%znt,grid%xlat,grid%xlong,                    &
1048                  zmid,z_at_w,grid%xland,grid%ash_fall,                                &
1049                  grid%h2oaj,grid%h2oai,grid%nu3,grid%ac3,grid%cor3,grid%asulf,        &
1050                  grid%ahno3,grid%anh3,grid%cvaro1,grid%cvaro2,grid%cvalk1,grid%cvole1,&
1051                  grid%cvapi1,grid%cvapi2,grid%cvlim1,grid%cvlim2,grid%dep_vel_o3,     &
1052                  grid%ddlen,grid%ddflx, &
1053                  emis_ant,ebu_in,                                                     &
1054                  config_flags%sf_urban_physics,numgas,current_month,dvel,grid%snowh,  &
1055                  grid%dustdrydep_1,grid%dustdrydep_2,grid%dustdrydep_3,               &
1056                  grid%dustdrydep_4,grid%dustdrydep_5,                                 &      
1057                  grid%depvelocity,                                                    &      
1058                  grid%dustgraset_1,grid%dustgraset_2,grid%dustgraset_3,               &
1059                  grid%dustgraset_4,grid%dustgraset_5,                                 &
1060                  grid%setvel_1,grid%setvel_2,grid%setvel_3,grid%setvel_4,             &
1061                  grid%setvel_5, imod,                                                 &                  
1062                  grid%is_CAMMGMP_used,                                                &
1063                  grid%dep_vel,grid%num_vert_mix,                                      &
1064                  ids,ide, jds,jde, kds,kde,                                           &
1065                  ims,ime, jms,jme, kms,kme,                                           &
1066                  its,ite,jts,jte,kts,kte)
1067                  
1068             if( chm_is_mozart ) then
1069                call mozcart_lbc_set( chem, num_chem, grid%id, &
1070                                      ims, ime, jms, jme, kms, kme,    & 
1071                                      its, ite, jts, jte, kts )
1072             end if
1074          end if
1075          
1076 ! accumulate vertical mixing tendencies
1077    DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1078       grid%vmix_ct(its:ite,kts:kte,jts:jte,nv) = grid%vmix_ct(its:ite,kts:kte,jts:jte,nv) + &
1079                                                         (chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv)) - &
1080                                                      chem_old(its:ite,kts:kte,jts:jte,nv))
1081    ENDDO
1083 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
1084        if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and. &
1085            (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and. &
1086            (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
1087           call wrf_debug(15,'calling chem_dbg after dry_deposition_driver')
1088           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
1089                dz8w,t_phy,p_phy,rho,chem,emis_ant,                                    &
1090                ids,ide, jds,jde, kds,kde,                                             &
1091                ims,ime, jms,jme, kms,kme,                                             &
1092                its,ite, jts,jte, kts, kte,                                            &
1093                config_flags%kemit,                                                    &
1094                grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                    &
1095                grid%ph_cl2,grid%ph_hocl,grid%ph_clno2,grid%ph_fmcl,                   &
1096                grid%ph_no3o2,                                                         &
1097                grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,      &
1098                grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,           &
1099                grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,                     &
1100                grid%ph_hcochest,grid%ph_ch3o2h,grid%ph_ch3coo2h,                      &
1101                grid%ph_ch3ono2,grid%ph_hcochob,grid%ph_n2o5,grid%ph_o2                &
1102                                                                                       )
1103        end if
1104 #endif
1106       end if
1108 ! cfrick - 2014 - WET DEPOSITION OF DUST FOLLOWING JUNG (2004)   -   mklose [03082015]
1109         if(config_flags%dustwd_onoff>0)then
1110           if(config_flags%mp_physics.ne.2 .and. config_flags%mp_physics.ne.10) then    ! mklose [03032015]
1111              write(msg,*)'CHEM_DRIVER - UoC wet deposition is not yet implemented for this & 
1112                           & microphysics option, mp_physics=', config_flags%mp_physics, & 
1113                           & ' and dustwd_onoff=', config_flags%dustwd_onoff
1114              call wrf_error_fatal( msg )
1115           endif
1117            call wrf_debug(15,'UoC dust wet deposition')
1118            call uoc_dustwd_driver(grid%precr,chem,p_phy,t_phy,                  &
1119                                   ids,ide, jds,jde, kds,kde,                    &
1120                                   ims,ime, jms,jme, kms,kme,                    &
1121                                   its,ite, jts,jte, kts,kte,                    &
1122                                   dtstepc,                                      &
1123                                   grid%dustwd_1, grid%dustwd_2,                 &
1124                                   grid%dustwd_3, grid%dustwd_4,                 &
1125                                   grid%dustwd_5,                                &
1126                                   grid%wetdep_1, grid%wetdep_2,                 &
1127                                   grid%wetdep_3, grid%wetdep_4,                 &
1128                                   grid%wetdep_5,                                &
1129                                   grid%dustwdload_1, grid%dustwdload_2,         &
1130                                   grid%dustwdload_3, grid%dustwdload_4,         &
1131                                   grid%dustwdload_5,                            &
1132                                   rri, dz8w, epsilc                             )
1133         endif
1135 !   convective transport/wet deposition
1138       if( config_flags%cu_physics>0 .and. config_flags%chem_conv_tr>0 &
1139            .and. config_flags%cu_physics/=kfcupscheme                 &
1140            .and. config_flags%cu_physics/=gfscheme ) then                !BSINGH - For WRFCuP scheme
1141         call wrf_debug(15,'calling conv transport for chemical species')
1142         if(config_flags%chem_opt >0 )then
1143         ! save old concentrations for convective tendencies
1144         DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1145            chem_old(its:ite,kts:kte,jts:jte,nv) = chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv))
1146         ENDDO
1147         call grelldrvct(grid%DT,ktau,grid%DX,                                         &
1148              rho,grid%RAINCV_B,chem,                                                  &
1149              U_phy,V_phy,t_phy,moist,dz8w,                                            &
1150              p_phy,XLV,CP,G,r_v,                                                      &
1151              z_at_w,grid%cu_co_ten,                                                   &
1152              grid%wd_no3_cu,grid%wd_so4_cu,                                 &
1153              grid%wd_nh4_cu,grid%wd_oa_cu,                                            &
1154              grid%wd_so2_cu, grid%wd_sulf_cu, grid%wd_hno3_cu, grid%wd_nh3_cu,        &
1155              grid%wd_cvasoa_cu, grid%wd_cvbsoa_cu, grid%wd_asoa_cu, grid%wd_bsoa_cu,  &
1156              grid%k22_shallow,grid%kbcon_shallow,grid%ktop_shallow,grid%xmb_shallow,  &
1157              config_flags%ishallow,num_moist,numgas,num_chem,config_flags%chem_opt,0, &
1158              config_flags%conv_tr_wetscav,config_flags%conv_tr_aqchem,                &
1159              ids,ide, jds,jde, kds,kde,                                               &
1160              ims,ime, jms,jme, kms,kme,                                               &
1161              its,ite,jts,jte,kts,k_end)
1162              if( chm_is_mozart ) then
1163                      call mozcart_lbc_set( chem, num_chem, grid%id, &
1164                                  ims, ime, jms, jme, kms, kme,    &
1165                                  its, ite, jts, jte, kts )
1166              end if
1167         ! accumulate convective tendencies
1168         DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1169           grid%conv_ct(its:ite,kts:kte,jts:jte,nv) = grid%conv_ct(its:ite,kts:kte,jts:jte,nv) + &
1170                                                             (chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv)) - &
1171                                                          chem_old(its:ite,kts:kte,jts:jte,nv))
1172         ENDDO
1173         endif
1174         if (config_flags%tracer_opt > 0)then
1175         call wrf_debug(15,'calling conv transport for tracers')
1176         call grelldrvct(grid%DT,ktau,grid%DX,                    &
1177              rho,grid%RAINCV_B,tracer,                               &
1178              U_phy,V_phy,t_phy,moist,dz8w,                                            &
1179              p_phy,XLV,CP,G,r_v,                                                      &
1180              z_at_w, grid%cu_co_ten,                                                  &
1181              grid%wd_no3_cu,grid%wd_so4_cu,                                 &
1182              grid%wd_nh4_cu,grid%wd_oa_cu,                                            &
1183              grid%wd_so2_cu, grid%wd_sulf_cu, grid%wd_hno3_cu, grid%wd_nh3_cu,        &
1184              grid%wd_cvasoa_cu, grid%wd_cvbsoa_cu, grid%wd_asoa_cu, grid%wd_bsoa_cu,  &
1185              grid%k22_shallow,grid%kbcon_shallow,grid%ktop_shallow,grid%xmb_shallow,  &
1186              config_flags%ishallow,num_moist,0,num_tracer,0,config_flags%tracer_opt,  &
1187              config_flags%conv_tr_wetscav,config_flags%conv_tr_aqchem,                &
1188              ids,ide, jds,jde, kds,kde,                                               &
1189              ims,ime, jms,jme, kms,kme,                                               &
1190              its,ite,jts,jte,kts,k_end)
1193           end if
1194      end if
1199      n2o5_het(its:ite,kts:kte,jts:jte)=0.
1200 ! Calculate rate of n2o5 hydrolysis 
1201      call wrf_debug(15,'calling calc_het_n2o5')
1203      write(msg,'(''chem_driver('',i2.2,''): Calling dchm_tstep_init'')') grid%id
1204      call wrf_debug( 200,trim(msg) )
1205      call trajectory_dchm_tstep_init( grid, do_chemstep )
1208 ! For the chemistry tracer mode, only emissions and vertical mixing are done.
1209 ! So, finish any remaining tiles and then skip to the end of chem_driver.
1211      if( do_chemstep .and.                           &
1212           config_flags%chem_opt /= CHEM_TRACER .and. &
1213           config_flags%chem_opt /= CHEM_VASH .and. &
1214           config_flags%chem_opt /= CHEM_VOLC .and. &
1215           config_flags%chem_opt /= CHEM_VOLC_4BIN .and. &
1216           config_flags%chem_opt /= DUST .and. &
1217           config_flags%chem_opt /= CHEM_TRACE2 .and. &
1218           config_flags%chem_opt /= CO2_TRACER  .and. &
1219           config_flags%chem_opt /= GHG_TRACER ) then
1222 ! chemical mechanisms
1225 ! save old concentrations for chemistry tendencies
1226    DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1227       chem_old(its:ite,kts:kte,jts:jte,nv) = chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv))
1228    ENDDO
1229         if ( cam_mam_aerosols ) &
1230            del_h2so4_gasprod(:,:,:) = chem(:,:,:,p_sulf)
1232         if(config_flags%gaschem_onoff>0)then
1234           call mechanism_driver(grid%id,curr_secs,ktau,grid%dt,grid%ktauc,dtstepc,config_flags, &
1235               grid%gmt,ijulian,t_phy,moist,p8w,t8w,grid%gd_cldfr,                     &
1236               p_phy,chem,rho,dz8w,grid%dx,g,                                          &
1237               zmid,z_at_w,grid%xlat,grid%xlong,                                       &
1238               vdrog3,vcsulf_old,vcso2_old,vch2o2_old,grid%ttday,grid%tcosz,           &
1239               grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                     &
1240               grid%ph_cl2,grid%ph_hocl,grid%ph_clno2,grid%ph_fmcl,                    &
1241               grid%ph_no3o2,                                                          &
1242               grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,       &
1243               grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,            &
1244               grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,grid%ph_hcochest,     &
1245               grid%ph_ch3o2h,grid%ph_ch3coo2h,grid%ph_ch3ono2,grid%ph_hcochob,        &
1246               grid%ph_n2o5,grid%ph_o2,grid%backg_oh,grid%backg_h2o2,grid%backg_no3,   &
1247               grid%addt,grid%addx,grid%addc,grid%etep,                                &
1248               grid%oltp,grid%olip,grid%cslp,grid%limp,grid%hc5p,grid%hc8p,grid%tolp,  &
1249               grid%xylp,grid%apip,grid%isop,grid%hc3p,grid%ethp,grid%o3p,grid%tco3,   &
1250               grid%mo2,grid%o1d,grid%olnn,grid%rpho,grid%xo2,                         &
1251               grid%ketp,grid%olnd,                                                    &
1252               ids,ide, jds,jde, kds,kde,                                              &
1253               ims,ime, jms,jme, kms,kme,                                              &
1254               its,ite,jts,jte,kts,kte        )
1257 #if ( WRF_KPP == 1 )
1259          if( config_flags%chem_opt == SAPRC99_MOSAIC_4BIN_VBS2_KPP.or. &
1260               config_flags%chem_opt ==SAPRC99_MOSAIC_8BIN_VBS2_AQ_KPP ) then
1261           do k=kts,kte
1262            do i=its,ite
1263             do j=jts,jte
1264               chem(i,k,j,p_psd1)=0.0
1265               chem(i,k,j,p_psd2)=0.0
1266             enddo
1267            enddo
1268           enddo
1269          endif
1271    CALL wrf_debug(15,'calling kpp_mechanism_driver')
1273 CALL kpp_mechanism_driver (chem,                                                      &
1274    grid%id,dtstepc,config_flags,                                                      &
1275    p_phy,t_phy,rho,                                                                   &
1276    moist,aero_srf_area,                                                               &
1277    vdrog3, ldrog, vdrog3_vbs, ldrog_vbs,                                              &
1279 #include "call_to_kpp_mech_drive.inc"
1281    ids,ide, jds,jde, kds,kde,                                                         &
1282    ims,ime, jms,jme, kms,kme,                                                         &
1283    its,ite,jts,jte,kts,kte,grid%id,num_irr_diag,irr_rates)
1284           if( chm_is_mozart ) then
1285              call mozcart_lbc_set( chem, num_chem, grid%id, &
1286                                    ims, ime, jms, jme, kms, kme,    &
1287                                     its, ite, jts, jte, kts )
1288           end if
1289     IF( grid%irr_opt == 1 ) then
1290       select case( config_flags%chem_opt )
1291         case( mozcart_kpp )
1292           do n=param_first_scalar,num_irr_diag
1293             do j=jts,jte
1294               do k=kts,kte
1295                 irr_diag_mozcart(its:ite,k,j,n) = irr_diag_mozcart(its:ite,k,j,n) &
1296                  + dtstepc*irr_rates(its:ite,k,j,n)*1.e6/(dens2con_a*rho(its:ite,k,j))
1297               enddo
1298             enddo
1299           enddo
1300         case( t1_mozcart_kpp )
1301           do n=param_first_scalar,num_irr_diag
1302             do j=jts,jte
1303               do k=kts,kte
1304                 irr_diag_t1_mozcart(its:ite,k,j,n) = irr_diag_t1_mozcart(its:ite,k,j,n) &
1305                  + dtstepc*irr_rates(its:ite,k,j,n)*1.e6/(dens2con_a*rho(its:ite,k,j))
1306               enddo
1307             enddo
1308           enddo
1309         case( mozart_mosaic_4bin_kpp )
1310           do n=param_first_scalar,num_irr_diag
1311             do j=jts,jte
1312               do k=kts,kte
1313                 irr_diag_mozart_mosaic_4bin(its:ite,k,j,n) = &
1314                       irr_diag_mozart_mosaic_4bin(its:ite,k,j,n) &
1315                       + dtstepc*irr_rates(its:ite,k,j,n)*1.e6/(dens2con_a*rho(its:ite,k,j))
1316               enddo
1317             enddo
1318           enddo
1319         case( mozart_mosaic_4bin_aq_kpp )
1320           do n=param_first_scalar,num_irr_diag
1321             do j=jts,jte
1322               do k=kts,kte
1323                 irr_diag_mozart_mosaic_4bin_aq(its:ite,k,j,n) = &
1324                       irr_diag_mozart_mosaic_4bin_aq(its:ite,k,j,n) &
1325                       + dtstepc*irr_rates(its:ite,k,j,n)*1.e6/(dens2con_a*rho(its:ite,k,j))
1326               enddo
1327             enddo
1328           enddo
1329       end select
1330     ENDIF
1331    
1332     if(config_flags%chem_opt == 301 ) then
1333        chem(its:ite,kts:kte,jts:jte,p_sulf)=vcsulf_old(its:ite,kts:kte,jts:jte)
1334        chem(its:ite,kts:kte,jts:jte,p_so2)=vcso2_old(its:ite,kts:kte,jts:jte)
1335 !      chem(its:ite,kts:kte,jts:jte,p_h2o2)=vch2o2_old(its:ite,kts:kte,jts:jte)
1336    endif
1338    write(msg,'(''chem_driver('',i2.2,''): Calling dchm_tstep_set'')') grid%id
1339    call wrf_debug( 200,trim(msg) )
1340    call trajectory_dchm_tstep_set( grid )
1342    IF(config_flags%conv_tr_aqchem == 0 ) THEN
1343       so2so4_selecta: SELECT CASE(config_flags%chem_opt)
1344       CASE (RADM2SORG,RADM2SORG_KPP,RACMSORG_KPP,RACM_SOA_VBS_KPP,RACM_SOA_VBS_HET_KPP)
1345          CALL wrf_debug(15,'gocart so2-so4 conversion')
1346          CALL  so2so4(0,chem,p_so2,p_sulf,p_h2o2,p_QC,T_PHY,MOIST,           &
1347               grid%qc_cu, grid%gd_cldfr, config_flags%cu_diag,                   &
1348               NUM_CHEM,NUM_MOIST,                                                &
1349               ids,ide, jds,jde, kds,kde,                                         &
1350               ims,ime, jms,jme, kms,kme,                                         &
1351               its,ite, jts,jte, kts,kte                                          )
1352       CASE DEFAULT
1353          CALL wrf_debug(15,'no gocart so2-so4 conversion')
1354       END SELECT so2so4_selecta
1355    else IF(config_flags%conv_tr_aqchem == 1 ) THEN
1356       so2so4_selectb: SELECT CASE(config_flags%chem_opt)
1357       CASE (RADM2SORG,RADM2SORG_KPP,RACMSORG_KPP,RACM_SOA_VBS_KPP,RACM_SOA_VBS_HET_KPP)
1358          CALL wrf_debug(15,'gocart so2-so4 conversion')
1359          CALL  so2so4(1,chem,p_so2,p_sulf,p_h2o2,p_QC,T_PHY,MOIST,           &
1360               grid%qc_cu, grid%gd_cldfr, config_flags%cu_diag,                 &
1361               NUM_CHEM,NUM_MOIST,                                                &
1362               ids,ide, jds,jde, kds,kde,                                         &
1363               ims,ime, jms,jme, kms,kme,                                         &
1364               its,ite, jts,jte, kts,kte                                          )
1365       CASE DEFAULT
1366          CALL wrf_debug(15,'no gocart so2-so4 conversion')
1367       END SELECT so2so4_selectb
1369    ENDIF
1370 #else
1371     if(config_flags%chem_opt == 301 ) then
1372        chem(its:ite,kts:kte,jts:jte,p_sulf)=vcsulf_old(its:ite,kts:kte,jts:jte)
1373        chem(its:ite,kts:kte,jts:jte,p_so2)=vcso2_old(its:ite,kts:kte,jts:jte)
1374 !      chem(its:ite,kts:kte,jts:jte,p_h2o2)=vch2o2_old(its:ite,kts:kte,jts:jte)
1375    endif
1377    IF(config_flags%conv_tr_aqchem == 0 ) THEN
1378       so2so4_selectc: SELECT CASE(config_flags%chem_opt)
1379       CASE (RADM2SORG)
1380          CALL wrf_debug(15,'gocart so2-so4 conversion')
1381          CALL  so2so4(0,chem,p_so2,p_sulf,p_h2o2,p_QC,T_PHY,MOIST,           &
1382               grid%qc_cu, grid%gd_cldfr, config_flags%cu_diag,                 &
1383               NUM_CHEM,NUM_MOIST,                                                &
1384               ids,ide, jds,jde, kds,kde,                                         &
1385               ims,ime, jms,jme, kms,kme,                                         &
1386               its,ite, jts,jte, kts,kte                                          )
1387       CASE DEFAULT
1388          CALL wrf_debug(15,'no gocart so2-so4 conversion')
1389       END SELECT so2so4_selectc
1390    else IF(config_flags%conv_tr_aqchem == 1 ) THEN
1391       so2so4_selectd: SELECT CASE(config_flags%chem_opt)
1392       CASE (RADM2SORG)
1393          CALL wrf_debug(15,'gocart so2-so4 conversion')
1394          CALL  so2so4(1,chem,p_so2,p_sulf,p_h2o2,p_QC,T_PHY,MOIST,           &
1395               grid%qc_cu, grid%gd_cldfr, config_flags%cu_diag,                   &
1396               NUM_CHEM,NUM_MOIST,                                                &
1397               ids,ide, jds,jde, kds,kde,                                         &
1398               ims,ime, jms,jme, kms,kme,                                         &
1399               its,ite, jts,jte, kts,kte                                          )
1400       CASE DEFAULT
1401          CALL wrf_debug(15,'no gocart so2-so4 conversion')
1402       END SELECT so2so4_selectd
1403    ENDIF
1406 #endif
1407 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
1408        if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and.                          &
1409            (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and.                          &
1410            (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
1411           call wrf_debug(15,'calling chem_dbg after mechanism_driver')
1412           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
1413                dz8w,t_phy,p_phy,rho,chem,emis_ant,                                    &
1414                ids,ide, jds,jde, kds,kde,                                             &
1415                ims,ime, jms,jme, kms,kme,                                             &
1416                its,ite, jts,jte, kts,kte,                                             &
1417                config_flags%kemit,                                                    &
1418                grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                    &
1419                grid%ph_cl2,grid%ph_hocl,grid%ph_clno2,grid%ph_fmcl,                   &
1420                grid%ph_no3o2,                                                         &
1421                grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,      &
1422                grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,           &
1423                grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,                     &
1424                grid%ph_hcochest,grid%ph_ch3o2h,grid%ph_ch3coo2h,                      &
1425                grid%ph_ch3ono2,grid%ph_hcochob,grid%ph_n2o5,grid%ph_o2                &
1426                                                                                       )
1427        end if
1428 #endif
1429         endif !gaschem_onoff
1431 ! accumulate gas phase chemistry tendencies
1432    DO nv=PARAM_FIRST_SCALAR,num_chem_ct
1433       grid%chem_ct(its:ite,kts:kte,jts:jte,nv) = grid%chem_ct(its:ite,kts:kte,jts:jte,nv) + &
1434                                                         (chem(its:ite,kts:kte,jts:jte,chem_ct_indices(nv)) - &
1435                                                      chem_old(its:ite,kts:kte,jts:jte,nv))
1436    ENDDO
1437 if ( cam_mam_aerosols ) &
1438            del_h2so4_gasprod(:,:,:) = chem(:,:,:,p_sulf) - del_h2so4_gasprod(:,:,:)
1440 ! initialize gas_aqfrac if either cldchem or wetscav is on
1441         if ( (config_flags%cldchem_onoff > 0) .or.                                    &
1442              (config_flags%wetscav_onoff > 0) ) then
1443             gas_aqfrac(its:ite,kts:kte,jts:jte,:) = 0.0
1444         end if
1445         
1446         !BSINGH -  The following call is for wet deposition of gases for MAM
1447         !aerosol model ONLY. In other WRF aerosol models wet deposition of gasses
1448         !is handled after dry deposition. This call doesn't belong here as there should 
1449         !be a formal driver for gas wet deposition and this call should be called
1450         !from that driver (something for future mods)
1452         !**IMPORTANT-WARNING**: We are going to operate on only 6 gases (following MOSAIC -
1453         !module_mosaic_wetscav.F). MOSAIC actually operates upon 7 gases but MAM
1454         !doesn't have msa, therefore MAM will operate on only 6 gases        
1455         !so2,h2o2,h2so4,hno3,hcl,nh3 (msa not included)
1456         
1457         !For adding more gasses, look into phys/module_cam_support.F
1458         if (((config_flags%chem_opt == CBMZ_CAM_MAM3_NOAQ) .or. &
1459              (config_flags%chem_opt == CBMZ_CAM_MAM3_AQ  ) .or. &
1460              (config_flags%chem_opt == CBMZ_CAM_MAM7_NOAQ) .or. &
1461              (config_flags%chem_opt == CBMZ_CAM_MAM7_AQ  )).and. &
1462              (config_flags%cu_physics == CAMZMSCHEME)) then
1463            
1464            call cam_mam_gas_wetdep_driver(                &
1465                 !Intent in-outs
1466                 chem,                                     &
1467                 !Intent ins
1468                 dtstepc, config_flags, grid%ht,grid%XLAT, &
1469                 grid%nevapr3d, grid%rprdsh,               &
1470                 grid%rprddp3d, grid%prain3d, grid%z,      & 
1471                 p_phy, t_phy, grid%alt, moist, scalar,    &
1472                 ids,ide, jds,jde, kds,kde,                &
1473                 ims,ime, jms,jme, kms,kme,                &
1474                 its,ite, jts,jte, kts,kte                 )
1475            
1476         endif
1477         !BSINGH - For WRFCuP scheme
1478         
1479         !
1480         !   now do convective cloud processing for cup cumulus scheme
1481         !
1482         !   issues for future work
1483         !      convective cloud processing should be done at same place (in the code)
1484         !         and same time frequency for cup and grell schemes,
1485         !         but doing it here is ok for initial testing
1486         !      there should be a "convchem_driver" or "ctrans_driver" routine
1487         !         to interface to the various options
1488         !
1489         if( config_flags%chem_conv_tr>0 .and. &
1490              config_flags%cu_physics==kfcupscheme ) then
1491            
1492            call chem_cup_driver(                                                        &
1493                 grid%id, ktau, grid%ktauc, grid%dt, dtstepc, config_flags,              &
1494                 t_phy, p_phy, rho, rri, dz8w, zmid, z_at_w,                             &
1495 #if (NMM_CORE == 1)
1496                 moist_trans, grid%cldfra, grid%ph_no2,                                  &
1497 #endif
1498 #if (EM_CORE == 1)
1499                 moist, grid%cldfra, grid%ph_no2,                                        &
1500 #endif
1501                 chem, grid%chem_cupflag,                                                &
1502                 grid%cupflag, grid%shall, grid%tcloud_cup, grid%nca, grid%wact_cup,     &
1503                 grid%cldfra_cup, grid%updfra_cup, grid%qc_ic_cup, grid%qc_iu_cup,       &
1504                 grid%mfup_cup, grid%mfup_ent_cup, grid%mfdn_cup, grid%mfdn_ent_cup,     &
1505                 grid%fcvt_qc_to_pr_cup, grid%fcvt_qc_to_qi_cup, grid%fcvt_qi_to_pr_cup, &
1506                 !BSINGH(12/03/2013): Commenting out *_ams_* variables and replacing them with 
1507                 !*_ic_cup* variables
1508                 !grid%so4_a_ams_ic_cup, grid%so4_cw_ams_ic_cup,                         &
1509                 !grid%nh4_a_ams_ic_cup, grid%nh4_cw_ams_ic_cup,                         &
1510                 !grid%no3_a_ams_ic_cup, grid%no3_cw_ams_ic_cup,                         &
1511                 !grid%oa_a_ams_ic_cup,  grid%oa_cw_ams_ic_cup,                          &
1512                 grid%co_a_ic_cup,       grid%hno3_a_ic_cup,                             &
1513                 grid%so4_a_1to4_ic_cup, grid%so4_cw_1to4_ic_cup,                        &
1514                 grid%nh4_a_1to4_ic_cup, grid%nh4_cw_1to4_ic_cup,                        &
1515                 grid%no3_a_1to4_ic_cup, grid%no3_cw_1to4_ic_cup,                        &
1516                 grid%oa_a_1to4_ic_cup,  grid%oa_cw_1to4_ic_cup,                         &
1517                 grid%oin_a_1to4_ic_cup, grid%oin_cw_1to4_ic_cup,                        &
1518                 grid%bc_a_1to4_ic_cup,  grid%bc_cw_1to4_ic_cup,                         &
1519                 grid%na_a_1to4_ic_cup,  grid%na_cw_1to4_ic_cup,                         &
1520                 grid%cl_a_1to4_ic_cup,  grid%cl_cw_1to4_ic_cup,                         &
1521                 grid%water_1to4_ic_cup,                                                 &
1522                 grid%so4_a_5to6_ic_cup, grid%so4_cw_5to6_ic_cup,                        &
1523                 grid%nh4_a_5to6_ic_cup, grid%nh4_cw_5to6_ic_cup,                        &
1524                 grid%no3_a_5to6_ic_cup, grid%no3_cw_5to6_ic_cup,                        &
1525                 grid%oa_a_5to6_ic_cup,  grid%oa_cw_5to6_ic_cup,                         &
1526                 grid%oin_a_5to6_ic_cup, grid%oin_cw_5to6_ic_cup,                        &
1527                 grid%bc_a_5to6_ic_cup,  grid%bc_cw_5to6_ic_cup,                         &
1528                 grid%na_a_5to6_ic_cup,  grid%na_cw_5to6_ic_cup,                         &
1529                 grid%cl_a_5to6_ic_cup,  grid%cl_cw_5to6_ic_cup,                         &
1530                 grid%water_5to6_ic_cup,                                                 & 
1531                 ids,ide, jds,jde, kds,kde,                                              &
1532                 ims,ime, jms,jme, kms,kme,                                              &
1533                 its,ite, jts,jte, kts,kte                                               )
1534            
1535         else
1536            grid%chem_cupflag = 0
1537         endif
1538         
1539         
1540         !BSINGH - Ends
1541         
1542         
1545 !   now do cloud chemistry
1547         if (config_flags%cldchem_onoff > 0) then
1549         call cloudchem_driver(                                                        &
1550                grid%id, ktau, grid%ktauc, grid%dt, dtstepc, config_flags,             &
1551                t_phy, p_phy, rho, rri, dz8w,                                          &
1552                p8w,grid%prain3d,scalar,grid%dvmrdt_sv13d,grid%dvmrcwdt_sv13d,         &
1553                grid%f_ice_phy,grid%f_rain_phy,grid%cldfrai, grid%cldfral,             &
1554                moist, grid%cldfra, grid%cldfra_mp_all, grid%ph_no2,                   &
1555                chem, gas_aqfrac, numgas_mam,grid%is_CAMMGMP_used,                     &
1556                grid%ph_cw,                                                            &
1557                ids,ide, jds,jde, kds,kde,                                             &
1558                ims,ime, jms,jme, kms,kme,                                             &
1559                its,ite, jts,jte, kts,kte                                              )
1561        endif
1565 !   now do aerosols
1567         if(config_flags%aerchem_onoff>0)then
1569         call aerosols_driver (grid%id,curr_secs,ktau,grid%dt,grid%ktauc,              &
1570              config_flags,dtstepc,grid%dx,                                            &
1571               rri,t_phy,moist,grid%aerwrf,p8w,t8w,                                    &
1572               p_phy,chem,rho,dz8w, rh,                                                & 
1573               zmid,z_at_w,pbl_h,grid%cldfra,grid%cldfra_mp_all,grid%vbs_nbin,         &
1574               grid%gamn2o5,grid%cn2o5,grid%kn2o5,grid%yclno2,grid%snu,grid%sac,       &
1575               grid%h2oaj,grid%h2oai,grid%nu3,grid%ac3,grid%cor3,grid%asulf,           &
1576               grid%ahno3,grid%anh3,grid%cvaro1,grid%cvaro2,grid%cvalk1,grid%cvole1,   &
1577               grid%cvapi1,grid%cvapi2,grid%cvlim1,grid%cvlim2,vcsulf_old,             &
1578               vdrog3,vdrog3_vbs,grid%br_rto,grid%dgnum4d,grid%dgnumwet4d,wetdens_ap,  &
1579               del_h2so4_gasprod,grid%dvmrdt_sv13d,grid%dvmrcwdt_sv13d,                &
1580               grid%is_CAMMGMP_used,                                                   &
1581               grid%ph_aer01, grid%ph_aer02, grid%ph_aer03, grid%ph_aer04,             &
1582               ids,ide, jds,jde, kds,kde,                                              &
1583               ims,ime, jms,jme, kms,kme,                                              &
1584               its,ite,jts,jte,kts,kte                                                 )
1586 #if (defined(CHEM_DBG_I) && defined(CHEM_DBG_J) && defined(CHEM_DBG_K))
1587        if( (its <= CHEM_DBG_I .and. ite >= CHEM_DBG_I) .and.                          &
1588            (jts <= CHEM_DBG_J .and. jte >= CHEM_DBG_J) .and.                          &
1589            (k_start <= CHEM_DBG_K .and. k_end >= CHEM_DBG_K)  ) then
1590           call wrf_debug(15,'calling chem_dbg after aerosols_driver')
1591           call chem_dbg(CHEM_DBG_I,CHEM_DBG_J,CHEM_DBG_K,grid%dt,ktau,                &
1592                dz8w,t_phy,p_phy,rho,chem,emis_ant,                                    &
1593                ids,ide, jds,jde, kds,kde,                                             &
1594                ims,ime, jms,jme, kms,kme,                                             &
1595                its,ite, jts,jte, kts, kte,                                            &
1596                config_flags%kemit,                                                    &
1597                grid%ph_macr,grid%ph_o31d,grid%ph_o33p,grid%ph_no2,                    &
1598                grid%ph_cl2,grid%ph_hocl,grid%ph_clno2,grid%ph_fmcl,                   &
1599                grid%ph_no3o2,                                                         &
1600                grid%ph_no3o,grid%ph_hno2,grid%ph_hno3,grid%ph_hno4,grid%ph_h2o2,      &
1601                grid%ph_ch2or,grid%ph_ch2om,grid%ph_ch3cho,grid%ph_ch3coch3,           &
1602                grid%ph_ch3coc2h5,grid%ph_hcocho,grid%ph_ch3cocho,                     &
1603                grid%ph_hcochest,grid%ph_ch3o2h,grid%ph_ch3coo2h,                      &
1604                grid%ph_ch3ono2,grid%ph_hcochob,grid%ph_n2o5,grid%ph_o2                &
1605                                                                                       )
1606        end if
1607 #endif
1608        endif
1611         if (config_flags%wetscav_onoff > 0) then
1612         call wetscav_driver (grid%id,ktau,grid%dt,grid%ktauc,config_flags,dtstepc,    &
1613               rri,t_phy,moist,p8w,t8w,                                                &
1614               grid%dx, grid%dy,                                                       &
1615               p_phy,chem,rho,grid%cldfra,grid%cldfra2,                                &
1616               wetscav_frcing(ims,kms,jms,p_rainprod),                                 &
1617               wetscav_frcing(ims,kms,jms,p_evapprod),                                 &
1618               grid%hno3_col_mdel,                                                     &                
1619               grid%qlsink,grid%precr,grid%preci,grid%precs,grid%precg,                &
1620               grid%wdflx,                                                             &
1621               gas_aqfrac, numgas_mam,dz8w,                                            &
1622               grid%h2oaj,grid%h2oai,grid%nu3,grid%ac3,grid%cor3,                      &
1623               grid%asulf,grid%ahno3,grid%anh3,grid%cvaro1,grid%cvaro2,                &
1624               grid%cvalk1,grid%cvole1,grid%cvapi1,grid%cvapi2,                        &
1625               grid%cvlim1,grid%cvlim2,                                                &
1626               grid%wd_no3_sc, grid%wd_so4_sc, grid%wd_nh4_sc,grid%wd_oa_sc,           &
1627               grid%wd_so2_sc, grid%wd_sulf_sc, grid%wd_hno3_sc, grid%wd_nh3_sc,       &
1628               grid%wd_cvasoa_sc, grid%wd_cvbsoa_sc, grid%wd_asoa_sc, grid%wd_bsoa_sc, &
1629               grid%qv_b4mp, grid%qc_b4mp, grid%qi_b4mp, grid%qs_b4mp,                 &
1630 !======================================================================================
1631 !Variables required for CAM_MAM_WETSCAV
1632               grid%p_hyd,scalar,grid%dgnum4d,grid%dgnumwet4d,grid%dlf,grid%dlf2,      &
1633               grid%qme3d,grid%prain3d,grid%nevapr3d,grid%rate1ord_cw2pr_st3d,         &
1634               grid%shfrc3d,grid%cmfmc,grid%cmfmc2,grid%evapcsh,grid%icwmrsh,          &
1635               grid%rprdsh,grid%evapcdp3d,grid%icwmrdp3d,grid%rprddp3d,grid%fracis3d,  &               
1636               grid%f_ice_phy,grid%f_rain_phy,grid%cldfrai,grid%cldfral,               &
1637               grid%cldfra_mp_all,grid%is_CAMMGMP_used,                                &
1638 !======================================================================================
1639               ids,ide, jds,jde, kds,kde,                                              &
1640               ims,ime, jms,jme, kms,kme,                                              &
1641               its,ite, jts,jte, kts,kte)                                              
1642        
1643        endif
1645        !BSINGH - zm_conv_tend_2 call HAS to be called after wetscavenging call as it uses 'fracis3d'
1646        !variable from wetscav driver.
1647        if (((config_flags%chem_opt == CBMZ_CAM_MAM3_NOAQ) .or. &
1648             (config_flags%chem_opt == CBMZ_CAM_MAM3_AQ  ) .or. &
1649             (config_flags%chem_opt == CBMZ_CAM_MAM7_NOAQ) .or. &
1650             (config_flags%chem_opt == CBMZ_CAM_MAM7_AQ  )).and. &
1651             (config_flags%cu_physics == CAMZMSCHEME)) then
1652           !BSINGH - CALL ZM cumulus scheme for convective transport of the aerosols
1653           call zm_conv_tend_2(grid%itimestep, grid%dt, p8w, grid%fracis3d, grid%dp3d,    &
1654                grid%du3d, grid%ed3d, grid%eu3d, grid%md3d, grid%mu3d, grid%dsubcld2d,    &
1655                grid%ideep2d, grid%jt2d, grid%maxg2d, grid%lengath2d, moist, scalar, chem,&
1656                ids,ide, jds,jde, kds,kde, &
1657                ims,ime, jms,jme, kms,kme, &
1658                its,ite, jts,jte, kts,kte)     
1659        endif
1661 end if !Chemistry time step check
1663 ! do so2 to sulf conversion for full volc case
1664 ! since we do not have h2o2 as a variable, pass in p_h2o2 as zero
1665 ! will have to use backgrund value
1667         if(config_flags%chem_opt == CHEM_VOLC)then
1668           CALL wrf_debug(15,'gocart so2-so4 conversion')
1669           CALL  so2so4(0,chem,p_so2,p_sulf,p_h2o2,p_QC,T_PHY,MOIST,         &
1670          grid%qc_cu, grid%gd_cldfr, config_flags%cu_diag,                   &
1671          NUM_CHEM,NUM_MOIST,                                                &
1672          ids,ide, jds,jde, kds,kde,                                         &
1673          ims,ime, jms,jme, kms,kme,                                         &
1674          its,ite, jts,jte, kts,kte                                          )
1675         endif
1677 !   now do wet removal; first LS if there is no explicit aqeous phase
1679         if(config_flags%wetscav_onoff<0)then
1680            call wrf_debug(15,'calculate LS wet deposition')
1681            call wetdep_ls(grid%dt,chem,grid%rainncv,moist,rho,num_moist, &
1682                 num_chem,numgas,dz8w,vvel,grid%chem_opt,                 &
1683                 ids,ide, jds,jde, kds,kde,                               &
1684                 ims,ime, jms,jme, kms,kme,                               &
1685                 its,ite, jts,jte, kts,kte                                )
1686          endif
1688 ! Sum up the aerosol mass for radiation and diagnostic purposes. Unlike
1689 ! aerosol_driver, which is called every dtchem, this must be done every
1690 ! time step because of emissions and deposition.
1692       call sum_pm_driver ( config_flags,                                              &
1693            rri, chem, grid%h2oaj, grid%h2oai,                                         &
1694            grid%pm2_5_dry, grid%pm2_5_water, grid%pm2_5_dry_ec, grid%pm10,            &
1695            grid%tsoa,grid%asoa,grid%bsoa,                                             &
1696            grid%hoa_a01,grid%hoa_a02,grid%hoa_a03,grid%hoa_a04,                       &
1697            grid%hoa_a05,     grid%hoa_a06,     grid%hoa_a07,      grid%hoa_a08,       & !BSINGH(12/03/2013) Added 4 more bins for each species for SAPRC 8 bin version
1698            grid%bboa_a01,grid%bboa_a02,grid%bboa_a03,grid%bboa_a04,                   &
1699            grid%bboa_a05,    grid%bboa_a06,    grid%bboa_a07,     grid%bboa_a08,      &
1700            grid%soa_a01,grid%soa_a02,grid%soa_a03,grid%soa_a04,                       &
1701            grid%soa_a05,     grid%soa_a06,     grid%soa_a07,      grid%soa_a08,       &
1702            grid%bbsoa_a01,grid%bbsoa_a02,grid%bbsoa_a03,grid%bbsoa_a04,               &
1703            grid%bbsoa_a05,   grid%bbsoa_a06,   grid%bbsoa_a07,    grid%bbsoa_a08,     &
1704            grid%hsoa_a01,grid%hsoa_a02,grid%hsoa_a03,grid%hsoa_a04,                   &
1705            grid%hsoa_a05,    grid%hsoa_a06,    grid%hsoa_a07,     grid%hsoa_a08,      &
1706            grid%biog_a01,grid%biog_a02,grid%biog_a03,grid%biog_a04,                   &
1707            grid%biog_a05,    grid%biog_a06,    grid%biog_a07,     grid%biog_a08,      &
1708            grid%asmpsoa_a01,grid%asmpsoa_a02,grid%asmpsoa_a03,grid%asmpsoa_a04,       &
1709            grid%arosoa_a01,grid%arosoa_a02,grid%arosoa_a03,grid%arosoa_a04,           &
1710            grid%arosoa_a05,  grid%arosoa_a06,  grid%arosoa_a07,   grid%arosoa_a08,    &
1711            grid%totoa_a01,grid%totoa_a02,grid%totoa_a03,grid%totoa_a04,               &
1712            grid%totoa_a05,   grid%totoa_a06,   grid%totoa_a07,    grid%totoa_a08,     &
1713            grid%hsoa_c,grid%hsoa_o,grid%bbsoa_c,grid%bbsoa_o,                         &
1714            grid%biog_v1,grid%biog_v2,grid%biog_v3,grid%biog_v4,                       &
1715            grid%ant_v1,grid%ant_v2,grid%ant_v3,grid%ant_v4,                           &
1716            grid%smpa_v1,grid%smpbb_v1,                                                &
1717             !BSINGH(12/03/2013) - Added cw aerosols(for VBS)
1718            grid%hoa_cw01,    grid%hoa_cw02,    grid%hoa_cw03,    grid%hoa_cw04,       &
1719            grid%hoa_cw05,    grid%hoa_cw06,    grid%hoa_cw07,    grid%hoa_cw08,       &
1720            grid%bboa_cw01,   grid%bboa_cw02,   grid%bboa_cw03,   grid%bboa_cw04,      &
1721            grid%bboa_cw05,   grid%bboa_cw06,   grid%bboa_cw07,   grid%bboa_cw08,      &
1722            grid%soa_cw01,    grid%soa_cw02,    grid%soa_cw03,    grid%soa_cw04,       &
1723            grid%soa_cw05,    grid%soa_cw06,    grid%soa_cw07,    grid%soa_cw08,       &
1724            grid%bbsoa_cw01,  grid%bbsoa_cw02,  grid%bbsoa_cw03,  grid%bbsoa_cw04,     &
1725            grid%bbsoa_cw05,  grid%bbsoa_cw06,  grid%bbsoa_cw07,  grid%bbsoa_cw08,     &
1726            grid%biog_cw01,   grid%biog_cw02,   grid%biog_cw03,   grid%biog_cw04,      &
1727            grid%biog_cw05,   grid%biog_cw06,   grid%biog_cw07,   grid%biog_cw08,      &
1728            grid%hsoa_cw01,   grid%hsoa_cw02,   grid%hsoa_cw03,   grid%hsoa_cw04,      &
1729            grid%hsoa_cw05,   grid%hsoa_cw06,   grid%hsoa_cw07,   grid%hsoa_cw08,      &
1730            grid%arosoa_cw01, grid%arosoa_cw02, grid%arosoa_cw03, grid%arosoa_cw04,    &
1731            grid%arosoa_cw05, grid%arosoa_cw06, grid%arosoa_cw07, grid%arosoa_cw08,    &
1732            grid%totoa_cw01,  grid%totoa_cw02,  grid%totoa_cw03,  grid%totoa_cw04,     &
1733            grid%totoa_cw05,  grid%totoa_cw06,  grid%totoa_cw07,  grid%totoa_cw08,     &        
1734            grid%hsoa_cw_c,   grid%hsoa_cw_o,   grid%bbsoa_cw_c,  grid%bbsoa_cw_o,     &
1735            grid%biog_cw_v1,                                                           &
1736            grid%ant_cw_v1,                                                            &
1737            !BSINGH(12/03/2013) -ENDS
1738              ids,ide, jds,jde, kds,kde,                                               &
1739              ims,ime, jms,jme, kms,kme,                                               &
1740              its,ite, jts,jte, kts,kte             )
1741              
1742      call dust_load_driver ( config_flags,                                            &
1743            rri, chem, dz8w, grid%dustload_1, grid%dustload_2, grid%dustload_3,        &
1744            grid%dustload_4, grid%dustload_5,                                          &
1745            ids,ide, jds,jde, kds,kde,                                                 &
1746            ims,ime, jms,jme, kms,kme,                                                 &
1747            its,ite, jts,jte, kts, kte                                                 )
1750 ! Fill top level to prevent spurious interpolation results (no extrapolation)
1751       do nv=1,num_chem
1752          do j=jts,jte
1753             do i=its,ite
1754                   chem(i,k_end,j,nv)=chem(i,kte,j,nv)
1755             enddo
1756          enddo
1757       enddo
1758       call wrf_debug(15,'done tileloop in chem_driver')
1759    if( grid%OPT_PARS_OUT == 1) then
1760       call wrf_debug(15,'calculate optical output stuff')
1761       call aer_opt_out(TAUAER300=grid%tauaer1, TAUAER400=grid%tauaer2                            & 
1762      &        ,TAUAER600=grid%tauaer3, TAUAER999=grid%tauaer4                                    &
1763      &        ,GAER300=grid%gaer1, GAER400=grid%gaer2, GAER600=grid%gaer3, GAER999=grid%gaer4    &
1764      &        ,WAER300=grid%waer1, WAER400=grid%waer2, WAER600=grid%waer3, WAER999=grid%waer4    &
1765               ,ext_coeff=grid%ext_coef,bscat_coeff=grid%bscat_coef,asym_par=grid%asym_par        &
1766               ,num_ext_coef=num_ext_coef,num_bscat_coef=num_bscat_coef,num_asym_par=num_asym_par &
1767      &        ,dz8w=dz8w                                                                         &
1768      &        ,IDS=ids,IDE=ide, JDS=jds,JDE=jde, KDS=kds,KDE=kde                                 &
1769      &        ,IMS=ims,IME=ime, JMS=jms,JME=jme, KMS=kms,KME=kme                                 &
1770      &        ,its=its,ite=ite,jts=jts,jte=jte,kts=kts, kte=kte                                  )
1773    endif
1774    tracer2: SELECT CASE(config_flags%tracer_opt)
1775     CASE (TRACER_TEST1, TRACER_TEST2, TRACER_TEST3)
1776        CALL wrf_debug(15,'tracer mode: reset some tracers')
1777        call set_tracer(grid%dt,ktau,pbl_h,tracer,t_phy,         &
1778                         config_flags%tracer_opt,num_tracer,     &
1779                        zmid,grid%ht,ids,ide, jds,jde, kds,kde,  & ! domain dims
1780                                ims,ime, jms,jme, kms,kme,       & ! memory dims
1781                                its,ite, jts,jte, kts,kte )
1782    END SELECT tracer2
1784 !-- set upper boundary condition
1785     if( config_flags%have_bcs_upper )then
1786        call wrf_debug(15,'set upper boundary condition')
1787        call tropopause_driver( grid%id, grid%dt, current_date_char,         &
1788                                t_phy, p_phy, p8w, zmid, z_at_w,             &
1789                                grid%tropo_lev, grid%tropo_p,  grid%tropo_z, &
1790                                ids, ide, jds, jde, kds, kde,                &
1791                                ims, ime, jms, jme, kms, kme,                &
1792                                its, ite, jts, jte, kts, kte                 )
1793        call upper_bc_driver  ( grid%id, grid%dt, current_date_char, &
1794                                chem, p_phy, p8w, grid%tropo_lev,    &
1795                                ids,ide, jds,jde, kds,kde,           &
1796                                ims,ime, jms,jme, kms,kme,           &
1797                                its,ite, jts,jte, kts,kte            )
1798     endif
1800    END DO chem_tile_loop_1
1802 !-- Work around for dgnum and dgnumwet not being written to restart files.
1803    
1804    grid%dgnum_a1(its:ite, kts:kte, jts:jte) = grid%dgnum4d(its:ite, kts:kte, jts:jte, 1)
1805    grid%dgnum_a2(its:ite, kts:kte, jts:jte) = grid%dgnum4d(its:ite, kts:kte, jts:jte, 2)
1806    grid%dgnum_a3(its:ite, kts:kte, jts:jte) = grid%dgnum4d(its:ite, kts:kte, jts:jte, 3)
1807    
1808    grid%dgnumwet_a1(its:ite, kts:kte, jts:jte) = grid%dgnumwet4d(its:ite, kts:kte, jts:jte, 1)
1809    grid%dgnumwet_a2(its:ite, kts:kte, jts:jte) = grid%dgnumwet4d(its:ite, kts:kte, jts:jte, 2)
1810    grid%dgnumwet_a3(its:ite, kts:kte, jts:jte) = grid%dgnumwet4d(its:ite, kts:kte, jts:jte, 3)
1812    IF( allocated( irr_rates ) ) THEN
1813      deallocate( irr_rates )
1814    ENDIF
1817     END subroutine chem_driver