1 MODULE module_emissions_anthropogenics
2 !WRF:MODEL_LAYER:CHEMICS
6 ! currently this only adds in the emissions...
7 ! this may be done differently for different chemical mechanisms
8 ! in the future. aerosols are already added somewhere else....
10 subroutine add_anthropogenics(id,dtstep,dz8w,config_flags,rho_phy,alt, &
11 chem, emis_ant, emis_aircraft, &
12 ids,ide, jds,jde, kds,kde, &
13 ims,ime, jms,jme, kms,kme, &
14 its,ite, jts,jte, kts,kte )
15 !----------------------------------------------------------------------
17 USE module_state_description
19 USE module_data_sorgam, only : mw_so4_aer
20 USE module_model_constants, only : mwdry
24 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
26 INTEGER, INTENT(IN ) :: id, &
27 ids,ide, jds,jde, kds,kde, &
28 ims,ime, jms,jme, kms,kme, &
29 its,ite, jts,jte, kts,kte
30 REAL, INTENT(IN ) :: &
33 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
34 INTENT(INOUT ) :: chem
38 ! REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
39 REAL, DIMENSION( ims:ime, kms:config_flags%kemit, jms:jme,num_emis_ant ), &
42 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
43 INTENT(IN ) :: rho_phy, &
51 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , &
55 ! stuff for aircraft emissions
57 REAL, DIMENSION( ims:ime, kms:config_flags%kemit_aircraft, jms:jme, num_emis_aircraft), &
58 OPTIONAL,INTENT(IN ) :: emis_aircraft
60 real, parameter :: voc_fac = .04*28./250.
62 real :: conv_rho(its:ite)
63 real :: conv_rho_aer(its:ite)
66 character(len=256) :: msg
68 !--- deposition and emissions stuff
73 ! .. Intrinsic Functions ..
75 call wrf_debug(15,'add_anthropogenics')
77 is_moz_chm = config_flags%chem_opt == MOZART_KPP .or. &
78 config_flags%chem_opt == MOZCART_KPP .or. &
79 config_flags%chem_opt == T1_MOZCART_KPP .or. &
80 config_flags%chem_opt == MOZART_MOSAIC_4BIN_KPP .or. &
81 config_flags%chem_opt == MOZART_MOSAIC_4BIN_AQ_KPP
86 k_loop: DO k=kts,min(config_flags%kemit,kte)
87 conv_rho(its:ite) = 4.828e-4/rho_phy(its:ite,k,j)*dtstep/(dz8w(its:ite,k,j)*60.)
88 conv_rho_aer(its:ite) = alt(its:ite,k,j)*dtstep/dz8w(its:ite,k,j)
89 chem(its:ite,k,j,p_so2) = chem(its:ite,k,j,p_so2) + emis_ant(its:ite,k,j,p_e_so2)*conv_rho(its:ite)
90 chem(its:ite,k,j,p_co) = chem(its:ite,k,j,p_co) + emis_ant(its:ite,k,j,p_e_co)*conv_rho(its:ite)
91 chem(its:ite,k,j,p_no) = chem(its:ite,k,j,p_no) + emis_ant(its:ite,k,j,p_e_no)*conv_rho(its:ite)
92 chem(its:ite,k,j,p_no2) = chem(its:ite,k,j,p_no2) + emis_ant(its:ite,k,j,p_e_no2)*conv_rho(its:ite)
93 chem(its:ite,k,j,p_nh3) = chem(its:ite,k,j,p_nh3) + emis_ant(its:ite,k,j,p_e_nh3)*conv_rho(its:ite)
94 chem(its:ite,k,j,p_hcl) = chem(its:ite,k,j,p_hcl) + emis_ant(its:ite,k,j,p_e_hcl)*conv_rho(its:ite)
95 chem(its:ite,k,j,p_ch3cl) = chem(its:ite,k,j,p_ch3cl) + emis_ant(its:ite,k,j,p_e_ch3cl)*conv_rho(its:ite)
96 is_mozart:if( is_moz_chm ) then
97 chem(its:ite,k,j,p_bigalk) = chem(its:ite,k,j,p_bigalk) + emis_ant(its:ite,k,j,p_e_bigalk)*conv_rho(its:ite)
98 chem(its:ite,k,j,p_bigene) = chem(its:ite,k,j,p_bigene) + emis_ant(its:ite,k,j,p_e_bigene)*conv_rho(its:ite)
99 chem(its:ite,k,j,p_c2h4) = chem(its:ite,k,j,p_c2h4) + emis_ant(its:ite,k,j,p_e_c2h4)*conv_rho(its:ite)
100 chem(its:ite,k,j,p_c2h5oh) = chem(its:ite,k,j,p_c2h5oh) + emis_ant(its:ite,k,j,p_e_c2h5oh)*conv_rho(its:ite)
101 chem(its:ite,k,j,p_c2h6) = chem(its:ite,k,j,p_c2h6) + emis_ant(its:ite,k,j,p_e_c2h6)*conv_rho(its:ite)
102 chem(its:ite,k,j,p_c3h6) = chem(its:ite,k,j,p_c3h6) + emis_ant(its:ite,k,j,p_e_c3h6)*conv_rho(its:ite)
103 chem(its:ite,k,j,p_c3h8) = chem(its:ite,k,j,p_c3h8) + emis_ant(its:ite,k,j,p_e_c3h8)*conv_rho(its:ite)
104 chem(its:ite,k,j,p_hcho) = chem(its:ite,k,j,p_hcho) + emis_ant(its:ite,k,j,p_e_ch2o)*conv_rho(its:ite)
105 chem(its:ite,k,j,p_ald) = chem(its:ite,k,j,p_ald) + emis_ant(its:ite,k,j,p_e_ch3cho)*conv_rho(its:ite)
106 chem(its:ite,k,j,p_acet) = chem(its:ite,k,j,p_acet) + emis_ant(its:ite,k,j,p_e_ch3coch3)*conv_rho(its:ite)
107 chem(its:ite,k,j,p_ch3oh) = chem(its:ite,k,j,p_ch3oh) + emis_ant(its:ite,k,j,p_e_ch3oh)*conv_rho(its:ite)
108 chem(its:ite,k,j,p_mek) = chem(its:ite,k,j,p_mek) + emis_ant(its:ite,k,j,p_e_mek)*conv_rho(its:ite)
109 chem(its:ite,k,j,p_tol) = chem(its:ite,k,j,p_tol) + emis_ant(its:ite,k,j,p_e_toluene)*conv_rho(its:ite)
110 chem(its:ite,k,j,p_isopr) = chem(its:ite,k,j,p_isopr) + emis_ant(its:ite,k,j,p_e_isop)*conv_rho(its:ite)
111 chem(its:ite,k,j,p_sulf) = chem(its:ite,k,j,p_sulf) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_sulf)*mwdry/mw_so4_aer*1.e-3
112 if( config_flags%emiss_opt == mozmem ) then
113 if( p_gly >= param_first_scalar ) then
114 chem(its:ite,k,j,p_gly) = chem(its:ite,k,j,p_gly) + emis_ant(its:ite,k,j,p_e_gly)*conv_rho(its:ite)
116 if( p_mgly >= param_first_scalar ) then
117 chem(its:ite,k,j,p_mgly) = chem(its:ite,k,j,p_mgly) + emis_ant(its:ite,k,j,p_e_mgly)*conv_rho(its:ite)
119 if( p_macr >= param_first_scalar ) then
120 chem(its:ite,k,j,p_macr) = chem(its:ite,k,j,p_macr) + emis_ant(its:ite,k,j,p_e_macr)*conv_rho(its:ite)
122 if( p_mvk >= param_first_scalar ) then
123 chem(its:ite,k,j,p_mvk) = chem(its:ite,k,j,p_mvk) + emis_ant(its:ite,k,j,p_e_mvk)*conv_rho(its:ite)
125 if( p_c2h2 >= param_first_scalar ) then
126 chem(its:ite,k,j,p_c2h2) = chem(its:ite,k,j,p_c2h2) + emis_ant(its:ite,k,j,p_e_c2h2)*conv_rho(its:ite)
128 if( p_benzene >= param_first_scalar ) then
129 chem(its:ite,k,j,p_benzene) = chem(its:ite,k,j,p_benzene) + emis_ant(its:ite,k,j,p_e_benzene)*conv_rho(its:ite)
131 if( p_xyl >= param_first_scalar ) then
132 chem(its:ite,k,j,p_xyl) = chem(its:ite,k,j,p_xyl) + emis_ant(its:ite,k,j,p_e_xylene)*conv_rho(its:ite)
134 if( p_apin >= param_first_scalar ) then
135 chem(its:ite,k,j,p_apin) = chem(its:ite,k,j,p_apin) + emis_ant(its:ite,k,j,p_e_apin)*conv_rho(its:ite)
138 if( config_flags%emiss_opt == mozem .or. config_flags%emiss_opt == mozcem ) then
139 if( p_c10h16 >= param_first_scalar ) then
140 chem(its:ite,k,j,p_c10h16) = chem(its:ite,k,j,p_c10h16) + emis_ant(its:ite,k,j,p_e_c10h16)*conv_rho(its:ite)
143 if( config_flags%chem_opt == MOZCART_KPP .or. config_flags%chem_opt == T1_MOZCART_KPP ) then
144 if( config_flags%emiss_opt == mozcem .or. config_flags%emiss_opt == mozc_t1_em ) then
145 chem(its:ite,k,j,p_p10) = chem(its:ite,k,j,p_p10) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_pm_10)
146 chem(its:ite,k,j,p_p25) = chem(its:ite,k,j,p_p25) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_pm_25)
147 chem(its:ite,k,j,p_oc1) = chem(its:ite,k,j,p_oc1) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_oc)
148 chem(its:ite,k,j,p_bc1) = chem(its:ite,k,j,p_bc1) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_bc)
150 if( config_flags%chem_opt == T1_MOZCART_KPP .and. config_flags%emiss_opt == mozc_t1_em ) then
151 chem(its:ite,k,j,p_gly) = chem(its:ite,k,j,p_gly) + emis_ant(its:ite,k,j,p_e_gly)*conv_rho(its:ite)
152 chem(its:ite,k,j,p_macr) = chem(its:ite,k,j,p_macr) + emis_ant(its:ite,k,j,p_e_macr)*conv_rho(its:ite)
153 chem(its:ite,k,j,p_mgly) = chem(its:ite,k,j,p_mgly) + emis_ant(its:ite,k,j,p_e_mgly)*conv_rho(its:ite)
154 chem(its:ite,k,j,p_c2h2) = chem(its:ite,k,j,p_c2h2) + emis_ant(its:ite,k,j,p_e_c2h2)*conv_rho(its:ite)
155 chem(its:ite,k,j,p_hcooh) = chem(its:ite,k,j,p_hcooh) + emis_ant(its:ite,k,j,p_e_hcooh)*conv_rho(its:ite)
156 chem(its:ite,k,j,p_mvk) = chem(its:ite,k,j,p_mvk) + emis_ant(its:ite,k,j,p_e_mvk)*conv_rho(its:ite)
157 chem(its:ite,k,j,p_benzene) = chem(its:ite,k,j,p_benzene) + emis_ant(its:ite,k,j,p_e_benzene)*conv_rho(its:ite)
158 chem(its:ite,k,j,p_xylenes) = chem(its:ite,k,j,p_xylenes) + emis_ant(its:ite,k,j,p_e_xylene)*conv_rho(its:ite)
160 elseif( config_flags%chem_opt == MOZART_MOSAIC_4BIN_KPP .or. config_flags%chem_opt == MOZART_MOSAIC_4BIN_AQ_KPP ) then
161 chem(its:ite,k,j,p_hcooh) = chem(its:ite,k,j,p_hcooh) + emis_ant(its:ite,k,j,p_e_hcooh)*conv_rho(its:ite)
162 chem(its:ite,k,j,p_hono) = chem(its:ite,k,j,p_hono) + emis_ant(its:ite,k,j,p_e_hono)*conv_rho(its:ite)
163 if( config_flags%chem_opt == MOZART_MOSAIC_4BIN_KPP ) then
164 ! emissions should be CO_A and CO_BB with yields instead of VOC_A and VOC_BB
165 ! chem(i,k,j,p_voca) = chem(i,k,j,p_voca) + conv_rho*emis_ant(i,k,j,p_e_voca)
166 ! chem(i,k,j,p_vocbb) = chem(i,k,j,p_vocbb) + conv_rho*emis_ant(i,k,j,p_e_vocbb)
167 ! chem(i,k,j,p_voca) = chem(i,k,j,p_voca) + conv_rho_aer*emis_ant(i,k,j,p_e_co_a)*0.04*28./250.
168 ! chem(i,k,j,p_vocbb) = chem(i,k,j,p_vocbb) + conv_rho_aer*emis_ant(i,k,j,p_e_co_bb)*0.04*28./250.
169 chem(its:ite,k,j,p_voca) = chem(its:ite,k,j,p_voca) + conv_rho(its:ite)*emis_ant(its:ite,k,j,p_e_co_a)*voc_fac
170 chem(its:ite,k,j,p_vocbb) = chem(its:ite,k,j,p_vocbb) + conv_rho(its:ite)*emis_ant(its:ite,k,j,p_e_co_bb)*voc_fac
173 else if( config_flags%chem_opt == CRIMECH_KPP .or. &
174 config_flags%chem_opt == CRI_MOSAIC_8BIN_AQ_KPP .or. &
175 config_flags%chem_opt == CRI_MOSAIC_4BIN_AQ_KPP) then
177 chem(i,k,j,p_c2h6) = chem(i,k,j,p_c2h6) + emis_ant(i,k,j,p_e_c2h6)*conv_rho(i)
178 chem(i,k,j,p_c2h4) = chem(i,k,j,p_c2h4) + emis_ant(i,k,j,p_e_c2h4)*conv_rho(i)
179 chem(i,k,j,p_c5h8) = chem(i,k,j,p_c5h8) + emis_ant(i,k,j,p_e_c5h8)*conv_rho(i)
180 chem(i,k,j,p_tm123b) = chem(i,k,j,p_tm123b) + emis_ant(i,k,j,p_e_tm123b)*conv_rho(i)
181 chem(i,k,j,p_tm124b) = chem(i,k,j,p_tm124b) + emis_ant(i,k,j,p_e_tm124b)*conv_rho(i)
182 chem(i,k,j,p_tm135b) = chem(i,k,j,p_tm135b) + emis_ant(i,k,j,p_e_tm135b)*conv_rho(i)
183 chem(i,k,j,p_oethtol) = chem(i,k,j,p_oethtol) + emis_ant(i,k,j,p_e_oethtol)*conv_rho(i)
184 chem(i,k,j,p_methtol) = chem(i,k,j,p_methtol) + emis_ant(i,k,j,p_e_methtol)*conv_rho(i)
185 chem(i,k,j,p_pethtol) = chem(i,k,j,p_pethtol) + emis_ant(i,k,j,p_e_pethtol)*conv_rho(i)
186 chem(i,k,j,p_dime35eb) = chem(i,k,j,p_dime35eb) + emis_ant(i,k,j,p_e_dime35eb)*conv_rho(i)
187 chem(i,k,j,p_hcho) = chem(i,k,j,p_hcho) + emis_ant(i,k,j,p_e_hcho)*conv_rho(i)
188 chem(i,k,j,p_ch3cho) = chem(i,k,j,p_ch3cho) + emis_ant(i,k,j,p_e_ch3cho)*conv_rho(i)
189 chem(i,k,j,p_c2h5cho) = chem(i,k,j,p_c2h5cho) + emis_ant(i,k,j,p_e_c2h5cho)*conv_rho(i)
190 chem(i,k,j,p_ket) = chem(i,k,j,p_ket) + emis_ant(i,k,j,p_e_ket)*conv_rho(i)
191 chem(i,k,j,p_mek) = chem(i,k,j,p_mek) + emis_ant(i,k,j,p_e_mek)*conv_rho(i)
192 chem(i,k,j,p_ch3oh) = chem(i,k,j,p_ch3oh) + emis_ant(i,k,j,p_e_ch3oh)*conv_rho(i)
193 chem(i,k,j,p_c2h5oh) = chem(i,k,j,p_c2h5oh) + emis_ant(i,k,j,p_e_c2h5oh)*conv_rho(i)
194 chem(i,k,j,p_c3h6) = chem(i,k,j,p_c3h6) + emis_ant(i,k,j,p_e_c3h6)*conv_rho(i)
195 chem(i,k,j,p_c2h2) = chem(i,k,j,p_c2h2) + emis_ant(i,k,j,p_e_c2h2)*conv_rho(i)
196 chem(i,k,j,p_benzene) = chem(i,k,j,p_benzene) + emis_ant(i,k,j,p_e_benzene)*conv_rho(i)
197 chem(i,k,j,p_nc4h10) = chem(i,k,j,p_nc4h10) + emis_ant(i,k,j,p_e_nc4h10)*conv_rho(i)
198 chem(i,k,j,p_toluene) = chem(i,k,j,p_toluene) + emis_ant(i,k,j,p_e_toluene)*conv_rho(i)
199 chem(i,k,j,p_oxyl) = chem(i,k,j,p_oxyl) + emis_ant(i,k,j,p_e_oxyl)*conv_rho(i)
200 chem(i,k,j,p_c3h8) = chem(i,k,j,p_c3h8) + emis_ant(i,k,j,p_e_c3h8)*conv_rho(i)
201 chem(i,k,j,p_tbut2ene) = chem(i,k,j,p_tbut2ene) + emis_ant(i,k,j,p_e_tbut2ene)*conv_rho(i)
202 chem(i,k,j,p_ch3co2h) = chem(i,k,j,p_ch3co2h) + emis_ant(i,k,j,p_e_ch3co2h)*conv_rho(i)
206 chem(i,k,j,p_csl) = chem(i,k,j,p_csl) + emis_ant(i,k,j,p_e_csl)*conv_rho(i)
207 chem(i,k,j,p_iso) = chem(i,k,j,p_iso) + emis_ant(i,k,j,p_e_iso)*conv_rho(i)
208 chem(i,k,j,p_ald) = chem(i,k,j,p_ald) + emis_ant(i,k,j,p_e_ald)*conv_rho(i)
209 chem(i,k,j,p_hcho) = chem(i,k,j,p_hcho) + emis_ant(i,k,j,p_e_hcho)*conv_rho(i)
210 chem(i,k,j,p_ora2) = chem(i,k,j,p_ora2) + emis_ant(i,k,j,p_e_ora2)*conv_rho(i)
211 chem(i,k,j,p_hc3) = chem(i,k,j,p_hc3) + emis_ant(i,k,j,p_e_hc3)*conv_rho(i)
212 chem(i,k,j,p_hc5) = chem(i,k,j,p_hc5) + emis_ant(i,k,j,p_e_hc5)*conv_rho(i)
213 chem(i,k,j,p_hc8) = chem(i,k,j,p_hc8) + emis_ant(i,k,j,p_e_hc8)*conv_rho(i)
214 chem(i,k,j,p_eth) = chem(i,k,j,p_eth) + emis_ant(i,k,j,p_e_eth)*conv_rho(i)
215 chem(i,k,j,p_olt) = chem(i,k,j,p_olt) + emis_ant(i,k,j,p_e_olt)*conv_rho(i)
216 chem(i,k,j,p_oli) = chem(i,k,j,p_oli) + emis_ant(i,k,j,p_e_oli)*conv_rho(i)
217 chem(i,k,j,p_tol) = chem(i,k,j,p_tol) + emis_ant(i,k,j,p_e_tol)*conv_rho(i)
218 chem(i,k,j,p_xyl) = chem(i,k,j,p_xyl) + emis_ant(i,k,j,p_e_xyl)*conv_rho(i)
219 chem(i,k,j,p_ket) = chem(i,k,j,p_ket) + emis_ant(i,k,j,p_e_ket)*conv_rho(i)
221 if(p_ch4 >= param_first_scalar) then
222 chem(its:ite,k,j,p_ch4) = chem(its:ite,k,j,p_ch4) + emis_ant(its:ite,k,j,p_e_ch4)*conv_rho(its:ite)
224 if(p_ol2 >= param_first_scalar) then
225 chem(its:ite,k,j,p_ol2) = chem(its:ite,k,j,p_ol2) + emis_ant(its:ite,k,j,p_e_ol2)*conv_rho(its:ite)
227 if(p_ete >= param_first_scalar) then
228 chem(its:ite,k,j,p_ete) = chem(its:ite,k,j,p_ete) + emis_ant(its:ite,k,j,p_e_ol2)*conv_rho(its:ite)
230 if(p_e_hono >= param_first_scalar .and. p_hono >= param_first_scalar) then
231 chem(its:ite,k,j,p_hono) = chem(its:ite,k,j,p_hono) + emis_ant(its:ite,k,j,p_e_hono)*conv_rho(its:ite)
233 if(p_e_terp >= param_first_scalar .and. p_api >= param_first_scalar .and. p_lim >= param_first_scalar) then
234 chem(its:ite,k,j,p_api) = chem(its:ite,k,j,p_api) + 0.5*emis_ant(its:ite,k,j,p_e_terp)*conv_rho(its:ite)
235 chem(its:ite,k,j,p_lim) = chem(its:ite,k,j,p_lim) + 0.5*emis_ant(its:ite,k,j,p_e_terp)*conv_rho(its:ite)
237 if(p_e_co2 >= param_first_scalar .and. p_co2 >= param_first_scalar) then
238 chem(its:ite,k,j,p_co2) = chem(its:ite,k,j,p_co2) + emis_ant(its:ite,k,j,p_e_co2)*conv_rho(its:ite)
240 if( config_flags%chem_opt == GOCARTRACM_KPP .or. config_flags%chem_opt == GOCARTRADM2 ) then
241 chem(its:ite,k,j,p_p10) = chem(its:ite,k,j,p_p10) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_pm_10)
242 chem(its:ite,k,j,p_p25) = chem(its:ite,k,j,p_p25) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_pm_25)
243 chem(its:ite,k,j,p_oc1) = chem(its:ite,k,j,p_oc1) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_oc)
244 chem(its:ite,k,j,p_bc1) = chem(its:ite,k,j,p_bc1) + conv_rho_aer(its:ite)*emis_ant(its:ite,k,j,p_e_bc)
250 ! add aircraft emissions
251 if (config_flags%aircraft_emiss_opt == 1 ) then
253 do k=kts,min(config_flags%kemit_aircraft,kte)
254 conv_rho(its:ite)=4.828e-4/rho_phy(its:ite,k,j)*dtstep/(dz8w(its:ite,k,j)*60.)
255 if( p_no >= param_first_scalar ) then
256 chem(its:ite,k,j,p_no) = chem(its:ite,k,j,p_no) + emis_aircraft(its:ite,k,j,p_eac_no) *conv_rho(its:ite)
258 if( p_co >= param_first_scalar ) then
259 chem(its:ite,k,j,p_co) = chem(its:ite,k,j,p_co) + emis_aircraft(its:ite,k,j,p_eac_co) *conv_rho(its:ite)
261 if( p_so2 >= param_first_scalar ) then
262 chem(its:ite,k,j,p_so2) = chem(its:ite,k,j,p_so2) + emis_aircraft(its:ite,k,j,p_eac_so2)*conv_rho(its:ite)
264 if( p_ch4 >= param_first_scalar ) then
265 chem(its:ite,k,j,p_ch4) = chem(its:ite,k,j,p_ch4) + emis_aircraft(its:ite,k,j,p_eac_ch4)*conv_rho(its:ite)
271 END subroutine add_anthropogenics
274 subroutine add_biogenics(id,dtstep,dz8w,config_flags,rho_phy,chem, &
276 ebio_iso,ebio_oli,ebio_api,ebio_lim,ebio_xyl, &
277 ebio_hc3,ebio_ete,ebio_olt,ebio_ket,ebio_ald, &
278 ebio_hcho,ebio_eth,ebio_ora2,ebio_co,ebio_nr,ebio_no, &
279 ebio_sesq,ebio_mbo, &
280 ids,ide, jds,jde, kds,kde, &
281 ims,ime, jms,jme, kms,kme, &
282 its,ite, jts,jte, kts,kte )
284 USE module_state_description
285 USE module_data_radm2
286 USE module_aerosols_sorgam
288 INTEGER, INTENT(IN ) :: id,ne_area, &
289 ids,ide, jds,jde, kds,kde, &
290 ims,ime, jms,jme, kms,kme, &
291 its,ite, jts,jte, kts,kte
292 REAL, INTENT(IN ) :: &
294 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
295 INTENT(INOUT ) :: chem
296 REAL, DIMENSION( ims:ime, jms:jme,ne_area ), &
299 REAL, DIMENSION( ims:ime, jms:jme ), &
301 ebio_iso,ebio_oli,ebio_api,ebio_lim,ebio_xyl, &
302 ebio_hc3,ebio_ete,ebio_olt,ebio_ket,ebio_ald, &
303 ebio_hcho,ebio_eth,ebio_ora2,ebio_co,ebio_nr,ebio_no, &
306 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , &
310 real :: conv_rho(its:ite)
311 real :: con2ppm(its:ite,jts:jte)
312 !--- deposition and emissions stuff
314 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
316 bioem_select: SELECT CASE(config_flags%bio_emiss_opt)
318 CALL wrf_debug(15,'adding biogenic emissions: Gunther1')
320 conv_rho(its:ite) = dtstep/(dz8w(its:ite,kts,j)*60.)
321 chem(its:ite,kts,j,p_iso)=chem(its:ite,kts,j,p_iso) + e_bio(its:ite,j,p_iso-1)*conv_rho(its:ite)
322 chem(its:ite,kts,j,p_oli)=chem(its:ite,kts,j,p_oli) + e_bio(its:ite,j,p_oli-1)*conv_rho(its:ite)
323 chem(its:ite,kts,j,p_xyl)=chem(its:ite,kts,j,p_xyl) + e_bio(its:ite,j,p_xyl-1)*conv_rho(its:ite)
324 chem(its:ite,kts,j,p_olt)=chem(its:ite,kts,j,p_olt) + e_bio(its:ite,j,p_olt-1)*conv_rho(its:ite)
325 chem(its:ite,kts,j,p_ket)=chem(its:ite,kts,j,p_ket) + e_bio(its:ite,j,p_ket-1)*conv_rho(its:ite)
326 chem(its:ite,kts,j,p_ald)=chem(its:ite,kts,j,p_ald) + e_bio(its:ite,j,p_ald-1)*conv_rho(its:ite)
327 chem(its:ite,kts,j,p_hcho)=chem(its:ite,kts,j,p_hcho) + e_bio(its:ite,j,p_hcho-1)*conv_rho(its:ite)
328 chem(its:ite,kts,j,p_eth)=chem(its:ite,kts,j,p_eth) + e_bio(its:ite,j,p_eth-1)*conv_rho(its:ite)
329 chem(its:ite,kts,j,p_ora2)=chem(its:ite,kts,j,p_ora2) + e_bio(its:ite,j,p_ora2-1)*conv_rho(its:ite)
330 chem(its:ite,kts,j,p_co)=chem(its:ite,kts,j,p_co) + e_bio(its:ite,j,p_co-1)*conv_rho(its:ite)
331 chem(its:ite,kts,j,p_no)=chem(its:ite,kts,j,p_no) + e_bio(its:ite,j,p_no-1)*conv_rho(its:ite)
335 if(p_hc3 >= param_first_scalar) then
336 chem(its:ite,kts,j,p_hc3)=chem(its:ite,kts,j,p_hc3)+ e_bio(its:ite,j,p_hc3-1)*conv_rho(its:ite)
338 if(p_ol2 >= param_first_scalar) then
339 chem(its:ite,kts,j,p_ol2)=chem(its:ite,kts,j,p_ol2)+ e_bio(its:ite,j,p_ol2-1)*conv_rho(its:ite)
342 if(p_c5h8 >= param_first_scalar) then
343 chem(its:ite,kts,j,p_c5h8)=chem(its:ite,kts,j,p_c5h8)+ e_bio(its:ite,j,liso)*conv_rho(its:ite)
345 if(p_oxyl >= param_first_scalar) then
346 chem(its:ite,kts,j,p_oxyl)=chem(its:ite,kts,j,p_oxyl)+ e_bio(its:ite,j,lxyl)*conv_rho(its:ite)
348 if(p_c3h8 >= param_first_scalar) then
349 chem(its:ite,kts,j,p_c3h8)=chem(its:ite,kts,j,p_c3h8)+ e_bio(its:ite,j,lhc3)*conv_rho(its:ite)
351 if(p_ket >= param_first_scalar) then
352 chem(its:ite,kts,j,p_ket)=chem(its:ite,kts,j,p_ket)+ e_bio(its:ite,j,lket)*conv_rho(its:ite)
354 if(p_ch3cho >= param_first_scalar) then
355 chem(its:ite,kts,j,p_ch3cho)=chem(its:ite,kts,j,p_ch3cho)+ e_bio(its:ite,j,lald)*conv_rho(its:ite)
357 if(p_apinene >= param_first_scalar) then
358 chem(its:ite,kts,j,p_apinene)=chem(its:ite,kts,j,p_apinene)+ 0.667*e_bio(its:ite,j,loli)*conv_rho(its:ite)
360 if(p_bpinene >= param_first_scalar) then
361 chem(its:ite,kts,j,p_bpinene)=chem(its:ite,kts,j,p_bpinene)+ .333*(e_bio(its:ite,j,loli) + e_bio(its:ite,j,lolt))*conv_rho(its:ite)
363 if(p_hcooh >= param_first_scalar) then
364 chem(its:ite,kts,j,p_hcooh)=chem(its:ite,kts,j,p_hcooh)+ e_bio(its:ite,j,lora1)*conv_rho(its:ite)
366 if(p_ch3co2h >= param_first_scalar) then
367 chem(its:ite,kts,j,p_ch3co2h)=chem(its:ite,kts,j,p_ch3co2h)+ e_bio(its:ite,j,lora2)*conv_rho(its:ite)
370 !BSINGH - Added for CBMZ
374 if(p_par .gt. 1) then
375 chem(its:ite,kts,j,p_par) = chem(its:ite,kts,j,p_par) &
376 + (e_bio(its:ite,j,p_ald-1)*0.4 &
377 + e_bio(its:ite,j,p_ket-1)*0.9 &
378 + e_bio(its:ite,j,p_olt-1)*1.8 &
379 + e_bio(its:ite,j,p_ora2-1))*conv_rho(its:ite)
381 if(p_hc3 .gt. 1) then
382 chem(its:ite,kts,j,p_par) = chem(its:ite,kts,j,p_par) + e_bio(its:ite,j,p_hc3-1)*4.8*conv_rho(its:ite)
387 CALL wrf_debug(100,'adding biogenic emissions: beis3.1.4')
389 conv_rho(its:ite)=4.828e-4/rho_phy(its:ite,kts,j)*dtstep/(dz8w(its:ite,kts,j)*60.)
390 chem(its:ite,kts,j,p_iso)=chem(its:ite,kts,j,p_iso)+ ebio_iso(its:ite,j)*conv_rho(its:ite)
391 chem(its:ite,kts,j,p_oli)=chem(its:ite,kts,j,p_oli)+ ebio_oli(its:ite,j)*conv_rho(its:ite)
392 chem(its:ite,kts,j,p_xyl)=chem(its:ite,kts,j,p_xyl)+ ebio_xyl(its:ite,j)*conv_rho(its:ite)
393 chem(its:ite,kts,j,p_hc3)=chem(its:ite,kts,j,p_hc3)+ ebio_hc3(its:ite,j)*conv_rho(its:ite)
394 chem(its:ite,kts,j,p_olt)=chem(its:ite,kts,j,p_olt)+ ebio_olt(its:ite,j)*conv_rho(its:ite)
395 chem(its:ite,kts,j,p_ket)=chem(its:ite,kts,j,p_ket)+ ebio_ket(its:ite,j)*conv_rho(its:ite)
396 chem(its:ite,kts,j,p_ald)=chem(its:ite,kts,j,p_ald)+ ebio_ald(its:ite,j)*conv_rho(its:ite)
397 chem(its:ite,kts,j,p_hcho)=chem(its:ite,kts,j,p_hcho)+ ebio_hcho(its:ite,j)*conv_rho(its:ite)
398 chem(its:ite,kts,j,p_eth)=chem(its:ite,kts,j,p_eth)+ ebio_eth(its:ite,j)*conv_rho(its:ite)
399 chem(its:ite,kts,j,p_ora2)=chem(its:ite,kts,j,p_ora2)+ ebio_ora2(its:ite,j)*conv_rho(its:ite)
400 chem(its:ite,kts,j,p_co)=chem(its:ite,kts,j,p_co)+ ebio_co(its:ite,j)*conv_rho(its:ite)
401 chem(its:ite,kts,j,p_no)=chem(its:ite,kts,j,p_no)+ ebio_no(its:ite,j)*conv_rho(its:ite)
405 if(p_ol2 >= param_first_scalar) then
406 chem(its:ite,kts,j,p_ol2)=chem(its:ite,kts,j,p_ol2)+ ebio_ete(its:ite,j)*conv_rho(its:ite)
411 if(p_api >= param_first_scalar) then
412 chem(its:ite,kts,j,p_api)=chem(its:ite,kts,j,p_api)+ ebio_api(its:ite,j)*conv_rho(its:ite)
414 if(p_lim >= param_first_scalar) then
415 chem(its:ite,kts,j,p_lim)=chem(its:ite,kts,j,p_lim)+ ebio_lim(its:ite,j)*conv_rho(its:ite)
417 if(p_ete >= param_first_scalar) then
418 chem(its:ite,kts,j,p_ete)=chem(its:ite,kts,j,p_ete)+ ebio_ete(its:ite,j)*conv_rho(its:ite)
421 ! SESQ and MBO are added, it works for RACM_SOA_VBS_KPP option
423 if(p_sesq >= param_first_scalar) then
424 chem(its:ite,kts,j,p_sesq)=chem(its:ite,kts,j,p_sesq)+ ebio_sesq(its:ite,j)*conv_rho(its:ite)
426 if(p_mbo >= param_first_scalar) then
427 chem(its:ite,kts,j,p_mbo)=chem(its:ite,kts,j,p_mbo)+ ebio_mbo(its:ite,j)*conv_rho(its:ite)
430 CASE (MEGAN2,MEGAN2_CLM)
432 con2ppm(its:ite,j) = dtstep/(dz8w(its:ite,kts,j)*60.)
435 if( any( e_bio(its:ite,jts:jte,n) /= 0. ) ) then
437 chem(its:ite,kts,j,n+1) = chem(its:ite,kts,j,n+1) &
438 + con2ppm(its:ite,j)*e_bio(its:ite,j,n)
445 END SELECT bioem_select
446 END subroutine add_biogenics
449 END MODULE module_emissions_anthropogenics