Line data Source code
1 : module mo_usrrxt
2 :
3 : use shr_kind_mod, only : r8 => shr_kind_r8
4 : use cam_logfile, only : iulog
5 : use ppgrid, only : pver, pcols
6 : use cam_abortutils, only : endrun
7 :
8 : implicit none
9 :
10 : private
11 : public :: usrrxt, usrrxt_inti, usrrxt_hrates
12 :
13 : save
14 :
15 : integer :: usr_O_O2_ndx
16 : integer :: usr_HO2_HO2_ndx
17 : integer :: usr_N2O5_M_ndx
18 : integer :: usr_HNO3_OH_ndx
19 : integer :: usr_HO2NO2_M_ndx
20 : integer :: usr_N2O5_aer_ndx
21 : integer :: usr_NO3_aer_ndx
22 : integer :: usr_NO2_aer_ndx
23 : integer :: usr_CO_OH_ndx
24 : integer :: usr_CO_OH_a_ndx
25 : integer :: usr_CO_OH_b_ndx
26 : integer :: usr_PAN_M_ndx
27 : integer :: usr_CH3COCH3_OH_ndx
28 : integer :: usr_MCO3_NO2_ndx
29 : integer :: usr_MPAN_M_ndx
30 : integer :: usr_XOOH_OH_ndx
31 : integer :: usr_SO2_OH_ndx
32 : integer :: usr_DMS_OH_ndx
33 : integer :: usr_HO2_aer_ndx
34 : integer :: usr_GLYOXAL_aer_ndx
35 :
36 : integer :: tag_NO2_NO3_ndx
37 : integer :: tag_NO2_OH_ndx
38 : integer :: tag_NO2_HO2_ndx
39 : integer :: tag_C2H4_OH_ndx
40 : integer :: tag_C3H6_OH_ndx
41 : integer :: tag_CH3CO3_NO2_ndx
42 :
43 : !lke-TS1
44 : integer :: usr_PBZNIT_M_ndx
45 : integer :: tag_ACBZO2_NO2_ndx
46 : integer :: usr_ISOPNITA_aer_ndx
47 : integer :: usr_ISOPNITB_aer_ndx
48 : integer :: usr_ONITR_aer_ndx
49 : integer :: usr_HONITR_aer_ndx
50 : integer :: usr_TERPNIT_aer_ndx
51 : integer :: usr_NTERPOOH_aer_ndx
52 : integer :: usr_NC4CHO_aer_ndx
53 : integer :: usr_NC4CH2OH_aer_ndx
54 : !TS2
55 : integer :: usr_ISOPZD1O2_ndx
56 : integer :: usr_ISOPZD4O2_ndx
57 : integer :: usr_ISOPFDN_aer_ndx
58 : integer :: usr_ISOPFNP_aer_ndx
59 : integer :: usr_ISOPN2B_aer_ndx
60 : integer :: usr_ISOPN1D_aer_ndx
61 : integer :: usr_ISOPN4D_aer_ndx
62 : integer :: usr_INOOHD_aer_ndx
63 : integer :: usr_INHEB_aer_ndx
64 : integer :: usr_INHED_aer_ndx
65 : integer :: usr_MACRN_aer_ndx
66 : integer :: usr_ISOPHFP_aer_ndx
67 : integer :: usr_IEPOX_aer_ndx
68 : integer :: usr_DHPMPAL_aer_ndx
69 : integer :: usr_ICHE_aer_ndx
70 : integer :: usr_ISOPFNC_aer_ndx
71 : integer :: usr_ISOPFDNC_aer_ndx
72 : integer :: usr_ISOPB1O2_NOn_ndx
73 : integer :: usr_ISOPB1O2_NOa_ndx
74 : integer :: usr_ISOPB4O2_NOn_ndx
75 : integer :: usr_ISOPB4O2_NOa_ndx
76 : integer :: usr_ISOPED1O2_NOn_ndx
77 : integer :: usr_ISOPED1O2_NOa_ndx
78 : integer :: usr_ISOPED4O2_NOn_ndx
79 : integer :: usr_ISOPED4O2_NOa_ndx
80 : integer :: usr_ISOPZD1O2_NOn_ndx
81 : integer :: usr_ISOPZD1O2_NOa_ndx
82 : integer :: usr_ISOPZD4O2_NOn_ndx
83 : integer :: usr_ISOPZD4O2_NOa_ndx
84 : integer :: usr_ISOPNO3_NOn_ndx
85 : integer :: usr_ISOPNO3_NOa_ndx
86 : integer :: usr_MVKO2_NOn_ndx
87 : integer :: usr_MVKO2_NOa_ndx
88 : integer :: usr_MACRO2_NOn_ndx
89 : integer :: usr_MACRO2_NOa_ndx
90 : integer :: usr_IEPOXOO_NOn_ndx
91 : integer :: usr_IEPOXOO_NOa_ndx
92 : integer :: usr_ISOPN1DO2_NOn_ndx
93 : integer :: usr_ISOPN1DO2_NOa_ndx
94 : integer :: usr_ISOPN2BO2_NOn_ndx
95 : integer :: usr_ISOPN2BO2_NOa_ndx
96 : integer :: usr_ISOPN3BO2_NOn_ndx
97 : integer :: usr_ISOPN3BO2_NOa_ndx
98 : integer :: usr_ISOPN4DO2_NOn_ndx
99 : integer :: usr_ISOPN4DO2_NOa_ndx
100 : integer :: usr_ISOPNBNO3O2_NOn_ndx
101 : integer :: usr_ISOPNBNO3O2_NOa_ndx
102 : integer :: usr_ISOPNOOHBO2_NOn_ndx
103 : integer :: usr_ISOPNOOHBO2_NOa_ndx
104 : integer :: usr_ISOPNOOHDO2_NOn_ndx
105 : integer :: usr_ISOPNOOHDO2_NOa_ndx
106 : integer :: usr_NC4CHOO2_NOn_ndx
107 : integer :: usr_NC4CHOO2_NOa_ndx
108 : integer :: tag_MCO3_NO2_ndx
109 : integer :: tag_TERPACO3_NO2_ndx
110 : integer :: usr_TERPAPAN_M_ndx
111 : integer :: tag_TERPA2CO3_NO2_ndx
112 : integer :: usr_TERPA2PAN_M_ndx
113 : integer :: tag_TERPA3CO3_NO2_ndx
114 : integer :: usr_TERPA3PAN_M_ndx
115 : integer :: usr_TERPNT_aer_ndx
116 : integer :: usr_TERPNT1_aer_ndx
117 : integer :: usr_TERPNPT_aer_ndx
118 : integer :: usr_TERPNPT1_aer_ndx
119 : integer :: usr_TERPFDN_aer_ndx
120 : integer :: usr_SQTN_aer_ndx
121 : integer :: usr_TERPHFN_aer_ndx
122 : integer :: usr_TERPDHDP_aer_ndx
123 : integer :: usr_TERPACID_aer_ndx
124 : !
125 : integer :: usr_OA_O2_NDX
126 : integer :: usr_XNO2NO3_M_ndx
127 : integer :: usr_NO2XNO3_M_ndx
128 : integer :: usr_XHNO3_OH_ndx
129 : integer :: usr_XHO2NO2_M_ndx
130 : integer :: usr_XNO2NO3_aer_ndx
131 : integer :: usr_NO2XNO3_aer_ndx
132 : integer :: usr_XNO3_aer_ndx
133 : integer :: usr_XNO2_aer_ndx
134 : integer :: usr_XPAN_M_ndx
135 : integer :: usr_XMPAN_M_ndx
136 : integer :: usr_MCO3_XNO2_ndx
137 :
138 : integer :: usr_C2O3_NO2_ndx
139 : integer :: usr_C2H4_OH_ndx
140 : integer :: usr_XO2N_HO2_ndx
141 : integer :: usr_C2O3_XNO2_ndx
142 :
143 : integer :: tag_XO2N_NO_ndx
144 : integer :: tag_XO2_HO2_ndx
145 : integer :: tag_XO2_NO_ndx
146 :
147 : integer :: usr_O_O_ndx
148 : integer :: usr_CL2O2_M_ndx
149 : integer :: usr_SO3_H2O_ndx
150 : integer :: tag_CLO_CLO_M_ndx
151 :
152 : integer :: ion1_ndx, ion2_ndx, ion3_ndx, ion11_ndx
153 : integer :: elec1_ndx, elec2_ndx, elec3_ndx
154 : integer :: elec4_ndx, elec5_ndx, elec6_ndx
155 : integer :: het1_ndx
156 :
157 : integer, parameter :: nean = 3
158 : integer :: ean_ndx(nean)
159 : integer, parameter :: nrpe = 5
160 : integer :: rpe_ndx(nrpe)
161 : integer, parameter :: npir = 16
162 : integer :: pir_ndx(npir)
163 : integer, parameter :: nedn = 2
164 : integer :: edn_ndx(nedn)
165 : integer, parameter :: nnir = 13
166 : integer :: nir_ndx(nnir)
167 : integer, parameter :: niira = 112
168 : integer :: iira_ndx(niira)
169 : integer, parameter :: niirb = 14
170 : integer :: iirb_ndx(niirb)
171 :
172 : integer :: usr_clm_h2o_m_ndx, usr_clm_hcl_m_ndx
173 : integer :: usr_oh_co_ndx, het_no2_h2o_ndx, usr_oh_dms_ndx, aq_so2_h2o2_ndx, aq_so2_o3_ndx
174 :
175 : integer :: h2o_ndx
176 : !
177 : ! jfl
178 : !
179 : integer, parameter :: num_strat_tau = 22
180 : integer :: usr_strat_tau_ndx(num_strat_tau)
181 : !
182 : !lke++
183 : integer :: usr_COhc_OH_ndx
184 : integer :: usr_COme_OH_ndx
185 : integer :: usr_CO01_OH_ndx
186 : integer :: usr_CO02_OH_ndx
187 : integer :: usr_CO03_OH_ndx
188 : integer :: usr_CO04_OH_ndx
189 : integer :: usr_CO05_OH_ndx
190 : integer :: usr_CO06_OH_ndx
191 : integer :: usr_CO07_OH_ndx
192 : integer :: usr_CO08_OH_ndx
193 : integer :: usr_CO09_OH_ndx
194 : integer :: usr_CO10_OH_ndx
195 : integer :: usr_CO11_OH_ndx
196 : integer :: usr_CO12_OH_ndx
197 : integer :: usr_CO13_OH_ndx
198 : integer :: usr_CO14_OH_ndx
199 : integer :: usr_CO15_OH_ndx
200 : integer :: usr_CO16_OH_ndx
201 : integer :: usr_CO17_OH_ndx
202 : integer :: usr_CO18_OH_ndx
203 : integer :: usr_CO19_OH_ndx
204 : integer :: usr_CO20_OH_ndx
205 : integer :: usr_CO21_OH_ndx
206 : integer :: usr_CO22_OH_ndx
207 : integer :: usr_CO23_OH_ndx
208 : integer :: usr_CO24_OH_ndx
209 : integer :: usr_CO25_OH_ndx
210 : integer :: usr_CO26_OH_ndx
211 : integer :: usr_CO27_OH_ndx
212 : integer :: usr_CO28_OH_ndx
213 : integer :: usr_CO29_OH_ndx
214 : integer :: usr_CO30_OH_ndx
215 : integer :: usr_CO31_OH_ndx
216 : integer :: usr_CO32_OH_ndx
217 : integer :: usr_CO33_OH_ndx
218 : integer :: usr_CO34_OH_ndx
219 : integer :: usr_CO35_OH_ndx
220 : integer :: usr_CO36_OH_ndx
221 : integer :: usr_CO37_OH_ndx
222 : integer :: usr_CO38_OH_ndx
223 : integer :: usr_CO39_OH_ndx
224 : integer :: usr_CO40_OH_ndx
225 : integer :: usr_CO41_OH_ndx
226 : integer :: usr_CO42_OH_ndx
227 : !lke--
228 :
229 : real(r8), parameter :: t0 = 300._r8 ! K
230 : real(r8), parameter :: trlim2 = 17._r8/3._r8 ! K
231 : real(r8), parameter :: trlim3 = 15._r8/3._r8 ! K
232 :
233 : logical :: has_ion_rxts, has_d_chem
234 :
235 : contains
236 :
237 1536 : subroutine usrrxt_inti
238 : !-----------------------------------------------------------------
239 : ! ... intialize the user reaction constants module
240 : !-----------------------------------------------------------------
241 :
242 : use mo_chem_utls, only : get_rxt_ndx, get_spc_ndx
243 : use spmd_utils, only : masterproc
244 :
245 : implicit none
246 :
247 : character(len=4) :: xchar
248 : character(len=32) :: rxtname
249 : integer :: i
250 :
251 : !
252 : ! full tropospheric chemistry
253 : !
254 1536 : usr_O_O2_ndx = get_rxt_ndx( 'usr_O_O2' )
255 1536 : usr_HO2_HO2_ndx = get_rxt_ndx( 'usr_HO2_HO2' )
256 1536 : usr_N2O5_M_ndx = get_rxt_ndx( 'usr_N2O5_M' )
257 1536 : usr_HNO3_OH_ndx = get_rxt_ndx( 'usr_HNO3_OH' )
258 1536 : usr_HO2NO2_M_ndx = get_rxt_ndx( 'usr_HO2NO2_M' )
259 1536 : usr_N2O5_aer_ndx = get_rxt_ndx( 'usr_N2O5_aer' )
260 1536 : usr_NO3_aer_ndx = get_rxt_ndx( 'usr_NO3_aer' )
261 1536 : usr_NO2_aer_ndx = get_rxt_ndx( 'usr_NO2_aer' )
262 1536 : usr_CO_OH_ndx = get_rxt_ndx( 'usr_CO_OH' )
263 1536 : usr_CO_OH_a_ndx = get_rxt_ndx( 'usr_CO_OH_a' )
264 1536 : usr_CO_OH_b_ndx = get_rxt_ndx( 'usr_CO_OH_b' )
265 1536 : usr_PAN_M_ndx = get_rxt_ndx( 'usr_PAN_M' )
266 1536 : usr_CH3COCH3_OH_ndx = get_rxt_ndx( 'usr_CH3COCH3_OH' )
267 1536 : usr_MCO3_NO2_ndx = get_rxt_ndx( 'usr_MCO3_NO2' )
268 1536 : tag_MCO3_NO2_ndx = get_rxt_ndx( 'tag_MCO3_NO2' )
269 1536 : if( tag_MCO3_NO2_ndx<0 .and. usr_MCO3_NO2_ndx>0 ) then
270 0 : tag_MCO3_NO2_ndx = usr_MCO3_NO2_ndx
271 : endif
272 :
273 1536 : usr_MPAN_M_ndx = get_rxt_ndx( 'usr_MPAN_M' )
274 1536 : usr_XOOH_OH_ndx = get_rxt_ndx( 'usr_XOOH_OH' )
275 1536 : usr_SO2_OH_ndx = get_rxt_ndx( 'usr_SO2_OH' )
276 1536 : usr_DMS_OH_ndx = get_rxt_ndx( 'usr_DMS_OH' )
277 1536 : usr_HO2_aer_ndx = get_rxt_ndx( 'usr_HO2_aer' )
278 1536 : usr_GLYOXAL_aer_ndx = get_rxt_ndx( 'usr_GLYOXAL_aer' )
279 : !
280 1536 : tag_NO2_NO3_ndx = get_rxt_ndx( 'tag_NO2_NO3' )
281 1536 : tag_NO2_OH_ndx = get_rxt_ndx( 'tag_NO2_OH' )
282 1536 : tag_NO2_HO2_ndx = get_rxt_ndx( 'tag_NO2_HO2' )
283 1536 : tag_C2H4_OH_ndx = get_rxt_ndx( 'tag_C2H4_OH' )
284 1536 : tag_C3H6_OH_ndx = get_rxt_ndx( 'tag_C3H6_OH' )
285 1536 : tag_CH3CO3_NO2_ndx = get_rxt_ndx( 'tag_CH3CO3_NO2' )
286 : !lke-TS1
287 1536 : usr_PBZNIT_M_ndx = get_rxt_ndx( 'usr_PBZNIT_M' )
288 1536 : tag_ACBZO2_NO2_ndx = get_rxt_ndx( 'tag_ACBZO2_NO2' )
289 1536 : usr_ISOPNITA_aer_ndx = get_rxt_ndx( 'usr_ISOPNITA_aer' )
290 1536 : usr_ISOPNITB_aer_ndx = get_rxt_ndx( 'usr_ISOPNITB_aer' )
291 1536 : usr_ONITR_aer_ndx = get_rxt_ndx( 'usr_ONITR_aer' )
292 1536 : usr_HONITR_aer_ndx = get_rxt_ndx( 'usr_HONITR_aer' )
293 1536 : usr_TERPNIT_aer_ndx = get_rxt_ndx( 'usr_TERPNIT_aer' )
294 1536 : usr_NTERPOOH_aer_ndx = get_rxt_ndx( 'usr_NTERPOOH_aer' )
295 1536 : usr_NC4CHO_aer_ndx = get_rxt_ndx( 'usr_NC4CHO_aer' )
296 1536 : usr_NC4CH2OH_aer_ndx = get_rxt_ndx( 'usr_NC4CH2OH_aer' )
297 : !TS2
298 1536 : usr_ISOPZD1O2_ndx = get_rxt_ndx( 'usr_ISOPZD1O2' )
299 1536 : usr_ISOPZD4O2_ndx = get_rxt_ndx( 'usr_ISOPZD4O2' )
300 1536 : usr_ISOPFDN_aer_ndx = get_rxt_ndx( 'usr_ISOPFDN_aer' )
301 1536 : usr_ISOPFNP_aer_ndx = get_rxt_ndx( 'usr_ISOPFNP_aer' )
302 1536 : usr_ISOPN2B_aer_ndx = get_rxt_ndx( 'usr_ISOPN2B_aer' )
303 1536 : usr_ISOPN1D_aer_ndx = get_rxt_ndx( 'usr_ISOPN1D_aer' )
304 1536 : usr_ISOPN4D_aer_ndx = get_rxt_ndx( 'usr_ISOPN4D_aer' )
305 1536 : usr_INOOHD_aer_ndx = get_rxt_ndx( 'usr_INOOHD_aer' )
306 1536 : usr_INHEB_aer_ndx = get_rxt_ndx( 'usr_INHEB_aer' )
307 1536 : usr_INHED_aer_ndx = get_rxt_ndx( 'usr_INHED_aer' )
308 1536 : usr_MACRN_aer_ndx = get_rxt_ndx( 'usr_MACRN_aer' )
309 1536 : usr_ISOPHFP_aer_ndx = get_rxt_ndx( 'usr_ISOPHFP_aer' )
310 1536 : usr_IEPOX_aer_ndx = get_rxt_ndx( 'usr_IEPOX_aer' )
311 1536 : usr_DHPMPAL_aer_ndx = get_rxt_ndx( 'usr_DHPMPAL_aer' )
312 1536 : usr_ICHE_aer_ndx = get_rxt_ndx( 'usr_ICHE_aer' )
313 1536 : usr_ISOPFNC_aer_ndx = get_rxt_ndx( 'usr_ISOPFNC_aer' )
314 1536 : usr_ISOPFDNC_aer_ndx = get_rxt_ndx( 'usr_ISOPFDNC_aer' )
315 1536 : usr_ISOPB1O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPB1O2_NOn' )
316 1536 : usr_ISOPB1O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPB1O2_NOa' )
317 1536 : usr_ISOPB4O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPB4O2_NOn' )
318 1536 : usr_ISOPB4O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPB4O2_NOa' )
319 1536 : usr_ISOPED1O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPED1O2_NOn' )
320 1536 : usr_ISOPED1O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPED1O2_NOa' )
321 1536 : usr_ISOPED4O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPED4O2_NOn' )
322 1536 : usr_ISOPED4O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPED4O2_NOa' )
323 1536 : usr_ISOPZD1O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPZD1O2_NOn' )
324 1536 : usr_ISOPZD1O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPZD1O2_NOa' )
325 1536 : usr_ISOPZD4O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPZD4O2_NOn' )
326 1536 : usr_ISOPZD4O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPZD4O2_NOa' )
327 1536 : usr_ISOPNO3_NOn_ndx = get_rxt_ndx( 'usr_ISOPNO3_NOn' )
328 1536 : usr_ISOPNO3_NOa_ndx = get_rxt_ndx( 'usr_ISOPNO3_NOa' )
329 1536 : usr_MVKO2_NOn_ndx = get_rxt_ndx( 'usr_MVKO2_NOn' )
330 1536 : usr_MVKO2_NOa_ndx = get_rxt_ndx( 'usr_MVKO2_NOa' )
331 1536 : usr_MACRO2_NOn_ndx = get_rxt_ndx( 'usr_MACRO2_NOn' )
332 1536 : usr_MACRO2_NOa_ndx = get_rxt_ndx( 'usr_MACRO2_NOa' )
333 1536 : usr_IEPOXOO_NOn_ndx = get_rxt_ndx( 'usr_IEPOXOO_NOn' )
334 1536 : usr_IEPOXOO_NOa_ndx = get_rxt_ndx( 'usr_IEPOXOO_NOa' )
335 1536 : usr_ISOPN1DO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPN1DO2_NOn' )
336 1536 : usr_ISOPN1DO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPN1DO2_NOa' )
337 1536 : usr_ISOPN2BO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPN2BO2_NOn' )
338 1536 : usr_ISOPN2BO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPN2BO2_NOa' )
339 1536 : usr_ISOPN3BO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPN3BO2_NOn' )
340 1536 : usr_ISOPN3BO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPN3BO2_NOa' )
341 1536 : usr_ISOPN4DO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPN4DO2_NOn' )
342 1536 : usr_ISOPN4DO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPN4DO2_NOa' )
343 1536 : usr_ISOPNBNO3O2_NOn_ndx = get_rxt_ndx( 'usr_ISOPNBNO3O2_NOn' )
344 1536 : usr_ISOPNBNO3O2_NOa_ndx = get_rxt_ndx( 'usr_ISOPNBNO3O2_NOa' )
345 1536 : usr_ISOPNOOHBO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPNOOHBO2_NOn' )
346 1536 : usr_ISOPNOOHBO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPNOOHBO2_NOa' )
347 1536 : usr_ISOPNOOHDO2_NOn_ndx = get_rxt_ndx( 'usr_ISOPNOOHDO2_NOn' )
348 1536 : usr_ISOPNOOHDO2_NOa_ndx = get_rxt_ndx( 'usr_ISOPNOOHDO2_NOa' )
349 1536 : usr_NC4CHOO2_NOn_ndx = get_rxt_ndx( 'usr_NC4CHOO2_NOn' )
350 1536 : usr_NC4CHOO2_NOa_ndx = get_rxt_ndx( 'usr_NC4CHOO2_NOa' )
351 1536 : tag_TERPACO3_NO2_ndx = get_rxt_ndx( 'tag_TERPACO3_NO2' )
352 1536 : usr_TERPAPAN_M_ndx = get_rxt_ndx( 'usr_TERPAPAN_M' )
353 1536 : tag_TERPA2CO3_NO2_ndx = get_rxt_ndx( 'tag_TERPA2CO3_NO2' )
354 1536 : usr_TERPA2PAN_M_ndx = get_rxt_ndx( 'usr_TERPA2PAN_M' )
355 1536 : tag_TERPA3CO3_NO2_ndx = get_rxt_ndx( 'tag_TERPA3CO3_NO2' )
356 1536 : usr_TERPA3PAN_M_ndx = get_rxt_ndx( 'usr_TERPA3PAN_M' )
357 1536 : usr_TERPNT_aer_ndx = get_rxt_ndx( 'usr_TERPNT_aer' )
358 1536 : usr_TERPNT1_aer_ndx = get_rxt_ndx( 'usr_TERPNT1_aer' )
359 1536 : usr_TERPNPT_aer_ndx = get_rxt_ndx( 'usr_TERPNPT_aer' )
360 1536 : usr_TERPNPT1_aer_ndx = get_rxt_ndx( 'usr_TERPNPT1_aer' )
361 1536 : usr_TERPFDN_aer_ndx = get_rxt_ndx( 'usr_TERPFDN_aer' )
362 1536 : usr_SQTN_aer_ndx = get_rxt_ndx( 'usr_SQTN_aer' )
363 1536 : usr_TERPHFN_aer_ndx = get_rxt_ndx( 'usr_TERPHFN_aer' )
364 1536 : usr_TERPDHDP_aer_ndx = get_rxt_ndx( 'usr_TERPDHDP_aer' )
365 1536 : usr_TERPACID_aer_ndx = get_rxt_ndx( 'usr_TERPACID_aer' )
366 : !
367 : ! additional reactions for O3A/XNO
368 : !
369 1536 : usr_OA_O2_ndx = get_rxt_ndx( 'usr_OA_O2' )
370 1536 : usr_XNO2NO3_M_ndx = get_rxt_ndx( 'usr_XNO2NO3_M' )
371 1536 : usr_NO2XNO3_M_ndx = get_rxt_ndx( 'usr_NO2XNO3_M' )
372 1536 : usr_XNO2NO3_aer_ndx = get_rxt_ndx( 'usr_XNO2NO3_aer' )
373 1536 : usr_NO2XNO3_aer_ndx = get_rxt_ndx( 'usr_NO2XNO3_aer' )
374 1536 : usr_XHNO3_OH_ndx = get_rxt_ndx( 'usr_XHNO3_OH' )
375 1536 : usr_XNO3_aer_ndx = get_rxt_ndx( 'usr_XNO3_aer' )
376 1536 : usr_XNO2_aer_ndx = get_rxt_ndx( 'usr_XNO2_aer' )
377 1536 : usr_MCO3_XNO2_ndx = get_rxt_ndx( 'usr_MCO3_XNO2' )
378 1536 : usr_XPAN_M_ndx = get_rxt_ndx( 'usr_XPAN_M' )
379 1536 : usr_XMPAN_M_ndx = get_rxt_ndx( 'usr_XMPAN_M' )
380 1536 : usr_XHO2NO2_M_ndx = get_rxt_ndx( 'usr_XHO2NO2_M' )
381 : !
382 : ! reduced hydrocarbon chemistry
383 : !
384 1536 : usr_C2O3_NO2_ndx = get_rxt_ndx( 'usr_C2O3_NO2' )
385 1536 : usr_C2H4_OH_ndx = get_rxt_ndx( 'usr_C2H4_OH' )
386 1536 : usr_XO2N_HO2_ndx = get_rxt_ndx( 'usr_XO2N_HO2' )
387 1536 : usr_C2O3_XNO2_ndx = get_rxt_ndx( 'usr_C2O3_XNO2' )
388 : !
389 1536 : tag_XO2N_NO_ndx = get_rxt_ndx( 'tag_XO2N_NO' )
390 1536 : tag_XO2_HO2_ndx = get_rxt_ndx( 'tag_XO2_HO2' )
391 1536 : tag_XO2_NO_ndx = get_rxt_ndx( 'tag_XO2_NO' )
392 : !
393 : ! stratospheric chemistry
394 : !
395 1536 : usr_O_O_ndx = get_rxt_ndx( 'usr_O_O' )
396 1536 : usr_CL2O2_M_ndx = get_rxt_ndx( 'usr_CL2O2_M' )
397 1536 : usr_SO3_H2O_ndx = get_rxt_ndx( 'usr_SO3_H2O' )
398 : !
399 1536 : tag_CLO_CLO_M_ndx = get_rxt_ndx( 'tag_CLO_CLO_M' )
400 1536 : if (tag_CLO_CLO_M_ndx<0) then ! for backwards compatibility
401 1536 : tag_CLO_CLO_M_ndx = get_rxt_ndx( 'tag_CLO_CLO' )
402 : endif
403 : !
404 : ! reactions to remove BAM aerosols in the stratosphere
405 : !
406 1536 : usr_strat_tau_ndx( 1) = get_rxt_ndx( 'usr_CB1_strat_tau' )
407 1536 : usr_strat_tau_ndx( 2) = get_rxt_ndx( 'usr_CB2_strat_tau' )
408 1536 : usr_strat_tau_ndx( 3) = get_rxt_ndx( 'usr_OC1_strat_tau' )
409 1536 : usr_strat_tau_ndx( 4) = get_rxt_ndx( 'usr_OC2_strat_tau' )
410 1536 : usr_strat_tau_ndx( 5) = get_rxt_ndx( 'usr_SO4_strat_tau' )
411 1536 : usr_strat_tau_ndx( 6) = get_rxt_ndx( 'usr_SOA_strat_tau' )
412 1536 : usr_strat_tau_ndx( 7) = get_rxt_ndx( 'usr_NH4_strat_tau' )
413 1536 : usr_strat_tau_ndx( 8) = get_rxt_ndx( 'usr_NH4NO3_strat_tau' )
414 1536 : usr_strat_tau_ndx( 9) = get_rxt_ndx( 'usr_SSLT01_strat_tau' )
415 1536 : usr_strat_tau_ndx(10) = get_rxt_ndx( 'usr_SSLT02_strat_tau' )
416 1536 : usr_strat_tau_ndx(11) = get_rxt_ndx( 'usr_SSLT03_strat_tau' )
417 1536 : usr_strat_tau_ndx(12) = get_rxt_ndx( 'usr_SSLT04_strat_tau' )
418 1536 : usr_strat_tau_ndx(13) = get_rxt_ndx( 'usr_DST01_strat_tau' )
419 1536 : usr_strat_tau_ndx(14) = get_rxt_ndx( 'usr_DST02_strat_tau' )
420 1536 : usr_strat_tau_ndx(15) = get_rxt_ndx( 'usr_DST03_strat_tau' )
421 1536 : usr_strat_tau_ndx(16) = get_rxt_ndx( 'usr_DST04_strat_tau' )
422 1536 : usr_strat_tau_ndx(17) = get_rxt_ndx( 'usr_SO2t_strat_tau' )
423 1536 : usr_strat_tau_ndx(18) = get_rxt_ndx( 'usr_SOAI_strat_tau' )
424 1536 : usr_strat_tau_ndx(19) = get_rxt_ndx( 'usr_SOAM_strat_tau' )
425 1536 : usr_strat_tau_ndx(20) = get_rxt_ndx( 'usr_SOAB_strat_tau' )
426 1536 : usr_strat_tau_ndx(21) = get_rxt_ndx( 'usr_SOAT_strat_tau' )
427 1536 : usr_strat_tau_ndx(22) = get_rxt_ndx( 'usr_SOAX_strat_tau' )
428 : !
429 : ! stratospheric aerosol chemistry
430 : !
431 1536 : het1_ndx = get_rxt_ndx( 'het1' )
432 : !
433 : ! ion chemistry
434 : !
435 1536 : ion1_ndx = get_rxt_ndx( 'ion_Op_O2' )
436 1536 : ion2_ndx = get_rxt_ndx( 'ion_Op_N2' )
437 1536 : ion3_ndx = get_rxt_ndx( 'ion_N2p_Oa' )
438 1536 : ion11_ndx = get_rxt_ndx( 'ion_N2p_Ob' )
439 :
440 1536 : elec1_ndx = get_rxt_ndx( 'elec1' )
441 1536 : elec2_ndx = get_rxt_ndx( 'elec2' )
442 1536 : elec3_ndx = get_rxt_ndx( 'elec3' )
443 :
444 6144 : do i = 1,nean
445 4608 : write (xchar,'(i4)') i
446 4608 : rxtname = 'ean'//trim(adjustl(xchar))
447 6144 : ean_ndx(i) = get_rxt_ndx(trim(rxtname))
448 : enddo
449 :
450 9216 : do i = 1,nrpe
451 7680 : write (xchar,'(i4)') i
452 7680 : rxtname = 'rpe'//trim(adjustl(xchar))
453 9216 : rpe_ndx(i) = get_rxt_ndx(trim(rxtname))
454 : enddo
455 :
456 26112 : do i = 1,npir
457 24576 : write (xchar,'(i4)') i
458 24576 : rxtname = 'pir'//trim(adjustl(xchar))
459 26112 : pir_ndx(i) = get_rxt_ndx(trim(rxtname))
460 : enddo
461 :
462 4608 : do i = 1,nedn
463 3072 : write (xchar,'(i4)') i
464 3072 : rxtname = 'edn'//trim(adjustl(xchar))
465 4608 : edn_ndx(i) = get_rxt_ndx(trim(rxtname))
466 : enddo
467 :
468 21504 : do i = 1,nnir
469 19968 : write (xchar,'(i4)') i
470 19968 : rxtname = 'nir'//trim(adjustl(xchar))
471 21504 : nir_ndx(i) = get_rxt_ndx(trim(rxtname))
472 : enddo
473 :
474 173568 : do i = 1,niira
475 172032 : write (xchar,'(i4)') i
476 172032 : rxtname = 'iira'//trim(adjustl(xchar))
477 173568 : iira_ndx(i) = get_rxt_ndx(trim(rxtname))
478 : enddo
479 :
480 23040 : do i = 1,niirb
481 21504 : write (xchar,'(i4)') i
482 21504 : rxtname = 'iirb'//trim(adjustl(xchar))
483 23040 : iirb_ndx(i) = get_rxt_ndx(trim(rxtname))
484 : enddo
485 :
486 1536 : usr_clm_h2o_m_ndx = get_rxt_ndx( 'usr_CLm_H2O_M' )
487 1536 : usr_clm_hcl_m_ndx = get_rxt_ndx( 'usr_CLm_HCL_M' )
488 :
489 1536 : elec4_ndx = get_rxt_ndx( 'Op2P_ea' )
490 1536 : elec5_ndx = get_rxt_ndx( 'Op2P_eb' )
491 1536 : elec6_ndx = get_rxt_ndx( 'Op2D_e' )
492 :
493 : has_ion_rxts = ion1_ndx>0 .and. ion2_ndx>0 .and. ion3_ndx>0 .and. elec1_ndx>0 &
494 1536 : .and. elec2_ndx>0 .and. elec3_ndx>0
495 :
496 : has_d_chem = &
497 : all(ean_ndx>0) .and. &
498 : all(rpe_ndx>0) .and. &
499 : all(pir_ndx>0) .and. &
500 : all(edn_ndx>0) .and. &
501 : all(nir_ndx>0) .and. &
502 : all(iira_ndx>0) .and. &
503 : all(iirb_ndx>0) .and. &
504 10752 : usr_clm_h2o_m_ndx>0 .and. usr_clm_hcl_m_ndx>0
505 :
506 1536 : h2o_ndx = get_spc_ndx( 'H2O' )
507 :
508 : !
509 : ! llnl super fast
510 : !
511 1536 : usr_oh_co_ndx = get_rxt_ndx( 'usr_oh_co' )
512 1536 : het_no2_h2o_ndx = get_rxt_ndx( 'het_no2_h2o' )
513 1536 : usr_oh_dms_ndx = get_rxt_ndx( 'usr_oh_dms' )
514 1536 : aq_so2_h2o2_ndx = get_rxt_ndx( 'aq_so2_h2o2' )
515 1536 : aq_so2_o3_ndx = get_rxt_ndx( 'aq_so2_o3' )
516 :
517 : !lke++
518 : ! CO tags
519 : !
520 1536 : usr_COhc_OH_ndx = get_rxt_ndx( 'usr_COhc_OH' )
521 1536 : usr_COme_OH_ndx = get_rxt_ndx( 'usr_COme_OH' )
522 1536 : usr_CO01_OH_ndx = get_rxt_ndx( 'usr_CO01_OH' )
523 1536 : usr_CO02_OH_ndx = get_rxt_ndx( 'usr_CO02_OH' )
524 1536 : usr_CO03_OH_ndx = get_rxt_ndx( 'usr_CO03_OH' )
525 1536 : usr_CO04_OH_ndx = get_rxt_ndx( 'usr_CO04_OH' )
526 1536 : usr_CO05_OH_ndx = get_rxt_ndx( 'usr_CO05_OH' )
527 1536 : usr_CO06_OH_ndx = get_rxt_ndx( 'usr_CO06_OH' )
528 1536 : usr_CO07_OH_ndx = get_rxt_ndx( 'usr_CO07_OH' )
529 1536 : usr_CO08_OH_ndx = get_rxt_ndx( 'usr_CO08_OH' )
530 1536 : usr_CO09_OH_ndx = get_rxt_ndx( 'usr_CO09_OH' )
531 1536 : usr_CO10_OH_ndx = get_rxt_ndx( 'usr_CO10_OH' )
532 1536 : usr_CO11_OH_ndx = get_rxt_ndx( 'usr_CO11_OH' )
533 1536 : usr_CO12_OH_ndx = get_rxt_ndx( 'usr_CO12_OH' )
534 1536 : usr_CO13_OH_ndx = get_rxt_ndx( 'usr_CO13_OH' )
535 1536 : usr_CO14_OH_ndx = get_rxt_ndx( 'usr_CO14_OH' )
536 1536 : usr_CO15_OH_ndx = get_rxt_ndx( 'usr_CO15_OH' )
537 1536 : usr_CO16_OH_ndx = get_rxt_ndx( 'usr_CO16_OH' )
538 1536 : usr_CO17_OH_ndx = get_rxt_ndx( 'usr_CO17_OH' )
539 1536 : usr_CO18_OH_ndx = get_rxt_ndx( 'usr_CO18_OH' )
540 1536 : usr_CO19_OH_ndx = get_rxt_ndx( 'usr_CO19_OH' )
541 1536 : usr_CO20_OH_ndx = get_rxt_ndx( 'usr_CO20_OH' )
542 1536 : usr_CO21_OH_ndx = get_rxt_ndx( 'usr_CO21_OH' )
543 1536 : usr_CO22_OH_ndx = get_rxt_ndx( 'usr_CO22_OH' )
544 1536 : usr_CO23_OH_ndx = get_rxt_ndx( 'usr_CO23_OH' )
545 1536 : usr_CO24_OH_ndx = get_rxt_ndx( 'usr_CO24_OH' )
546 1536 : usr_CO25_OH_ndx = get_rxt_ndx( 'usr_CO25_OH' )
547 1536 : usr_CO26_OH_ndx = get_rxt_ndx( 'usr_CO26_OH' )
548 1536 : usr_CO27_OH_ndx = get_rxt_ndx( 'usr_CO27_OH' )
549 1536 : usr_CO28_OH_ndx = get_rxt_ndx( 'usr_CO28_OH' )
550 1536 : usr_CO29_OH_ndx = get_rxt_ndx( 'usr_CO29_OH' )
551 1536 : usr_CO30_OH_ndx = get_rxt_ndx( 'usr_CO30_OH' )
552 1536 : usr_CO31_OH_ndx = get_rxt_ndx( 'usr_CO31_OH' )
553 1536 : usr_CO32_OH_ndx = get_rxt_ndx( 'usr_CO32_OH' )
554 1536 : usr_CO33_OH_ndx = get_rxt_ndx( 'usr_CO33_OH' )
555 1536 : usr_CO34_OH_ndx = get_rxt_ndx( 'usr_CO34_OH' )
556 1536 : usr_CO35_OH_ndx = get_rxt_ndx( 'usr_CO35_OH' )
557 1536 : usr_CO36_OH_ndx = get_rxt_ndx( 'usr_CO36_OH' )
558 1536 : usr_CO37_OH_ndx = get_rxt_ndx( 'usr_CO37_OH' )
559 1536 : usr_CO38_OH_ndx = get_rxt_ndx( 'usr_CO38_OH' )
560 1536 : usr_CO39_OH_ndx = get_rxt_ndx( 'usr_CO39_OH' )
561 1536 : usr_CO40_OH_ndx = get_rxt_ndx( 'usr_CO40_OH' )
562 1536 : usr_CO41_OH_ndx = get_rxt_ndx( 'usr_CO41_OH' )
563 1536 : usr_CO42_OH_ndx = get_rxt_ndx( 'usr_CO42_OH' )
564 : !lke--
565 :
566 1536 : if (masterproc) then
567 2 : write(iulog,*) ' '
568 2 : write(iulog,*) 'usrrxt_inti: diagnostics '
569 2 : write(iulog,'(10i5)') usr_O_O2_ndx,usr_HO2_HO2_ndx,tag_NO2_NO3_ndx,usr_N2O5_M_ndx,tag_NO2_OH_ndx,usr_HNO3_OH_ndx &
570 2 : ,tag_NO2_HO2_ndx,usr_HO2NO2_M_ndx,usr_N2O5_aer_ndx,usr_NO3_aer_ndx,usr_NO2_aer_ndx &
571 2 : ,usr_CO_OH_b_ndx,tag_C2H4_OH_ndx,tag_C3H6_OH_ndx,tag_CH3CO3_NO2_ndx,usr_PAN_M_ndx,usr_CH3COCH3_OH_ndx &
572 2 : ,usr_MCO3_NO2_ndx,usr_MPAN_M_ndx,usr_XOOH_OH_ndx,usr_SO2_OH_ndx,usr_DMS_OH_ndx,usr_HO2_aer_ndx &
573 2 : ,usr_GLYOXAL_aer_ndx,usr_ISOPNITA_aer_ndx,usr_ISOPNITB_aer_ndx,usr_ONITR_aer_ndx,usr_HONITR_aer_ndx &
574 2 : ,usr_TERPNIT_aer_ndx,usr_NTERPOOH_aer_ndx,usr_NC4CHO_aer_ndx,usr_NC4CH2OH_aer_ndx,usr_ISOPZD1O2_ndx &
575 2 : ,usr_ISOPZD4O2_ndx,usr_ISOPFDN_aer_ndx,usr_ISOPFNP_aer_ndx,usr_ISOPN2B_aer_ndx,usr_ISOPN1D_aer_ndx &
576 2 : ,usr_ISOPN4D_aer_ndx,usr_INOOHD_aer_ndx,usr_INHEB_aer_ndx,usr_INHED_aer_ndx,usr_MACRN_aer_ndx &
577 2 : ,usr_ISOPHFP_aer_ndx,usr_IEPOX_aer_ndx,usr_DHPMPAL_aer_ndx,usr_ISOPB1O2_NOn_ndx,usr_ISOPB1O2_NOa_ndx &
578 2 : ,usr_ISOPB4O2_NOn_ndx,usr_ISOPB4O2_NOa_ndx,usr_ISOPED1O2_NOn_ndx,usr_ISOPED1O2_NOa_ndx &
579 2 : ,usr_ISOPED4O2_NOn_ndx,usr_ISOPED4O2_NOa_ndx,usr_ISOPZD1O2_NOn_ndx,usr_ISOPZD1O2_NOa_ndx &
580 2 : ,usr_ISOPZD4O2_NOn_ndx,usr_ISOPZD4O2_NOa_ndx,usr_ISOPNO3_NOn_ndx,usr_ISOPNO3_NOa_ndx &
581 2 : ,usr_MVKO2_NOn_ndx,usr_MVKO2_NOa_ndx,usr_MACRO2_NOn_ndx,usr_MACRO2_NOa_ndx &
582 2 : ,usr_IEPOXOO_NOn_ndx,usr_IEPOXOO_NOa_ndx,usr_ISOPN1DO2_NOn_ndx,usr_ISOPN1DO2_NOa_ndx &
583 2 : ,usr_ISOPN2BO2_NOn_ndx,usr_ISOPN2BO2_NOa_ndx,usr_ISOPN3BO2_NOn_ndx,usr_ISOPN3BO2_NOa_ndx &
584 2 : ,usr_ISOPN4DO2_NOn_ndx,usr_ISOPN4DO2_NOa_ndx,usr_ISOPNBNO3O2_NOn_ndx,usr_ISOPNBNO3O2_NOa_ndx &
585 2 : ,usr_ISOPNOOHBO2_NOn_ndx,usr_ISOPNOOHBO2_NOa_ndx,usr_ISOPNOOHDO2_NOn_ndx,usr_ISOPNOOHDO2_NOa_ndx &
586 2 : ,usr_NC4CHOO2_NOn_ndx,usr_NC4CHOO2_NOa_ndx,tag_MCO3_NO2_ndx &
587 2 : ,usr_TERPNT_aer_ndx,tag_TERPA2CO3_NO2_ndx,usr_TERPA2PAN_M_ndx &
588 2 : ,usr_TERPNT1_aer_ndx,usr_TERPNPT_aer_ndx,usr_TERPNPT1_aer_ndx,usr_TERPFDN_aer_ndx,usr_SQTN_aer_ndx &
589 2 : ,usr_TERPHFN_aer_ndx,usr_TERPDHDP_aer_ndx,usr_TERPACID_aer_ndx,tag_TERPACO3_NO2_ndx &
590 2 : ,usr_TERPAPAN_M_ndx,tag_TERPA3CO3_NO2_ndx, usr_TERPA3PAN_M_ndx,usr_ICHE_aer_ndx,usr_ISOPFNC_aer_ndx &
591 4 : ,usr_ISOPFDNC_aer_ndx
592 :
593 : end if
594 :
595 1536 : end subroutine usrrxt_inti
596 :
597 117648 : subroutine usrrxt( state, rxt, temp, tempi, tempe, invariants, h2ovmr, &
598 58824 : pmid, m, sulfate, mmr, relhum, strato_sad, &
599 58824 : tropchemlev, dlat, ncol, sad_trop, reff_trop, cwat, mbar, pbuf )
600 :
601 : !-----------------------------------------------------------------
602 : ! ... set the user specified reaction rates
603 : !-----------------------------------------------------------------
604 :
605 : use mo_constants, only : pi, avo => avogadro, boltz_cgs, rgas
606 : use chem_mods, only : nfs, rxntot, gas_pcnst, inv_m_ndx=>indexm
607 : use mo_setinv, only : inv_o2_ndx=>o2_ndx, inv_h2o_ndx=>h2o_ndx
608 : use physics_buffer, only : physics_buffer_desc
609 : use physics_types, only : physics_state
610 : use carma_flags_mod, only : carma_hetchem_feedback
611 : use aero_model, only : aero_model_surfarea
612 : use rad_constituents,only : rad_cnst_get_info
613 :
614 : implicit none
615 :
616 : !-----------------------------------------------------------------
617 : ! ... dummy arguments
618 : !-----------------------------------------------------------------
619 : integer, intent(in) :: ncol
620 : integer, intent(in) :: tropchemlev(pcols) ! trop/strat reaction separation vertical index
621 : real(r8), intent(in) :: dlat(:) ! degrees latitude
622 : real(r8), intent(in) :: temp(pcols,pver) ! temperature (K); neutral temperature
623 : real(r8), intent(in) :: tempi(pcols,pver) ! ionic temperature (K); only used if ion chemistry
624 : real(r8), intent(in) :: tempe(pcols,pver) ! electronic temperature (K); only used if ion chemistry
625 : real(r8), intent(in) :: m(ncol,pver) ! total atm density (/cm^3)
626 : real(r8), intent(in) :: sulfate(ncol,pver) ! sulfate aerosol (mol/mol)
627 : real(r8), intent(in) :: strato_sad(pcols,pver) ! stratospheric aerosol sad (1/cm)
628 : real(r8), intent(in) :: h2ovmr(ncol,pver) ! water vapor (mol/mol)
629 : real(r8), intent(in) :: relhum(ncol,pver) ! relative humidity
630 : real(r8), intent(in) :: pmid(pcols,pver) ! midpoint pressure (Pa)
631 : real(r8), intent(in) :: invariants(ncol,pver,nfs) ! invariants density (/cm^3)
632 : real(r8), intent(in) :: mmr(pcols,pver,gas_pcnst) ! species concentrations (kg/kg)
633 : real(r8), intent(in) :: cwat(ncol,pver) !PJC Condensed Water (liquid+ice) (kg/kg)
634 : real(r8), intent(in) :: mbar(ncol,pver) !PJC Molar mass of air (g/mol)
635 : real(r8), intent(inout) :: rxt(ncol,pver,rxntot) ! gas phase rates
636 : real(r8), intent(out) :: sad_trop(pcols,pver) ! tropospheric surface area density (cm2/cm3)
637 : real(r8), intent(out) :: reff_trop(pcols,pver) ! tropospheric effective radius (cm)
638 : type(physics_state), intent(in) :: state ! Physics state variables
639 : type(physics_buffer_desc), pointer :: pbuf(:)
640 :
641 : !-----------------------------------------------------------------
642 : ! ... local variables
643 : !-----------------------------------------------------------------
644 :
645 : real(r8), parameter :: dg = 0.1_r8 ! mole diffusion =0.1 cm2/s (Dentener, 1993)
646 :
647 : !-----------------------------------------------------------------
648 : ! ... reaction probabilities for heterogeneous reactions
649 : !-----------------------------------------------------------------
650 : real(r8), parameter :: gamma_n2o5 = 0.02_r8 ! JPL19
651 : real(r8), parameter :: gamma_ho2 = 0.10_r8 ! Gaubert et al., https://doi.org/10.5194/acp-20-14617-2020
652 : real(r8), parameter :: gamma_no2 = 8.0e-6_r8 ! Liu et al., Environ.Sci.&Tech, 53, 3517, 2019 doi:10.1021/acs.est.8b06367
653 : real(r8), parameter :: gamma_no3 = 0.002_r8 ! JPL19
654 : real(r8), parameter :: gamma_glyoxal = 2.0e-4_r8 ! Washenfelder et al, JGR, 2011
655 : !TS1 species
656 : real(r8), parameter :: gamma_isopnita = 0.005_r8 ! from Fisher et al., ACP, 2016
657 : real(r8), parameter :: gamma_isopnitb = 0.005_r8 !
658 : real(r8), parameter :: gamma_onitr = 0.005_r8 !
659 : real(r8), parameter :: gamma_honitr = 0.005_r8 !
660 : real(r8), parameter :: gamma_terpnit = 0.01_r8 !
661 : real(r8), parameter :: gamma_nterpooh = 0.01_r8 !
662 : real(r8), parameter :: gamma_nc4cho = 0.02_r8 !
663 : real(r8), parameter :: gamma_nc4ch2oh = 0.005_r8 !
664 : !TS2 species
665 : real(r8), parameter :: gamma_isopfdn = 0.1_r8 ! Marais 2015 for C5-LVOC
666 : real(r8), parameter :: gamma_isopfnp = 0.1_r8 ! Marais 2015 for C5-LVOC
667 : real(r8), parameter :: gamma_isopn2b = 0.02_r8 ! All isoprene nitrates Wolfe
668 : real(r8), parameter :: gamma_isopn1d = 0.02_r8 !
669 : real(r8), parameter :: gamma_isopn4d = 0.02_r8 !
670 : real(r8), parameter :: gamma_inoohd = 0.02_r8 !
671 : real(r8), parameter :: gamma_inheb = 0.02_r8 !Marais 2015 for IEPOX
672 : real(r8), parameter :: gamma_inhed = 0.02_r8 !Marais 2015 for IEPOX
673 : real(r8), parameter :: gamma_macrn = 0.02_r8 !
674 : real(r8), parameter :: gamma_isophfp = 0.1_r8 !Marais 2015 for C5-LVOC
675 : real(r8), parameter :: gamma_iepox = 0.0042_r8 !Marais 2015 for IEPOX
676 : real(r8), parameter :: gamma_dhpmpal = 0.1_r8 !Marais 2015 for C5-LVOC
677 : real(r8), parameter :: gamma_iche = 0.0042_r8 !Marais 2015 for IEPOX
678 : real(r8), parameter :: gamma_isopfnc = 0.1_r8 ! Marais 2015 for C5-LVOC
679 : real(r8), parameter :: gamma_isopfdnc = 0.1_r8 ! Marais 2015 for C5-LVOC
680 : real(r8), parameter :: gamma_terpnt = 0.02_r8 !
681 : real(r8), parameter :: gamma_terpnt1 = 0.02_r8 !
682 : real(r8), parameter :: gamma_terpnpt = 0.02_r8 !
683 : real(r8), parameter :: gamma_terpnpt1 = 0.02_r8 !
684 : real(r8), parameter :: gamma_terpfdn = 0.1_r8 !
685 : real(r8), parameter :: gamma_sqtn = 0.1_r8 !
686 : real(r8), parameter :: gamma_terphfn = 0.1_r8 !
687 : real(r8), parameter :: gamma_terpdhdp = 0.1_r8 !
688 : real(r8), parameter :: gamma_terpacid = 0.01_r8 !
689 :
690 : integer :: i, k
691 : integer :: l
692 117648 : real(r8) :: tp(ncol) ! 300/t
693 117648 : real(r8) :: tinv(ncol) ! 1/t
694 117648 : real(r8) :: ko(ncol)
695 117648 : real(r8) :: term1(ncol)
696 117648 : real(r8) :: term2(ncol)
697 117648 : real(r8) :: kinf(ncol)
698 117648 : real(r8) :: fc(ncol)
699 117648 : real(r8) :: xr(ncol)
700 117648 : real(r8) :: sur(ncol)
701 117648 : real(r8) :: sqrt_t(ncol) ! sqrt( temp )
702 117648 : real(r8) :: sqrt_t_58(ncol) ! sqrt( temp / 58.)
703 117648 : real(r8) :: exp_fac(ncol) ! vector exponential
704 117648 : real(r8) :: lwc(ncol)
705 117648 : real(r8) :: ko_m(ncol)
706 117648 : real(r8) :: k0(ncol)
707 117648 : real(r8) :: kinf_m(ncol)
708 117648 : real(r8) :: o2(ncol)
709 : real(r8) :: c_n2o5, c_ho2, c_no2, c_no3, c_glyoxal
710 : !TS1 species
711 : real(r8) :: c_isopnita, c_isopnitb, c_onitr, c_honitr, c_terpnit, c_nterpooh
712 : real(r8) :: c_nc4cho, c_nc4ch2oh
713 : !T2 species
714 : real(r8) :: c_isopfdn, c_isopfnp, c_isopn2b, c_isopn1d, c_isopn4d, c_inoohd
715 : real(r8) :: c_inheb, c_inhed, c_macrn, c_isophfp, c_iepox, c_dhpmpal
716 : real(r8) :: c_iche, c_isopfnc, c_isopfdnc
717 : real(r8) :: c_terpnt, c_terpnt1, c_terpnpt, c_terpnpt1, c_terpfdn, c_sqtn, c_terphfn, c_terpdhdp, c_terpacid
718 :
719 : real(r8) :: amas
720 : !-----------------------------------------------------------------
721 : ! ... density of sulfate aerosol
722 : !-----------------------------------------------------------------
723 : real(r8), parameter :: gam1 = 0.04_r8 ! N2O5+SUL ->2HNO3
724 : real(r8), parameter :: wso4 = 98._r8
725 : real(r8), parameter :: den = 1.15_r8 ! each molecule of SO4(aer) density g/cm3
726 : !-------------------------------------------------
727 : ! ... volume of sulfate particles
728 : ! assuming mean rm
729 : ! continient 0.05um 0.07um 0.09um
730 : ! ocean 0.09um 0.25um 0.37um
731 : ! 0.16um Blake JGR,7195, 1995
732 : !-------------------------------------------------
733 : real(r8), parameter :: rm1 = 0.16_r8*1.e-4_r8 ! mean radii in cm
734 : real(r8), parameter :: fare = 4._r8*pi*rm1*rm1 ! each mean particle(r=0.1u) area cm2/cm3
735 :
736 : !-----------------------------------------------------------------------
737 : ! ... Aqueous phase sulfur quantities for SO2 + H2O2 and SO2 + O3
738 : !-----------------------------------------------------------------------
739 : real(r8), parameter :: HENRY298_H2O2 = 7.45e+04_r8
740 : real(r8), parameter :: H298_H2O2 = -1.45e+04_r8
741 : real(r8), parameter :: HENRY298_SO2 = 1.23e+00_r8
742 : real(r8), parameter :: H298_SO2 = -6.25e+03_r8
743 : real(r8), parameter :: K298_SO2_HSO3 = 1.3e-02_r8
744 : real(r8), parameter :: H298_SO2_HSO3 = -4.16e+03_r8
745 : real(r8), parameter :: R_CONC = 82.05e+00_r8 / avo
746 : real(r8), parameter :: R_CAL = rgas * 0.239006e+00_r8
747 : real(r8), parameter :: K_AQ = 7.57e+07_r8
748 : real(r8), parameter :: ER_AQ = 4.43e+03_r8
749 :
750 : real(r8), parameter :: HENRY298_O3 = 1.13e-02_r8
751 : real(r8), parameter :: H298_O3 = -5.04e+03_r8
752 : real(r8), parameter :: K298_HSO3_SO3 = 6.6e-08_r8
753 : real(r8), parameter :: H298_HSO3_SO3 = -2.23e+03_r8
754 : real(r8), parameter :: K0_AQ = 2.4e+04_r8
755 : real(r8), parameter :: ER0_AQ = 0.0e+00_r8
756 : real(r8), parameter :: K1_AQ = 3.7e+05_r8
757 : real(r8), parameter :: ER1_AQ = 5.53e+03_r8
758 : real(r8), parameter :: K2_AQ = 1.5e+09_r8
759 : real(r8), parameter :: ER2_AQ = 5.28e+03_r8
760 :
761 : real(r8), parameter :: pH = 4.5e+00_r8
762 :
763 58824 : real(r8), pointer :: sfc(:), dm_aer(:)
764 : integer :: ntot_amode, nbins
765 :
766 58824 : real(r8), pointer :: sfc_array(:,:,:), dm_array(:,:,:)
767 : !TS2
768 117648 : real(r8) :: aterm(ncol)
769 : real(r8) :: natom
770 : real(r8) :: nyield
771 : real(r8) :: acorr
772 : real(r8) :: exp_natom
773 : character(len=*), parameter :: subname = 'usrrxt'
774 :
775 : ! get info about the modal aerosols
776 : ! get ntot_amode
777 58824 : call rad_cnst_get_info(0, nmodes=ntot_amode)
778 58824 : call rad_cnst_get_info(0, nbins=nbins)
779 :
780 58824 : if (ntot_amode>0.and.nbins>0) then
781 0 : call endrun(subname // ':: ERROR running with MAM and CARMA simultaneously not supported.')
782 : end if
783 :
784 58824 : if (ntot_amode>0) then
785 235296 : allocate(sfc_array(pcols,pver,ntot_amode), dm_array(pcols,pver,ntot_amode) )
786 0 : else if (nbins>0) then
787 0 : allocate(sfc_array(pcols,pver,nbins), dm_array(pcols,pver,nbins) )
788 : else
789 0 : allocate(sfc_array(pcols,pver,5), dm_array(pcols,pver,5) )
790 : endif
791 :
792 372297096 : sfc_array(:,:,:) = 0._r8
793 372297096 : dm_array(:,:,:) = 0._r8
794 58824 : sad_trop(:,:) = 0._r8
795 58824 : reff_trop(:,:) = 0._r8
796 :
797 58824 : if( usr_NO2_aer_ndx > 0 .or. usr_NO3_aer_ndx > 0 .or. usr_N2O5_aer_ndx > 0 .or. usr_HO2_aer_ndx > 0 ) then
798 :
799 : ! sad_trop should be set outside of usrrxt ??
800 0 : if( carma_hetchem_feedback ) then
801 0 : sad_trop(:ncol,:pver)=strato_sad(:ncol,:pver)
802 : else
803 :
804 : call aero_model_surfarea( &
805 : state, mmr, rm1, relhum, pmid, temp, strato_sad, sulfate, m, tropchemlev, dlat, &
806 0 : het1_ndx, pbuf, ncol, sfc_array, dm_array, sad_trop, reff_trop )
807 :
808 : endif
809 : endif
810 :
811 5529456 : level_loop : do k = 1,pver
812 91346832 : tinv(:) = 1._r8 / temp(:ncol,k)
813 91346832 : tp(:) = 300._r8 * tinv(:)
814 91346832 : sqrt_t(:) = sqrt( temp(:ncol,k) )
815 91346832 : sqrt_t_58(:) = sqrt( temp(:ncol,k) / 58.0_r8 )
816 :
817 : !-----------------------------------------------------------------
818 : ! ... o + o2 + m --> o3 + m (JPL15-10)
819 : !-----------------------------------------------------------------
820 5470632 : if( usr_O_O2_ndx > 0 ) then
821 0 : rxt(:,k,usr_O_O2_ndx) = 6.e-34_r8 * tp(:)**2.4_r8
822 : end if
823 5470632 : if( usr_OA_O2_ndx > 0 ) then
824 0 : rxt(:,k,usr_OA_O2_ndx) = 6.e-34_r8 * tp(:)**2.4_r8
825 : end if
826 :
827 : !-----------------------------------------------------------------
828 : ! ... o + o + m -> o2 + m
829 : !-----------------------------------------------------------------
830 5470632 : if ( usr_O_O_ndx > 0 ) then
831 0 : rxt(:,k,usr_O_O_ndx) = 2.76e-34_r8 * exp( 720.0_r8*tinv(:) )
832 : end if
833 :
834 : !-----------------------------------------------------------------
835 : ! ... cl2o2 + m -> 2*clo + m (JPL15-10)
836 : !-----------------------------------------------------------------
837 5470632 : if ( usr_CL2O2_M_ndx > 0 ) then
838 0 : if ( tag_CLO_CLO_M_ndx > 0 ) then
839 0 : ko(:) = 2.16e-27_r8 * exp( 8537.0_r8* tinv(:) )
840 0 : rxt(:,k,usr_CL2O2_M_ndx) = rxt(:,k,tag_CLO_CLO_M_ndx)/ko(:)
841 : else
842 0 : rxt(:,k,usr_CL2O2_M_ndx) = 0._r8
843 : end if
844 : end if
845 :
846 : !-----------------------------------------------------------------
847 : ! ... so3 + 2*h2o --> h2so4 + h2o
848 : ! Note: this reaction proceeds by the 2 intermediate steps below
849 : ! so3 + h2o --> adduct
850 : ! adduct + h2o --> h2so4 + h2o
851 : ! (Lovejoy et al., JCP, pp. 19911-19916, 1996)
852 : ! The first order rate constant used here is recommended by JPL 2011.
853 : ! This rate involves the water vapor number density.
854 : !-----------------------------------------------------------------
855 :
856 5470632 : if ( usr_SO3_H2O_ndx > 0 ) then
857 0 : call comp_exp( exp_fac, 6540.0_r8*tinv(:), ncol )
858 0 : if( h2o_ndx > 0 ) then
859 0 : fc(:) = 8.5e-21_r8 * m(:,k) * h2ovmr(:,k) * exp_fac(:)
860 : else
861 0 : fc(:) = 8.5e-21_r8 * invariants(:,k,inv_h2o_ndx) * exp_fac(:)
862 : end if
863 0 : rxt(:,k,usr_SO3_H2O_ndx) = 1.0e-20_r8 * fc(:)
864 : end if
865 :
866 : !-----------------------------------------------------------------
867 : ! ... n2o5 + m --> no2 + no3 + m (JPL15-10)
868 : !-----------------------------------------------------------------
869 5470632 : if( usr_N2O5_M_ndx > 0 ) then
870 0 : if( tag_NO2_NO3_ndx > 0 ) then
871 0 : call comp_exp( exp_fac, -10840.0_r8*tinv, ncol )
872 0 : rxt(:,k,usr_N2O5_M_ndx) = rxt(:,k,tag_NO2_NO3_ndx) * 1.724138e26_r8 * exp_fac(:)
873 : else
874 0 : rxt(:,k,usr_N2O5_M_ndx) = 0._r8
875 : end if
876 : end if
877 5470632 : if( usr_XNO2NO3_M_ndx > 0 ) then
878 0 : if( tag_NO2_NO3_ndx > 0 ) then
879 0 : call comp_exp( exp_fac, -10840.0_r8*tinv, ncol )
880 0 : rxt(:,k,usr_XNO2NO3_M_ndx) = rxt(:,k,tag_NO2_NO3_ndx) *1.724138e26_r8 * exp_fac(:)
881 : else
882 0 : rxt(:,k,usr_XNO2NO3_M_ndx) = 0._r8
883 : end if
884 : end if
885 5470632 : if( usr_NO2XNO3_M_ndx > 0 ) then
886 0 : if( tag_NO2_NO3_ndx > 0 ) then
887 0 : call comp_exp( exp_fac, -10840.0_r8*tinv, ncol )
888 0 : rxt(:,k,usr_NO2XNO3_M_ndx) = rxt(:,k,tag_NO2_NO3_ndx) * 1.734138e26_r8 * exp_fac(:)
889 : else
890 0 : rxt(:,k,usr_NO2XNO3_M_ndx) = 0._r8
891 : end if
892 : end if
893 :
894 : !-----------------------------------------------------------------
895 : ! set rates for:
896 : ! ... hno3 + oh --> no3 + h2o
897 : ! ho2no2 + m --> ho2 + no2 + m
898 : !-----------------------------------------------------------------
899 5470632 : if( usr_HNO3_OH_ndx > 0 ) then
900 0 : call comp_exp( exp_fac, 1335._r8*tinv, ncol )
901 0 : ko(:) = m(:,k) * 6.5e-34_r8 * exp_fac(:)
902 0 : call comp_exp( exp_fac, 2199._r8*tinv, ncol )
903 0 : ko(:) = ko(:) / (1._r8 + ko(:)/(2.7e-17_r8*exp_fac(:)))
904 0 : call comp_exp( exp_fac, 460._r8*tinv, ncol )
905 0 : rxt(:,k,usr_HNO3_OH_ndx) = ko(:) + 2.4e-14_r8*exp_fac(:)
906 : end if
907 5470632 : if( usr_XHNO3_OH_ndx > 0 ) then
908 0 : call comp_exp( exp_fac, 1335._r8*tinv, ncol )
909 0 : ko(:) = m(:,k) * 6.5e-34_r8 * exp_fac(:)
910 0 : call comp_exp( exp_fac, 2199._r8*tinv, ncol )
911 0 : ko(:) = ko(:) / (1._r8 + ko(:)/(2.7e-17_r8*exp_fac(:)))
912 0 : call comp_exp( exp_fac, 460._r8*tinv, ncol )
913 0 : rxt(:,k,usr_XHNO3_OH_ndx) = ko(:) + 2.4e-14_r8*exp_fac(:)
914 : end if
915 5470632 : if( usr_HO2NO2_M_ndx > 0 ) then
916 0 : if( tag_NO2_HO2_ndx > 0 ) then
917 0 : call comp_exp( exp_fac, -10900._r8*tinv, ncol )
918 0 : rxt(:,k,usr_HO2NO2_M_ndx) = rxt(:,k,tag_NO2_HO2_ndx) * exp_fac(:) / 2.1e-27_r8
919 : else
920 0 : rxt(:,k,usr_HO2NO2_M_ndx) = 0._r8
921 : end if
922 : end if
923 5470632 : if( usr_XHO2NO2_M_ndx > 0 ) then
924 0 : if( tag_NO2_HO2_ndx > 0 ) then
925 0 : call comp_exp( exp_fac, -10900._r8*tinv, ncol )
926 0 : rxt(:,k,usr_XHO2NO2_M_ndx) = rxt(:,k,tag_NO2_HO2_ndx) * exp_fac(:) / 2.1e-27_r8
927 : else
928 0 : rxt(:,k,usr_XHO2NO2_M_ndx) = 0._r8
929 : end if
930 : end if
931 : !-----------------------------------------------------------------
932 : ! ... co + oh --> co2 + ho2 (new single reaction for combined branches [JPL19])
933 : !-----------------------------------------------------------------
934 5470632 : if( usr_CO_OH_ndx > 0 ) then
935 0 : ko (:) = 6.9e-33_r8 * ( 298._r8 / temp(:ncol,k) )**(2.1_r8)
936 0 : kinf(:) = 1.1e-12_r8 * ( 298._r8 / temp(:ncol,k) )**(-1.3_r8)
937 :
938 0 : term2(:) = (1 + (log10( ko(:)*m(:,k) / kinf(:) ))**2)**(-1)
939 :
940 0 : term1(:) = (kinf(:) * ko(:)*m(:,k)) / (kinf(:) + ko(:)*m(:,k)) * (0.6_r8)**term2(:)
941 :
942 0 : rxt(:ncol,k,usr_CO_OH_ndx) = term1(:) + 1.85e-13_r8 * exp(-65._r8/temp(:ncol,k)) * (1._r8 - term1(:)/kinf(:))
943 :
944 : end if
945 : !-----------------------------------------------------------------
946 : ! ... co + oh --> co2 + ho2 (combined branches - do not use with CO_OH_b)
947 : ! note: for mechanisms prior to Dec 2022
948 : !-----------------------------------------------------------------
949 5470632 : if( usr_CO_OH_a_ndx > 0 ) then
950 0 : rxt(:,k,usr_CO_OH_a_ndx) = 1.5e-13_r8 * &
951 0 : (1._r8 + 6.e-7_r8*boltz_cgs*m(:,k)*temp(:ncol,k))
952 : end if
953 : !-----------------------------------------------------------------
954 : ! ... co + oh --> co2 + h (second branch JPL15-10, with CO+OH+M)
955 : ! note: for mechanisms prior to Dec 2022
956 : !-----------------------------------------------------------------
957 5470632 : if( usr_CO_OH_b_ndx > 0 ) then
958 0 : kinf(:) = 2.1e+09_r8 * (temp(:ncol,k)/ t0)**(6.1_r8)
959 0 : ko (:) = 1.5e-13_r8
960 :
961 0 : term1(:) = ko(:) / ( (kinf(:) / m(:,k)) )
962 0 : term2(:) = ko(:) / (1._r8 + term1(:))
963 :
964 0 : term1(:) = log10( term1(:) )
965 0 : term1(:) = 1.0_r8 / (1.0_r8 + term1(:)*term1(:))
966 :
967 0 : rxt(:ncol,k,usr_CO_OH_b_ndx) = term2(:) * (0.6_r8)**term1(:)
968 : end if
969 :
970 : !-----------------------------------------------------------------
971 : ! ... ho2 + ho2 --> h2o2
972 : ! note: this rate involves the water vapor number density
973 : !-----------------------------------------------------------------
974 5470632 : if( usr_HO2_HO2_ndx > 0 ) then
975 :
976 91346832 : call comp_exp( exp_fac, 460._r8*tinv, ncol )
977 91346832 : ko(:) = 3.0e-13_r8 * exp_fac(:)
978 91346832 : call comp_exp( exp_fac, 920._r8*tinv, ncol )
979 91346832 : kinf(:) = 2.1e-33_r8 * m(:,k) * exp_fac(:)
980 91346832 : call comp_exp( exp_fac, 2200._r8*tinv, ncol )
981 :
982 5470632 : if( h2o_ndx > 0 ) then
983 91346832 : fc(:) = 1._r8 + 1.4e-21_r8 * m(:,k) * h2ovmr(:,k) * exp_fac(:)
984 : else
985 0 : fc(:) = 1._r8 + 1.4e-21_r8 * invariants(:,k,inv_h2o_ndx) * exp_fac(:)
986 : end if
987 91346832 : rxt(:,k,usr_HO2_HO2_ndx) = (ko(:) + kinf(:)) * fc(:)
988 :
989 : end if
990 :
991 : !-----------------------------------------------------------------
992 : ! ... mco3 + no2 -> mpan
993 : !-----------------------------------------------------------------
994 5470632 : if( usr_MCO3_NO2_ndx > 0 ) then
995 0 : rxt(:,k,usr_MCO3_NO2_ndx) = 1.1e-11_r8 * tp(:) / m(:,k)
996 : end if
997 5470632 : if( usr_MCO3_XNO2_ndx > 0 ) then
998 0 : rxt(:,k,usr_MCO3_XNO2_ndx) = 1.1e-11_r8 * tp(:) / m(:,k)
999 : end if
1000 :
1001 : !-----------------------------------------------------------------
1002 : ! ... pan + m --> ch3co3 + no2 + m (JPL15-10)
1003 : !-----------------------------------------------------------------
1004 91346832 : call comp_exp( exp_fac, -14000._r8*tinv, ncol )
1005 5470632 : if( usr_PAN_M_ndx > 0 ) then
1006 0 : if( tag_CH3CO3_NO2_ndx > 0 ) then
1007 0 : rxt(:,k,usr_PAN_M_ndx) = rxt(:,k,tag_CH3CO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1008 : else
1009 0 : rxt(:,k,usr_PAN_M_ndx) = 0._r8
1010 : end if
1011 : end if
1012 5470632 : if( usr_XPAN_M_ndx > 0 ) then
1013 0 : if( tag_CH3CO3_NO2_ndx > 0 ) then
1014 0 : rxt(:,k,usr_XPAN_M_ndx) = rxt(:,k,tag_CH3CO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1015 : else
1016 0 : rxt(:,k,usr_XPAN_M_ndx) = 0._r8
1017 : end if
1018 : end if
1019 :
1020 : !-----------------------------------------------------------------
1021 : ! ... mpan + m --> mco3 + no2 + m (JPL15-10)
1022 : !-----------------------------------------------------------------
1023 5470632 : if( usr_MPAN_M_ndx > 0 ) then
1024 0 : if( tag_MCO3_NO2_ndx > 0 ) then
1025 0 : rxt(:,k,usr_MPAN_M_ndx) = rxt(:,k,tag_MCO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1026 : else
1027 0 : rxt(:,k,usr_MPAN_M_ndx) = 0._r8
1028 : end if
1029 : end if
1030 5470632 : if( usr_XMPAN_M_ndx > 0 ) then
1031 0 : if( tag_MCO3_NO2_ndx > 0 ) then
1032 0 : rxt(:,k,usr_XMPAN_M_ndx) = rxt(:,k,tag_MCO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1033 : else
1034 0 : rxt(:,k,usr_XMPAN_M_ndx) = 0._r8
1035 : end if
1036 : end if
1037 :
1038 : !lke-TS1
1039 : !-----------------------------------------------------------------
1040 : ! ... pbznit + m --> acbzo2 + no2 + m
1041 : !-----------------------------------------------------------------
1042 5470632 : if( usr_PBZNIT_M_ndx > 0 ) then
1043 0 : if( tag_ACBZO2_NO2_ndx > 0 ) then
1044 0 : rxt(:,k,usr_PBZNIT_M_ndx) = rxt(:,k,tag_ACBZO2_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1045 : else
1046 0 : rxt(:,k,usr_PBZNIT_M_ndx) = 0._r8
1047 : end if
1048 : end if
1049 : !-----------------------------------------------------------------
1050 : ! ... TERPAPAN + m --> TERPACO3 + no2 + m
1051 : !-----------------------------------------------------------------
1052 5470632 : if( usr_TERPAPAN_M_ndx > 0 ) then
1053 0 : if( tag_TERPACO3_NO2_ndx > 0 ) then
1054 0 : rxt(:,k,usr_TERPAPAN_M_ndx) = rxt(:,k,tag_TERPACO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1055 : else
1056 0 : rxt(:,k,usr_TERPAPAN_M_ndx) = 0._r8
1057 : end if
1058 : end if
1059 : !-----------------------------------------------------------------
1060 : ! ... TERPA2PAN + m --> TERPA2CO3 + no2 + m
1061 : !-----------------------------------------------------------------
1062 5470632 : if( usr_TERPA2PAN_M_ndx > 0 ) then
1063 0 : if( tag_TERPA2CO3_NO2_ndx > 0 ) then
1064 0 : rxt(:,k,usr_TERPA2PAN_M_ndx) = rxt(:,k,tag_TERPA2CO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1065 : else
1066 0 : rxt(:,k,usr_TERPA2PAN_M_ndx) = 0._r8
1067 : end if
1068 : end if
1069 : !-----------------------------------------------------------------
1070 : ! ... TERPA3PAN + m --> TERPA3CO3 + no2 + m
1071 : !-----------------------------------------------------------------
1072 5470632 : if( usr_TERPA3PAN_M_ndx > 0 ) then
1073 0 : if( tag_TERPA3CO3_NO2_ndx > 0 ) then
1074 0 : rxt(:,k,usr_TERPA3PAN_M_ndx) = rxt(:,k,tag_TERPA3CO3_NO2_ndx) * 1.111e28_r8 * exp_fac(:)
1075 : else
1076 0 : rxt(:,k,usr_TERPA3PAN_M_ndx) = 0._r8
1077 : end if
1078 : end if
1079 : !-----------------------------------------------------------------
1080 : ! ... xooh + oh -> h2o + oh
1081 : !-----------------------------------------------------------------
1082 5470632 : if( usr_XOOH_OH_ndx > 0 ) then
1083 0 : call comp_exp( exp_fac, 253._r8*tinv, ncol )
1084 0 : rxt(:,k,usr_XOOH_OH_ndx) = temp(:ncol,k)**2._r8 * 7.69e-17_r8 * exp_fac(:)
1085 : end if
1086 :
1087 : !-----------------------------------------------------------------
1088 : ! ... ch3coch3 + oh -> ro2 + h2o
1089 : !-----------------------------------------------------------------
1090 5470632 : if( usr_CH3COCH3_OH_ndx > 0 ) then
1091 0 : call comp_exp( exp_fac, -2000._r8*tinv, ncol )
1092 0 : rxt(:,k,usr_CH3COCH3_OH_ndx) = 3.82e-11_r8 * exp_fac(:) + 1.33e-13_r8
1093 : end if
1094 :
1095 : !-----------------------------------------------------------------
1096 : ! ... DMS + OH --> .5 * SO2
1097 : ! JPL15-10 (use [O2] = 0.21*[M])
1098 : ! k = 8.2E-39 * exp(5376/T) * [O2] / (1 + 1.05E-5 *([O2]/[M]) * exp(3644/T))
1099 : !-----------------------------------------------------------------
1100 5470632 : if( usr_DMS_OH_ndx > 0 ) then
1101 91346832 : call comp_exp( exp_fac, 3644._r8*tinv, ncol )
1102 91346832 : ko(:) = 1._r8 + 1.05e-5_r8 * exp_fac * 0.21_r8
1103 91346832 : call comp_exp( exp_fac, 5376._r8*tinv, ncol )
1104 91346832 : rxt(:,k,usr_DMS_OH_ndx) = 8.2e-39_r8 * exp_fac * m(:,k) * 0.21_r8 / ko(:)
1105 : end if
1106 :
1107 : !-----------------------------------------------------------------
1108 : ! ... SO2 + OH --> SO4 (REFERENCE?? - not Liao)
1109 : !-----------------------------------------------------------------
1110 5470632 : if( usr_SO2_OH_ndx > 0 ) then
1111 0 : fc(:) = 3.0e-31_r8 *(300._r8*tinv(:))**3.3_r8
1112 0 : ko(:) = fc(:)*m(:,k)/(1._r8 + fc(:)*m(:,k)/1.5e-12_r8)
1113 0 : rxt(:,k,usr_SO2_OH_ndx) = ko(:)*.6_r8**(1._r8 + (log10(fc(:)*m(:,k)/1.5e-12_r8))**2._r8)**(-1._r8)
1114 : end if
1115 : !RHS TS2
1116 : !-----------------------------------------------------------------
1117 : ! ... ISOPZD1O2 --> HPALD etc. Wennberg 2018 for rate
1118 : !-----------------------------------------------------------------
1119 5470632 : if( usr_ISOPZD1O2_ndx > 0 ) then
1120 0 : call comp_exp( exp_fac, -12200._r8*tinv, ncol )
1121 0 : ko(:) = 5.05e15_r8 * exp_fac(:)
1122 0 : call comp_exp( exp_fac, 1.e8_r8*tinv**3._r8, ncol )
1123 0 : rxt(:,k,usr_ISOPZD1O2_ndx) = ko(:)*exp_fac(:)
1124 : end if
1125 : !-----------------------------------------------------------------
1126 : ! ... ISOPZD4O2 --> HPALD etc. Wennberg 2018 for rate
1127 : !-----------------------------------------------------------------
1128 5470632 : if( usr_ISOPZD4O2_ndx > 0 ) then
1129 0 : call comp_exp( exp_fac, -7160._r8*tinv, ncol )
1130 0 : ko(:) = 2.22e9_r8 * exp_fac(:)
1131 0 : call comp_exp( exp_fac, 1.e8_r8*tinv**3._r8, ncol )
1132 0 : rxt(:,k,usr_ISOPZD4O2_ndx) = ko(:)*exp_fac(:)
1133 : end if
1134 : !-----------------------------------------------------------------
1135 : ! ... ISOPB1O2_NOn Temp/Pressure Dependent Nitrate Yield
1136 : !-----------------------------------------------------------------
1137 5470632 : if( usr_ISOPB1O2_NOn_ndx > 0 ) then
1138 0 : nyield = (1._r8-0.14_r8)/0.14_r8
1139 0 : natom = 6.0_r8
1140 0 : exp_natom = exp( natom )
1141 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1142 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1143 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1144 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1145 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1146 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1147 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1148 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1149 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1150 0 : rxt(:,k,usr_ISOPB1O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1151 0 : rxt(:,k,usr_ISOPB1O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1152 : end if
1153 : !-----------------------------------------------------------------
1154 : ! ... ISOPB4O2_NOn Temp/Pressure Dependent Nitrate Yield
1155 : !-----------------------------------------------------------------
1156 5470632 : if( usr_ISOPB4O2_NOn_ndx > 0 ) then
1157 0 : nyield = (1._r8-0.13_r8)/0.13_r8
1158 0 : natom = 6.0_r8
1159 0 : exp_natom = exp( natom )
1160 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1161 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1162 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1163 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1164 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1165 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1166 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1167 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1168 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1169 0 : rxt(:,k,usr_ISOPB4O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1170 0 : rxt(:,k,usr_ISOPB4O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1171 : end if
1172 : !-----------------------------------------------------------------
1173 : ! ... ISOPED1O2_NOn Temp/Pressure Dependent Nitrate Yield
1174 : !-----------------------------------------------------------------
1175 5470632 : if( usr_ISOPED1O2_NOn_ndx > 0 ) then
1176 0 : nyield = (1._r8-0.12_r8)/0.12_r8
1177 0 : natom = 6.0_r8
1178 0 : exp_natom = exp( natom )
1179 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1180 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1181 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1182 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1183 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1184 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1185 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1186 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1187 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1188 0 : rxt(:,k,usr_ISOPED1O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1189 0 : rxt(:,k,usr_ISOPED1O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1190 : end if
1191 : !-----------------------------------------------------------------
1192 : ! ... ISOPED4O2_NOn Temp/Pressure Dependent Nitrate Yield
1193 : !-----------------------------------------------------------------
1194 5470632 : if( usr_ISOPED4O2_NOn_ndx > 0 ) then
1195 0 : nyield = (1._r8-0.12_r8)/0.12_r8
1196 0 : natom = 6.0_r8
1197 0 : exp_natom = exp( natom )
1198 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1199 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1200 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1201 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1202 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1203 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1204 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1205 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1206 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1207 0 : rxt(:,k,usr_ISOPED4O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1208 0 : rxt(:,k,usr_ISOPED4O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1209 : end if
1210 : !-----------------------------------------------------------------
1211 : ! ... ISOPZD1O2_NOn Temp/Pressure Dependent Nitrate Yield
1212 : !-----------------------------------------------------------------
1213 5470632 : if( usr_ISOPZD1O2_NOn_ndx > 0 ) then
1214 0 : nyield = (1._r8-0.12_r8)/0.12_r8
1215 0 : natom = 6.0_r8
1216 0 : exp_natom = exp( natom )
1217 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1218 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1219 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1220 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1221 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1222 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1223 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1224 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1225 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1226 0 : rxt(:,k,usr_ISOPZD1O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1227 0 : rxt(:,k,usr_ISOPZD1O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1228 : end if
1229 : !-----------------------------------------------------------------
1230 : ! ... ISOPZD4O2_NOn Temp/Pressure Dependent Nitrate Yield
1231 : !-----------------------------------------------------------------
1232 5470632 : if( usr_ISOPZD4O2_NOn_ndx > 0 ) then
1233 0 : nyield = (1._r8-0.12_r8)/0.12_r8
1234 0 : natom = 6.0_r8
1235 0 : exp_natom = exp( natom )
1236 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1237 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1238 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1239 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1240 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1241 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1242 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1243 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1244 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1245 0 : rxt(:,k,usr_ISOPZD4O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1246 0 : rxt(:,k,usr_ISOPZD4O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1247 : end if
1248 : !-----------------------------------------------------------------
1249 : ! ... ISOPNO3_NOn Temp/Pressure Dependent Nitrate Yield
1250 : !-----------------------------------------------------------------
1251 5470632 : if( usr_ISOPNO3_NOn_ndx > 0 ) then
1252 0 : nyield = (1._r8-0.135_r8)/0.135_r8
1253 0 : natom = 9.0_r8
1254 0 : exp_natom = exp( natom )
1255 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1256 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1257 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1258 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1259 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1260 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1261 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1262 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1263 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1264 0 : rxt(:,k,usr_ISOPNO3_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1265 0 : rxt(:,k,usr_ISOPNO3_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1266 : end if
1267 : !-----------------------------------------------------------------
1268 : ! ... MVKO2_NOn Temp/Pressure Dependent Nitrate Yield
1269 : !-----------------------------------------------------------------
1270 5470632 : if( usr_MVKO2_NOn_ndx > 0 ) then
1271 0 : nyield = (1._r8-0.04_r8)/0.04_r8
1272 0 : natom = 6.0_r8
1273 0 : exp_natom = exp( natom )
1274 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1275 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1276 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1277 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1278 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1279 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1280 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1281 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1282 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1283 0 : rxt(:,k,usr_MVKO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1284 0 : rxt(:,k,usr_MVKO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1285 : end if
1286 : !-----------------------------------------------------------------
1287 : ! ... MACRO2_NOn Temp/Pressure Dependent Nitrate Yield
1288 : !-----------------------------------------------------------------
1289 5470632 : if( usr_MACRO2_NOn_ndx > 0 ) then
1290 0 : nyield = (1._r8-0.06_r8)/0.06_r8
1291 0 : natom = 6.0_r8
1292 0 : exp_natom = exp( natom )
1293 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1294 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1295 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1296 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1297 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1298 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1299 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1300 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1301 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1302 0 : rxt(:,k,usr_MACRO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1303 0 : rxt(:,k,usr_MACRO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1304 : end if
1305 : !-----------------------------------------------------------------
1306 : ! ... IEPOXOO_NOn Temp/Pressure Dependent Nitrate Yield
1307 : !-----------------------------------------------------------------
1308 5470632 : if( usr_IEPOXOO_NOn_ndx > 0 ) then
1309 0 : nyield = (1._r8-0.025_r8)/0.025_r8
1310 0 : natom = 8.0_r8
1311 0 : exp_natom = exp( natom )
1312 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1313 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1314 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1315 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1316 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1317 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1318 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1319 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1320 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1321 0 : rxt(:,k,usr_IEPOXOO_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1322 0 : rxt(:,k,usr_IEPOXOO_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1323 : end if
1324 : !-----------------------------------------------------------------
1325 : ! ... ISOPN1DO2_NOn Temp/Pressure Dependent Nitrate Yield
1326 : !-----------------------------------------------------------------
1327 5470632 : if( usr_ISOPN1DO2_NOn_ndx > 0 ) then
1328 0 : nyield = (1._r8-0.084_r8)/0.084_r8
1329 0 : natom = 11.0_r8
1330 0 : exp_natom = exp( natom )
1331 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1332 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1333 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1334 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1335 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1336 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1337 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1338 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1339 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1340 0 : rxt(:,k,usr_ISOPN1DO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1341 0 : rxt(:,k,usr_ISOPN1DO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1342 : end if
1343 : !-----------------------------------------------------------------
1344 : ! ... ISOPN2BO2_NOn Temp/Pressure Dependent Nitrate Yield
1345 : !-----------------------------------------------------------------
1346 5470632 : if( usr_ISOPN2BO2_NOn_ndx > 0 ) then
1347 0 : nyield = (1._r8-0.065_r8)/0.065_r8
1348 0 : natom = 11.0_r8
1349 0 : exp_natom = exp( natom )
1350 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1351 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1352 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1353 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1354 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1355 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1356 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1357 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1358 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1359 0 : rxt(:,k,usr_ISOPN2BO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1360 0 : rxt(:,k,usr_ISOPN2BO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1361 : end if
1362 : !-----------------------------------------------------------------
1363 : ! ... ISOPN3BO2_NOn Temp/Pressure Dependent Nitrate Yield
1364 : !-----------------------------------------------------------------
1365 5470632 : if( usr_ISOPN3BO2_NOn_ndx > 0 ) then
1366 0 : nyield = (1._r8-0.053_r8)/0.053_r8
1367 0 : natom = 11.0_r8
1368 0 : exp_natom = exp( natom )
1369 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1370 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1371 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1372 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1373 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1374 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1375 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1376 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1377 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1378 0 : rxt(:,k,usr_ISOPN3BO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1379 0 : rxt(:,k,usr_ISOPN3BO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1380 : end if
1381 : !-----------------------------------------------------------------
1382 : ! ... ISOPN4DO2_NOn Temp/Pressure Dependent Nitrate Yield
1383 : !-----------------------------------------------------------------
1384 5470632 : if( usr_ISOPN4DO2_NOn_ndx > 0 ) then
1385 0 : nyield = (1._r8-0.165_r8)/0.165_r8
1386 0 : natom = 11.0_r8
1387 0 : exp_natom = exp( natom )
1388 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1389 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1390 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1391 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1392 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1393 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1394 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1395 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1396 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1397 0 : rxt(:,k,usr_ISOPN4DO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1398 0 : rxt(:,k,usr_ISOPN4DO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1399 : end if
1400 : !-----------------------------------------------------------------
1401 : ! ... ISOPNBNO3O2_NOn Temp/Pressure Dependent Nitrate Yield
1402 : !-----------------------------------------------------------------
1403 5470632 : if( usr_ISOPNBNO3O2_NOn_ndx > 0 ) then
1404 0 : nyield = (1._r8-0.203_r8)/0.203_r8
1405 0 : natom = 11.0_r8
1406 0 : exp_natom = exp( natom )
1407 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1408 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1409 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1410 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1411 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1412 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1413 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1414 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1415 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1416 0 : rxt(:,k,usr_ISOPNBNO3O2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1417 0 : rxt(:,k,usr_ISOPNBNO3O2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1418 : end if
1419 : !-----------------------------------------------------------------
1420 : ! ... ISOPNOOHBO2_NOn Temp/Pressure Dependent Nitrate Yield
1421 : !-----------------------------------------------------------------
1422 5470632 : if( usr_ISOPNOOHBO2_NOn_ndx > 0 ) then
1423 0 : nyield = (1._r8-0.141_r8)/0.141_r8
1424 0 : natom = 12.0_r8
1425 0 : exp_natom = exp( natom )
1426 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1427 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1428 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1429 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1430 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1431 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1432 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1433 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1434 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1435 0 : rxt(:,k,usr_ISOPNOOHBO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1436 0 : rxt(:,k,usr_ISOPNOOHBO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1437 : end if
1438 : !-----------------------------------------------------------------
1439 : ! ... ISOPNOOHDO2_NOn Temp/Pressure Dependent Nitrate Yield
1440 : !-----------------------------------------------------------------
1441 5470632 : if( usr_ISOPNOOHDO2_NOn_ndx > 0 ) then
1442 0 : nyield = (1._r8-0.045_r8)/0.045_r8
1443 0 : natom = 12.0_r8
1444 0 : exp_natom = exp( natom )
1445 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1446 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1447 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1448 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1449 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1450 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1451 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1452 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1453 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1454 0 : rxt(:,k,usr_ISOPNOOHDO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1455 0 : rxt(:,k,usr_ISOPNOOHDO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1456 : end if
1457 : !-----------------------------------------------------------------
1458 : ! ... NC4CHOO2_NOn Temp/Pressure Dependent Nitrate Yield
1459 : !-----------------------------------------------------------------
1460 5470632 : if( usr_NC4CHOO2_NOn_ndx > 0 ) then
1461 0 : nyield = (1._r8-0.021_r8)/0.021_r8
1462 0 : natom = 11.0_r8
1463 0 : exp_natom = exp( natom )
1464 : acorr = (2.0e-22_r8*exp_natom*2.45e19_r8)/(1._r8+((2.0e-22_r8* &
1465 : exp_natom*2.45e19_r8)/(0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))* &
1466 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*2.45e19_r8)/ &
1467 0 : (0.43_r8*(298._r8*(1._r8/293._r8))**8._r8)))**2._r8))
1468 : aterm(:) = (2.0e-22_r8*exp_natom*m(:,k))/(1._r8+((2.0e-22_r8* &
1469 : exp_natom*m(:,k))/(0.43_r8*(298._r8*tinv(:))**8._r8)))* &
1470 : 0.41_r8**(1._r8/(1._r8+(log10((2.0e-22_r8*exp_natom*m(:,k))/ &
1471 0 : (0.43_r8*(298._r8*tinv(:))**8._r8)))**2._r8))
1472 0 : call comp_exp( exp_fac, 360._r8*tinv, ncol )
1473 0 : rxt(:,k,usr_NC4CHOO2_NOn_ndx) = 2.7e-12_r8 * exp_fac(:)*aterm(:)/(aterm(:)+acorr*nyield)
1474 0 : rxt(:,k,usr_NC4CHOO2_NOa_ndx) = 2.7e-12_r8 * exp_fac(:)*acorr*nyield/(aterm(:)+acorr*nyield)
1475 : end if
1476 : !
1477 : ! reduced hydrocarbon scheme
1478 : !
1479 5470632 : if ( usr_C2O3_NO2_ndx > 0 ) then
1480 0 : ko(:) = 2.6e-28_r8 * m(:,k)
1481 0 : kinf(:) = 1.2e-11_r8
1482 0 : rxt(:,k,usr_C2O3_NO2_ndx) = (ko/(1._r8+ko/kinf)) * 0.6_r8**(1._r8/(1._r8+(log10(ko/kinf))**2))
1483 : end if
1484 5470632 : if ( usr_C2O3_XNO2_ndx > 0 ) then
1485 0 : ko(:) = 2.6e-28_r8 * m(:,k)
1486 0 : kinf(:) = 1.2e-11_r8
1487 0 : rxt(:,k,usr_C2O3_XNO2_ndx) = (ko/(1._r8+ko/kinf)) * 0.6_r8**(1._r8/(1._r8+(log10(ko/kinf))**2))
1488 : end if
1489 5470632 : if ( usr_C2H4_OH_ndx > 0 ) then
1490 0 : ko(:) = 1.0e-28_r8 * m(:,k)
1491 0 : kinf(:) = 8.8e-12_r8
1492 0 : rxt(:,k,usr_C2H4_OH_ndx) = (ko/(1._r8+ko/kinf)) * 0.6_r8**(1._r8/(1._r8+(log(ko/kinf))**2))
1493 : end if
1494 5470632 : if ( usr_XO2N_HO2_ndx > 0 ) then
1495 0 : rxt(:,k,usr_XO2N_HO2_ndx) = rxt(:,k,tag_XO2N_NO_ndx)*rxt(:,k,tag_XO2_HO2_ndx)/(rxt(:,k,tag_XO2_NO_ndx)+1.e-36_r8)
1496 : end if
1497 :
1498 : !
1499 : ! hydrolysis reactions on wetted aerosols
1500 : !
1501 : if( usr_NO2_aer_ndx > 0 .or. usr_NO3_aer_ndx > 0 .or. usr_N2O5_aer_ndx > 0 .or. usr_HO2_aer_ndx > 0 &
1502 5470632 : .or. usr_GLYOXAL_aer_ndx > 0 ) then
1503 :
1504 0 : long_loop : do i = 1,ncol
1505 :
1506 0 : sfc => sfc_array(i,k,:)
1507 0 : dm_aer => dm_array(i,k,:)
1508 :
1509 0 : c_n2o5 = 1.40e3_r8 * sqrt_t(i) ! mean molecular speed of n2o5
1510 0 : c_no3 = 1.85e3_r8 * sqrt_t(i) ! mean molecular speed of no3
1511 0 : c_no2 = 2.15e3_r8 * sqrt_t(i) ! mean molecular speed of no2
1512 0 : c_ho2 = 2.53e3_r8 * sqrt_t(i) ! mean molecular speed of ho2
1513 0 : c_glyoxal = 1.455e4_r8 * sqrt_t_58(i) ! mean molecular speed of ho2
1514 0 : c_isopnita = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of isopnita
1515 0 : c_isopnitb = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of isopnitb
1516 0 : c_onitr = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of onitr
1517 0 : c_honitr = 1.26e3_r8 * sqrt_t(i) ! mean molecular speed of honitr
1518 0 : c_terpnit = 0.992e3_r8 * sqrt_t(i) ! mean molecular speed of terpnit
1519 0 : c_nterpooh = 0.957e3_r8 * sqrt_t(i) ! mean molecular speed of nterpooh
1520 0 : c_nc4cho = 1.21e3_r8 * sqrt_t(i) ! mean molecular speed of nc4cho
1521 0 : c_nc4ch2oh = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of nc4ch2oh
1522 0 : c_isopfdn = 9.68e2_r8 * sqrt_t(i) ! mean molecular speed of isopfdn
1523 0 : c_isopfnp = 1.04e3_r8 * sqrt_t(i) ! mean molecular speed of isopfnp
1524 0 : c_isopn2b = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of isopn2
1525 0 : c_isopn1d = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of isopn1d
1526 0 : c_isopn4d = 1.20e3_r8 * sqrt_t(i) ! mean molecular speed of isopn4d
1527 0 : c_inoohd = 1.14e3_r8 * sqrt_t(i) ! mean molecular speed of inoohd
1528 0 : c_inheb = 1.14e3_r8 * sqrt_t(i) ! mean molecular speed of inheb
1529 0 : c_inhed = 1.14e3_r8 * sqrt_t(i) ! mean molecular speed of inhed
1530 0 : c_macrn = 1.19e3_r8 * sqrt_t(i) ! mean molecular speed of macrn
1531 0 : c_isophfp = 1.19e3_r8 * sqrt_t(i) ! mean molecular speed of isophfp
1532 0 : c_iepox = 1.34e3_r8 * sqrt_t(i) ! mean molecular speed of iepox
1533 0 : c_dhpmpal = 1.24e3_r8 * sqrt_t(i) ! mean molecular speed of dhpmpal
1534 0 : c_iche = 1.35e3_r8 * sqrt_t(i) ! mean molecular speed of iche
1535 0 : c_isopfnc = 1.04e3_r8 * sqrt_t(i) ! mean molecular speed of isopfnc
1536 0 : c_isopfdnc = 9.72e2_r8 * sqrt_t(i) ! mean molecular speed of isopfdnc
1537 0 : c_terpnt = 9.92e2_r8 * sqrt_t(i) ! mean molecular speed of terpnt
1538 0 : c_terpnt1 = 9.92e2_r8 * sqrt_t(i) ! mean molecular speed of terpnt1
1539 0 : c_terpnpt = 9.57e2_r8 * sqrt_t(i) ! mean molecular speed of terpnpt
1540 0 : c_terpnpt1 = 9.57e2_r8 * sqrt_t(i) ! mean molecular speed of terpnpt1
1541 0 : c_terpfdn = 8.48e2_r8 * sqrt_t(i) ! mean molecular speed of terpfdn
1542 0 : c_sqtn = 8.64e2_r8 * sqrt_t(i) ! mean molecular speed of sqtn
1543 0 : c_terphfn = 8.93e2_r8 * sqrt_t(i) ! mean molecular speed of terphfn
1544 0 : c_terpdhdp = 9.47e2_r8 * sqrt_t(i) ! mean molecular speed of terpdhdp
1545 0 : c_terpacid = 1.07e3_r8 * sqrt_t(i) ! mean molecular speed of terpacid
1546 : !-------------------------------------------------------------------------
1547 : ! Heterogeneous reaction rates for uptake of a gas on an aerosol:
1548 : ! rxt = sfc / ( (rad_aer/Dg_gas) + (4/(c_gas*gamma_gas)))
1549 : !-------------------------------------------------------------------------
1550 : !-------------------------------------------------------------------------
1551 : ! ... n2o5 -> 2 hno3 (on sulfate, nh4no3, oc2, soa)
1552 : !-------------------------------------------------------------------------
1553 0 : if( usr_N2O5_aer_ndx > 0 ) then
1554 0 : rxt(i,k,usr_N2O5_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_n2o5, gamma_n2o5 )
1555 : end if
1556 0 : if( usr_XNO2NO3_aer_ndx > 0 ) then
1557 0 : rxt(i,k,usr_XNO2NO3_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_n2o5, gamma_n2o5 )
1558 : end if
1559 0 : if( usr_NO2XNO3_aer_ndx > 0 ) then
1560 0 : rxt(i,k,usr_NO2XNO3_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_n2o5, gamma_n2o5 )
1561 : end if
1562 : !-------------------------------------------------------------------------
1563 : ! ... no3 -> hno3 (on sulfate, nh4no3, oc, soa)
1564 : !-------------------------------------------------------------------------
1565 0 : if( usr_NO3_aer_ndx > 0 ) then
1566 0 : rxt(i,k,usr_NO3_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_no3, gamma_no3 )
1567 : end if
1568 0 : if( usr_XNO3_aer_ndx > 0 ) then
1569 0 : rxt(i,k,usr_XNO3_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_no3, gamma_no3 )
1570 : end if
1571 : !-------------------------------------------------------------------------
1572 : ! ... no2 -> 0.5 * (ho+no+hno3) (on sulfate, nh4no3, oc2, soa)
1573 : !-------------------------------------------------------------------------
1574 0 : if( usr_NO2_aer_ndx > 0 ) then
1575 0 : rxt(i,k,usr_NO2_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_no2, gamma_no2 )
1576 : end if
1577 0 : if( usr_XNO2_aer_ndx > 0 ) then
1578 0 : rxt(i,k,usr_XNO2_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_no2, gamma_no2 )
1579 : end if
1580 :
1581 : !-------------------------------------------------------------------------
1582 : ! ... ho2 -> 0.5 * h2o2 (on sulfate, nh4no3, oc2, soa)
1583 : !-------------------------------------------------------------------------
1584 0 : if( usr_HO2_aer_ndx > 0 ) then
1585 0 : rxt(i,k,usr_HO2_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_ho2, gamma_ho2 )
1586 : end if
1587 : !-------------------------------------------------------------------------
1588 : ! ... glyoxal -> soag1 (on sulfate, nh4no3, oc2, soa)
1589 : ! first order uptake, Fuchs and Sutugin, 1971, dCg = 1/4 * gamma * ! A * |v_mol| * Cg * dt
1590 : !-------------------------------------------------------------------------
1591 0 : if( usr_GLYOXAL_aer_ndx > 0 ) then
1592 0 : rxt(i,k,usr_GLYOXAL_aer_ndx) = hetrxtrate_gly( sfc, c_glyoxal, gamma_glyoxal )
1593 : end if
1594 : !-------------------------------------------------------------------------
1595 : ! ... ISOPNITA -> HNO3 (on sulfate, nh4no3, oc2, soa)
1596 : !-------------------------------------------------------------------------
1597 0 : if( usr_ISOPNITA_aer_ndx > 0 ) then
1598 0 : rxt(i,k,usr_ISOPNITA_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopnita, gamma_isopnita )
1599 : end if
1600 : !-------------------------------------------------------------------------
1601 : ! ... ISOPNITB -> HNO3 (on sulfate, nh4no3, oc2, soa)
1602 : !-------------------------------------------------------------------------
1603 0 : if( usr_ISOPNITB_aer_ndx > 0 ) then
1604 0 : rxt(i,k,usr_ISOPNITB_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopnitb, gamma_isopnitb )
1605 : end if
1606 : !-------------------------------------------------------------------------
1607 : ! ... ONITR -> HNO3 (on sulfate, nh4no3, oc2, soa)
1608 : !-------------------------------------------------------------------------
1609 0 : if( usr_ONITR_aer_ndx > 0 ) then
1610 0 : rxt(i,k,usr_ONITR_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_onitr, gamma_onitr )
1611 : end if
1612 : !-------------------------------------------------------------------------
1613 : ! ... HONITR -> HNO3 (on sulfate, nh4no3, oc2, soa)
1614 : !-------------------------------------------------------------------------
1615 0 : if( usr_HONITR_aer_ndx > 0 ) then
1616 0 : rxt(i,k,usr_HONITR_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_honitr, gamma_honitr )
1617 : end if
1618 : !-------------------------------------------------------------------------
1619 : ! ... TERPNIT -> HNO3 (on sulfate, nh4no3, oc2, soa)
1620 : !-------------------------------------------------------------------------
1621 0 : if( usr_TERPNIT_aer_ndx > 0 ) then
1622 0 : rxt(i,k,usr_TERPNIT_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpnit, gamma_terpnit )
1623 : end if
1624 : !-------------------------------------------------------------------------
1625 : ! ... NTERPOOH -> HNO3 (on sulfate, nh4no3, oc2, soa)
1626 : !-------------------------------------------------------------------------
1627 0 : if( usr_NTERPOOH_aer_ndx > 0 ) then
1628 0 : rxt(i,k,usr_NTERPOOH_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_nterpooh, gamma_nterpooh )
1629 : end if
1630 : !-------------------------------------------------------------------------
1631 : ! ... NC4CHO -> HNO3 (on sulfate, nh4no3, oc2, soa)
1632 : !-------------------------------------------------------------------------
1633 0 : if( usr_NC4CHO_aer_ndx > 0 ) then
1634 0 : rxt(i,k,usr_NC4CHO_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_nc4cho, gamma_nc4cho )
1635 : end if
1636 : !-------------------------------------------------------------------------
1637 : ! ... NC4CH2OH -> HNO3 (on sulfate, nh4no3, oc2, soa)
1638 : !-------------------------------------------------------------------------
1639 0 : if( usr_NC4CH2OH_aer_ndx > 0 ) then
1640 0 : rxt(i,k,usr_NC4CH2OH_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_nc4ch2oh, gamma_nc4ch2oh )
1641 : end if
1642 : !-------------------------------------------------------------------------
1643 : ! ... ISOPFDN -> HNO3 (on sulfate, nh4no3, oc2, soa)
1644 : !-------------------------------------------------------------------------
1645 0 : if( usr_ISOPFDN_aer_ndx > 0 ) then
1646 0 : rxt(i,k,usr_ISOPFDN_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopfdn, gamma_isopfdn )
1647 : end if
1648 : !-------------------------------------------------------------------------
1649 : ! ... ISOPFNP -> (on sulfate, nh4no3, oc2, soa)
1650 : !-------------------------------------------------------------------------
1651 0 : if( usr_ISOPFNP_aer_ndx > 0 ) then
1652 0 : rxt(i,k,usr_ISOPFNP_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopfnp, gamma_isopfnp )
1653 : end if
1654 : !-------------------------------------------------------------------------
1655 : ! ... ISOPN2B -> HNO3 (on sulfate, nh4no3, oc2, soa)
1656 : !-------------------------------------------------------------------------
1657 0 : if( usr_ISOPN2B_aer_ndx > 0 ) then
1658 0 : rxt(i,k,usr_ISOPN2B_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopn2b, gamma_isopn2b )
1659 : end if
1660 : !-------------------------------------------------------------------------
1661 : ! ... ISOPN1D -> HNO3 (on sulfate, nh4no3, oc2, soa)
1662 : !-------------------------------------------------------------------------
1663 0 : if( usr_ISOPN1D_aer_ndx > 0 ) then
1664 0 : rxt(i,k,usr_ISOPN1D_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopn1d, gamma_isopn1d )
1665 : end if
1666 : !-------------------------------------------------------------------------
1667 : ! ... ISOPN4D -> HNO3 (on sulfate, nh4no3, oc2, soa)
1668 : !-------------------------------------------------------------------------
1669 0 : if( usr_ISOPN4D_aer_ndx > 0 ) then
1670 0 : rxt(i,k,usr_ISOPN4D_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopn4d, gamma_isopn4d )
1671 : end if
1672 : !-------------------------------------------------------------------------
1673 : ! ... INOOHD -> HNO3 (on sulfate, nh4no3, oc2, soa)
1674 : !-------------------------------------------------------------------------
1675 0 : if( usr_INOOHD_aer_ndx > 0 ) then
1676 0 : rxt(i,k,usr_INOOHD_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_inoohd, gamma_inoohd )
1677 : end if
1678 : !-------------------------------------------------------------------------
1679 : ! ... INHEB -> HNO3 (on sulfate, nh4no3, oc2, soa)
1680 : !-------------------------------------------------------------------------
1681 0 : if( usr_INHEB_aer_ndx > 0 ) then
1682 0 : rxt(i,k,usr_INHEB_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_inheb, gamma_inheb )
1683 : end if
1684 : !-------------------------------------------------------------------------
1685 : ! ... INHED -> HNO3 (on sulfate, nh4no3, oc2, soa)
1686 : !-------------------------------------------------------------------------
1687 0 : if( usr_INHED_aer_ndx > 0 ) then
1688 0 : rxt(i,k,usr_INHED_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_inhed, gamma_inhed )
1689 : end if
1690 : !-------------------------------------------------------------------------
1691 : ! ... MACRN -> HNO3 (on sulfate, nh4no3, oc2, soa)
1692 : !-------------------------------------------------------------------------
1693 0 : if( usr_MACRN_aer_ndx > 0 ) then
1694 0 : rxt(i,k,usr_MACRN_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_macrn, gamma_macrn )
1695 : end if
1696 : !-------------------------------------------------------------------------
1697 : ! ... ISOPHFP -> (on sulfate, nh4no3, oc2, soa)
1698 : !-------------------------------------------------------------------------
1699 0 : if( usr_ISOPHFP_aer_ndx > 0 ) then
1700 0 : rxt(i,k,usr_ISOPHFP_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isophfp, gamma_isophfp )
1701 : end if
1702 : !-------------------------------------------------------------------------
1703 : ! ... IEPOX -> (on sulfate, nh4no3, oc2, soa)
1704 : !-------------------------------------------------------------------------
1705 0 : if( usr_IEPOX_aer_ndx > 0 ) then
1706 0 : rxt(i,k,usr_IEPOX_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_iepox, gamma_iepox )
1707 : end if
1708 : !-------------------------------------------------------------------------
1709 : ! ... DHPMPAL -> (on sulfate, nh4no3, oc2, soa)
1710 : !-------------------------------------------------------------------------
1711 0 : if( usr_DHPMPAL_aer_ndx > 0 ) then
1712 0 : rxt(i,k,usr_DHPMPAL_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_dhpmpal, gamma_dhpmpal )
1713 : end if
1714 : !-------------------------------------------------------------------------
1715 : ! ... ICHE -> (on sulfate, nh4no3, oc2, soa)
1716 : !-------------------------------------------------------------------------
1717 0 : if( usr_ICHE_aer_ndx > 0 ) then
1718 0 : rxt(i,k,usr_ICHE_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_iche, gamma_iche )
1719 : end if
1720 : !-------------------------------------------------------------------------
1721 : ! ... ISOPFNC -> (on sulfate, nh4no3, oc2, soa)
1722 : !-------------------------------------------------------------------------
1723 0 : if( usr_ISOPFNC_aer_ndx > 0 ) then
1724 0 : rxt(i,k,usr_ISOPFNC_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopfnc, gamma_isopfnc )
1725 : end if
1726 : !-------------------------------------------------------------------------
1727 : ! ... ISOPFDNC -> (on sulfate, nh4no3, oc2, soa)
1728 : !-------------------------------------------------------------------------
1729 0 : if( usr_ISOPFDNC_aer_ndx > 0 ) then
1730 0 : rxt(i,k,usr_ISOPFDNC_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_isopfdnc, gamma_isopfdnc )
1731 : end if
1732 : !-------------------------------------------------------------------------
1733 : ! ... TERPNT -> HNO3 (on sulfate, nh4no3, oc2, soa)
1734 : !-------------------------------------------------------------------------
1735 0 : if( usr_TERPNT_aer_ndx > 0 ) then
1736 0 : rxt(i,k,usr_TERPNT_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpnt, gamma_terpnt )
1737 : end if
1738 : !-------------------------------------------------------------------------
1739 : ! ... TERPNT1 -> HNO3 (on sulfate, nh4no3, oc2, soa)
1740 : !-------------------------------------------------------------------------
1741 0 : if( usr_TERPNT1_aer_ndx > 0 ) then
1742 0 : rxt(i,k,usr_TERPNT1_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpnt1, gamma_terpnt1 )
1743 : end if
1744 : !-------------------------------------------------------------------------
1745 : ! ... TERPNPT -> HNO3 (on sulfate, nh4no3, oc2, soa)
1746 : !-------------------------------------------------------------------------
1747 0 : if( usr_TERPNPT_aer_ndx > 0 ) then
1748 0 : rxt(i,k,usr_TERPNPT_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpnpt, gamma_terpnpt )
1749 : end if
1750 : !-------------------------------------------------------------------------
1751 : ! ... TERPNPT1 -> HNO3 (on sulfate, nh4no3, oc2, soa)
1752 : !-------------------------------------------------------------------------
1753 0 : if( usr_TERPNPT1_aer_ndx > 0 ) then
1754 0 : rxt(i,k,usr_TERPNPT1_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpnpt1, gamma_terpnpt1 )
1755 : end if
1756 : !-------------------------------------------------------------------------
1757 : ! ... TERPFDN -> HNO3 (on sulfate, nh4no3, oc2, soa)
1758 : !-------------------------------------------------------------------------
1759 0 : if( usr_TERPFDN_aer_ndx > 0 ) then
1760 0 : rxt(i,k,usr_TERPFDN_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpfdn, gamma_terpfdn )
1761 : end if
1762 : !-------------------------------------------------------------------------
1763 : ! ... SQTN -> (on sulfate, nh4no3, oc2, soa)
1764 : !-------------------------------------------------------------------------
1765 0 : if( usr_SQTN_aer_ndx > 0 ) then
1766 0 : rxt(i,k,usr_SQTN_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_sqtn, gamma_sqtn )
1767 : end if
1768 : !-------------------------------------------------------------------------
1769 : ! ... TERPHFN -> (on sulfate, nh4no3, oc2, soa)
1770 : !-------------------------------------------------------------------------
1771 0 : if( usr_TERPHFN_aer_ndx > 0 ) then
1772 0 : rxt(i,k,usr_TERPHFN_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terphfn, gamma_terphfn )
1773 : end if
1774 : !-------------------------------------------------------------------------
1775 : ! ... TERPDHDP -> (on sulfate, nh4no3, oc2, soa)
1776 : !-------------------------------------------------------------------------
1777 0 : if( usr_TERPDHDP_aer_ndx > 0 ) then
1778 0 : rxt(i,k,usr_TERPDHDP_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpdhdp, gamma_terpdhdp )
1779 : end if
1780 : !-------------------------------------------------------------------------
1781 : ! ... TERPACID -> (on sulfate, nh4no3, oc2, soa)
1782 : !-------------------------------------------------------------------------
1783 0 : if( usr_TERPACID_aer_ndx > 0 ) then
1784 0 : rxt(i,k,usr_TERPACID_aer_ndx) = hetrxtrate( sfc, dm_aer, dg, c_terpacid, gamma_terpacid )
1785 : end if
1786 :
1787 : end do long_loop
1788 : end if
1789 :
1790 : ! LLNL super fast chem reaction rates
1791 :
1792 : !-----------------------------------------------------------------------
1793 : ! ... CO + OH --> CO2 + HO2
1794 : !-----------------------------------------------------------------------
1795 5470632 : if ( usr_oh_co_ndx > 0 ) then
1796 0 : ko(:) = 5.9e-33_r8 * tp(:)**1.4_r8
1797 0 : kinf(:) = 1.1e-12_r8 * (temp(:ncol,k) / 300._r8)**1.3_r8
1798 0 : ko_m(:) = ko(:) * m(:,k)
1799 0 : k0(:) = 1.5e-13_r8 * (temp(:ncol,k) / 300._r8)**0.6_r8
1800 0 : kinf_m(:) = (2.1e+09_r8 * (temp(:ncol,k) / 300._r8)**6.1_r8) / m(:,k)
1801 0 : rxt(:,k,usr_oh_co_ndx) = (ko_m(:)/(1._r8+(ko_m(:)/kinf(:)))) * &
1802 : 0.6_r8**(1._r8/(1._r8+(log10(ko_m(:)/kinf(:)))**2._r8)) + &
1803 : (k0(:)/(1._r8+(k0(:)/kinf_m(:)))) * &
1804 0 : 0.6_r8**(1._r8/(1._r8+(log10(k0(:)/kinf_m(:)))**2._r8))
1805 : endif
1806 : !-----------------------------------------------------------------------
1807 : ! ... NO2 + H2O --> 0.5 HONO + 0.5 HNO3
1808 : !-----------------------------------------------------------------------
1809 5470632 : if ( het_no2_h2o_ndx > 0 ) then
1810 0 : rxt(:,k,het_no2_h2o_ndx) = 4.0e-24_r8
1811 : endif
1812 : !-----------------------------------------------------------------------
1813 : ! ... DMS + OH --> 0.75 SO2 + 0.25 MSA
1814 : !-----------------------------------------------------------------------
1815 5470632 : if ( usr_oh_dms_ndx > 0 ) then
1816 0 : o2(:ncol) = invariants(:ncol,k,inv_o2_ndx)
1817 0 : rxt(:,k,usr_oh_dms_ndx) = 2.000e-10_r8 * exp(5820.0_r8 * tinv(:)) / &
1818 0 : ((2.000e29_r8 / o2(:)) + exp(6280.0_r8 * tinv(:)))
1819 : endif
1820 5470632 : if ( aq_so2_h2o2_ndx > 0 .or. aq_so2_o3_ndx > 0 ) then
1821 0 : lwc(:) = cwat(:ncol,k) * invariants(:ncol,k,inv_m_ndx) * mbar(:ncol,k) /avo !PJC convert kg/kg to g/cm3
1822 : !-----------------------------------------------------------------------
1823 : ! ... SO2 + H2O2 --> S(VI)
1824 : !-----------------------------------------------------------------------
1825 0 : if ( aq_so2_h2o2_ndx > 0 ) then
1826 0 : rxt(:,k,aq_so2_h2o2_ndx) = lwc(:) * 1.0e-03_r8 * avo * &
1827 : K_AQ * &
1828 :
1829 : exp(ER_AQ * ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1830 : HENRY298_SO2 * &
1831 : K298_SO2_HSO3 * &
1832 : HENRY298_H2O2 * &
1833 : exp(((H298_SO2 + H298_SO2_HSO3 + H298_H2O2) / R_CAL) * &
1834 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1835 0 : (R_CONC * temp(:ncol,k))**2.0e+00_r8 / &
1836 :
1837 0 : (1.0e+00_r8 + 13.0e+00_r8 * 10.0e+00_r8**(-pH))
1838 : endif
1839 : !-----------------------------------------------------------------------
1840 : ! ... SO2 + O3 --> S(VI)
1841 : !-----------------------------------------------------------------------
1842 0 : if (aq_so2_o3_ndx >0) then
1843 0 : rxt(:,k,aq_so2_o3_ndx) = lwc(:) * 1.0e-03_r8 * avo * &
1844 : HENRY298_SO2 * exp((H298_SO2 / R_CAL) * &
1845 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1846 : (K0_AQ * exp(ER0_AQ * &
1847 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) + &
1848 : K298_SO2_HSO3 * exp((H298_SO2_HSO3 / R_CAL) * &
1849 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1850 : (K1_AQ * exp(ER1_AQ * &
1851 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) / &
1852 : 10.0e+00_r8**(-pH) + K2_AQ * exp(ER2_AQ * &
1853 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1854 : K298_HSO3_SO3 * exp((H298_HSO3_SO3 / R_CAL) * &
1855 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) / &
1856 : (10.0e+00_r8**(-pH))**2.0e+00_r8) ) * &
1857 : HENRY298_O3 * exp((H298_O3 / R_CAL) * &
1858 : ((1.0e+00_r8 / 298.0e+00_r8) - tinv(:))) * &
1859 0 : (R_CONC * temp(:ncol,k))**2.0e+00_r8
1860 : endif
1861 : endif
1862 :
1863 5529456 : if ( has_d_chem ) then
1864 :
1865 0 : call comp_exp( exp_fac, -600._r8 * tinv, ncol )
1866 0 : rxt(:,k,ean_ndx(1)) = 1.e-31_r8 * tp(:) * exp_fac(:)
1867 0 : rxt(:,k,ean_ndx(2)) = 9.1e-12_r8 * tp(:)**(-1.46_r8)
1868 0 : call comp_exp( exp_fac, -193._r8 * tinv, ncol )
1869 0 : rxt(:,k,ean_ndx(3)) = (4.e-30_r8 * exp_fac(:)) * 0.21_r8
1870 :
1871 0 : rxt(:,k,rpe_ndx(1)) = 4.2e-6_r8 * tp(:)**0.5_r8
1872 0 : rxt(:,k,rpe_ndx(2)) = 6.3e-7_r8 * tp(:)**0.5_r8
1873 0 : rxt(:,k,rpe_ndx(3)) = 2.5e-6_r8 * tp(:)**0.1_r8
1874 0 : rxt(:,k,rpe_ndx(4)) = 2.48e-6_r8 * tp(:)**0.76_r8
1875 0 : rxt(:,k,rpe_ndx(5)) = 1.4e-6_r8 * tp(:)**0.4_r8
1876 :
1877 0 : rxt(:,k,pir_ndx(1)) = 4.e-30_r8 * tp(:)**2.93_r8
1878 0 : rxt(:,k,pir_ndx(2)) = 4.6e-27_r8 * tp(:)**4._r8
1879 :
1880 0 : call comp_exp( exp_fac, -15900._r8 * tinv, ncol )
1881 0 : rxt(:,k,pir_ndx(3)) = (2.5e-2_r8 * tp(:)**5._r8) * exp_fac(:)
1882 0 : rxt(:,k,pir_ndx(4)) = 2.3e-27_r8 * tp(:)**7.5_r8
1883 :
1884 0 : call comp_exp( exp_fac, -10272._r8 * tinv, ncol )
1885 0 : rxt(:,k,pir_ndx(5)) = (2.6e-3_r8 * tp(:)**8.5_r8) * exp_fac(:)
1886 0 : rxt(:,k,pir_ndx(6)) = 3.6e-27_r8 * tp(:)**8.1_r8
1887 :
1888 0 : call comp_exp( exp_fac, -9000._r8 * tinv, ncol )
1889 0 : rxt(:,k,pir_ndx(7)) = (1.5e-1_r8 * tp(:)**9.1_r8) * exp_fac(:)
1890 0 : rxt(:,k,pir_ndx(8)) = 4.6e-28_r8 * tp(:)**14._r8
1891 :
1892 0 : call comp_exp( exp_fac, -6400._r8 * tinv, ncol )
1893 0 : rxt(:,k,pir_ndx(9)) = (1.7e-3_r8 * tp(:)**15._r8) * exp_fac(:)
1894 0 : rxt(:,k,pir_ndx(10)) = 1.35e-28_r8 * tp(:)**2.83_r8
1895 :
1896 0 : rxt(:,k,pir_ndx(11)) = 1.e-27_r8 * (308._r8 * tinv(:))**4.7_r8
1897 0 : rxt(:,k,pir_ndx(12)) = rxt(:,k,pir_ndx(11))
1898 0 : rxt(:,k,pir_ndx(13)) = 1.4e-29_r8 * tp(:)**4._r8
1899 :
1900 0 : call comp_exp( exp_fac, -3872._r8 * tinv, ncol )
1901 0 : rxt(:,k,pir_ndx(14)) = (3.4e-7_r8 * tp(:)**5._r8) * exp_fac(:)
1902 :
1903 0 : rxt(:,k,pir_ndx(15)) = 3.0e-31_r8 * tp(:)**4.3_r8
1904 0 : call comp_exp( exp_fac, -2093._r8 * tinv, ncol )
1905 0 : rxt(:,k,pir_ndx(16)) = (1.5e-8_r8 * tp(:)**4.3_r8) * exp_fac(:)
1906 :
1907 0 : rxt(:,k,edn_ndx(1)) = 3.1e-10_r8 * tp(:)**0.83_r8
1908 0 : call comp_exp( exp_fac, -4990._r8 * tinv, ncol )
1909 0 : rxt(:,k,edn_ndx(2)) = (1.9e-12_r8 * tp(:)**(-1.5_r8)) * exp_fac(:)
1910 :
1911 0 : rxt(:,k,nir_ndx(1)) = 1.05e-12_r8 * tp(:)**2.15_r8
1912 0 : rxt(:,k,nir_ndx(2)) = 2.5e-11_r8 * tp(:)**0.79_r8
1913 0 : rxt(:,k,nir_ndx(3)) = 7.5e-11_r8 * tp(:)**0.79_r8
1914 0 : rxt(:,k,nir_ndx(4)) = rxt(:,k,nir_ndx(1))
1915 0 : rxt(:,k,nir_ndx(5)) = 1.3e-11_r8 * tp(:)**1.64_r8
1916 0 : rxt(:,k,nir_ndx(6)) = 3.3e-11_r8 * tp(:)**2.38_r8
1917 :
1918 0 : call comp_exp( exp_fac, -7300_r8 * tinv, ncol )
1919 0 : rxt(:,k,nir_ndx(7)) = (1.0e-3_r8 * tp(:)) * exp_fac(:)
1920 0 : call comp_exp( exp_fac, -7050_r8 * tinv, ncol )
1921 0 : rxt(:,k,nir_ndx(8)) = (7.2e-4_r8 * tp(:)) * exp_fac(:)
1922 0 : call comp_exp( exp_fac, -6800_r8 * tinv, ncol )
1923 0 : rxt(:,k,nir_ndx(9)) = (6.5e-3_r8 * tp(:)) * exp_fac(:)
1924 0 : call comp_exp( exp_fac, -7600_r8 * tinv, ncol )
1925 0 : rxt(:,k,nir_ndx(10)) = (5.7e-4_r8 * tp(:)) * exp_fac(:)
1926 :
1927 0 : call comp_exp( exp_fac, -7150_r8 * tinv, ncol )
1928 0 : rxt(:,k,nir_ndx(11)) = (1.5e-2_r8 * tp(:)) * exp_fac(:)
1929 :
1930 0 : call comp_exp( exp_fac, -13130_r8 * tinv, ncol )
1931 0 : rxt(:,k,nir_ndx(12)) = (6.0e-3_r8 * tp(:)) * exp_fac(:)
1932 0 : rxt(:,k,nir_ndx(13)) = 5.22e-28_r8 * tp(:)**2.62_r8
1933 :
1934 0 : rxt(:,k,iira_ndx(1)) = 6.0e-8_r8 * tp(:)**.5_r8
1935 0 : do i = 2,niira
1936 0 : rxt(:,k,iira_ndx(i)) = rxt(:,k,iira_ndx(i-1))
1937 : enddo
1938 :
1939 0 : rxt(:,k,iirb_ndx(1)) = 1.25e-25_r8 * tp(:)**4._r8
1940 0 : do i = 2,niirb
1941 0 : rxt(:,k,iirb_ndx(i)) = rxt(:,k,iirb_ndx(i-1))
1942 : enddo
1943 :
1944 0 : call comp_exp( exp_fac, -6600._r8 * tinv, ncol )
1945 0 : rxt(:,k,usr_clm_h2o_m_ndx) = 2.e-8_r8 * exp_fac(:)
1946 :
1947 0 : call comp_exp( exp_fac, -11926._r8 * tinv, ncol )
1948 0 : rxt(:,k,usr_clm_hcl_m_ndx) = tinv(:) * exp_fac(:)
1949 :
1950 : endif
1951 : end do level_loop
1952 :
1953 : !-----------------------------------------------------------------
1954 : ! ... the ionic rates
1955 : !-----------------------------------------------------------------
1956 58824 : if ( has_ion_rxts ) then
1957 0 : level_loop2 : do k = 1,pver
1958 0 : tp(:ncol) = (2._r8*tempi(:ncol,k) + temp(:ncol,k)) / ( 3._r8 * t0 )
1959 0 : tp(:) = max( min( tp(:),20._r8 ),1._r8 )
1960 0 : rxt(:,k,ion1_ndx) = 2.82e-11_r8 + tp(:)*(-7.74e-12_r8 + tp(:)*(1.073e-12_r8 &
1961 0 : + tp(:)*(-5.17e-14_r8 + 9.65e-16_r8*tp(:))))
1962 0 : tp(:ncol) = (.6363_r8*tempi(:ncol,k) + .3637_r8*temp(:ncol,k)) / t0
1963 0 : tp(:) = max( min( tp(:),trlim2 ),1._r8 )
1964 0 : rxt(:,k,ion2_ndx) = 1.533e-12_r8 + tp(:)*(-5.92e-13_r8 + tp(:)*8.6e-14_r8)
1965 0 : tp(:ncol) = 2._r8 * t0 /(tempi(:ncol,k) + temp(:ncol,k))
1966 0 : where( tp(:ncol) < trlim3 )
1967 0 : rxt(:,k,ion3_ndx) = 1.4e-10_r8 * tp(:)**.44_r8
1968 0 : rxt(:,k,ion11_ndx) = 1.e-11_r8 * tp(:)**.23_r8
1969 : elsewhere
1970 : rxt(:,k,ion3_ndx) = 5.2e-11_r8 / tp(:)**.2_r8
1971 : rxt(:,k,ion11_ndx) = 3.6e-12_r8 / tp(:)**.41_r8
1972 : end where
1973 0 : tp(:ncol) = t0 / tempe(:ncol,k)
1974 0 : rxt(:,k,elec1_ndx) = 4.e-7_r8 * tp(:)**.85_r8
1975 0 : rxt(:,k,elec3_ndx) = 1.8e-7_r8 * tp(:)**.39_r8
1976 0 : where( tp(:ncol) < 4._r8 )
1977 : rxt(:,k,elec2_ndx) = 2.7e-7_r8 * tp(:)**.7_r8
1978 : elsewhere
1979 0 : rxt(:,k,elec2_ndx) = 1.6e-7_r8 * tp(:)**.55_r8
1980 : end where
1981 : end do level_loop2
1982 : endif
1983 :
1984 : ! quenching of O+(2P) and O+(2D) by e to produce O+
1985 : ! See TABLE 1 of Roble (1995)
1986 : ! drm 2015-07-27
1987 58824 : if (elec4_ndx > 0 .and. elec5_ndx > 0 .and. elec6_ndx > 0) then
1988 0 : do k=1,pver
1989 0 : tp(:ncol) = sqrt(300._r8 / tempe(:ncol,k))
1990 0 : rxt(:,k,elec4_ndx) = 1.5e-7_r8 * tp(:)
1991 0 : rxt(:,k,elec5_ndx) = 4.0e-8_r8 * tp(:)
1992 0 : rxt(:,k,elec6_ndx) = 6.6e-8_r8 * tp(:)
1993 : end do
1994 : endif
1995 :
1996 : !-----------------------------------------------------------------
1997 : ! ... tropospheric "aerosol" rate constants
1998 : !-----------------------------------------------------------------
1999 58824 : if ( het1_ndx > 0 .AND. (.NOT. usr_N2O5_aer_ndx > 0) ) then
2000 : amas = 4._r8*pi*rm1**3*den/3._r8 ! each mean particle(r=0.1u) mass (g)
2001 0 : do k = 1,pver
2002 : !-------------------------------------------------------------------------
2003 : ! ... estimate humidity effect on aerosols (from Shettle and Fenn, 1979)
2004 : ! xr is a factor of the increase aerosol radii with hum (hum=0., factor=1)
2005 : !-------------------------------------------------------------------------
2006 0 : xr(:) = .999151_r8 + relhum(:ncol,k)*(1.90445_r8 + relhum(:ncol,k)*(-6.35204_r8 + relhum(:ncol,k)*5.32061_r8))
2007 : !-------------------------------------------------------------------------
2008 : ! ... estimate sulfate particles surface area (cm2/cm3) in each grid
2009 : !-------------------------------------------------------------------------
2010 0 : if ( carma_hetchem_feedback ) then
2011 0 : sur(:ncol) = strato_sad(:ncol,k)
2012 : else
2013 : sur(:) = sulfate(:,k)*m(:,k)/avo*wso4 & ! xform mixing ratio to g/cm3
2014 : / amas & ! xform g/cm3 to num particles/cm3
2015 : * fare & ! xform num particles/cm3 to cm2/cm3
2016 0 : * xr(:)*xr(:) ! humidity factor
2017 : endif
2018 : !-----------------------------------------------------------------
2019 : ! ... compute the "aerosol" reaction rates
2020 : !-----------------------------------------------------------------
2021 : ! k = gam * A * velo/4
2022 : !
2023 : ! where velo = sqrt[ 8*bk*T/pi/(w/av) ]
2024 : ! bk = 1.381e-16
2025 : ! av = 6.02e23
2026 : ! w = 108 (n2o5) HO2(33) CH2O (30) NH3(15)
2027 : !
2028 : ! so that velo = 1.40e3*sqrt(T) (n2o5) gama=0.1
2029 : ! so that velo = 2.53e3*sqrt(T) (HO2) gama>0.2
2030 : ! so that velo = 2.65e3*sqrt(T) (CH2O) gama>0.022
2031 : ! so that velo = 3.75e3*sqrt(T) (NH3) gama=0.4
2032 : !--------------------------------------------------------
2033 : !-----------------------------------------------------------------
2034 : ! ... use this n2o5 -> 2*hno3 only in troposphere
2035 : !-----------------------------------------------------------------
2036 0 : rxt(:,k,het1_ndx) = rxt(:,k,het1_ndx) &
2037 0 : +.25_r8 * gam1 * sur(:) * 1.40e3_r8 * sqrt( temp(:ncol,k) )
2038 : end do
2039 : end if
2040 :
2041 : !lke++
2042 : !-----------------------------------------------------------------
2043 : ! ... CO tags
2044 : !-----------------------------------------------------------------
2045 58824 : if( usr_CO_OH_b_ndx > 0 .and. usr_CO_OH_ndx < 0 ) then
2046 0 : usr_CO_OH_ndx = usr_CO_OH_b_ndx
2047 : end if
2048 58824 : if( usr_CO_OH_ndx > 0 ) then
2049 0 : if( usr_COhc_OH_ndx > 0 ) then
2050 0 : rxt(:ncol,:,usr_COhc_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2051 : end if
2052 0 : if( usr_COme_OH_ndx > 0 ) then
2053 0 : rxt(:ncol,:,usr_COme_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2054 : end if
2055 0 : if( usr_CO01_OH_ndx > 0 ) then
2056 0 : rxt(:ncol,:,usr_CO01_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2057 : end if
2058 0 : if( usr_CO02_OH_ndx > 0 ) then
2059 0 : rxt(:ncol,:,usr_CO02_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2060 : end if
2061 0 : if( usr_CO03_OH_ndx > 0 ) then
2062 0 : rxt(:ncol,:,usr_CO03_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2063 : end if
2064 0 : if( usr_CO04_OH_ndx > 0 ) then
2065 0 : rxt(:ncol,:,usr_CO04_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2066 : end if
2067 0 : if( usr_CO05_OH_ndx > 0 ) then
2068 0 : rxt(:ncol,:,usr_CO05_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2069 : end if
2070 0 : if( usr_CO06_OH_ndx > 0 ) then
2071 0 : rxt(:ncol,:,usr_CO06_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2072 : end if
2073 0 : if( usr_CO07_OH_ndx > 0 ) then
2074 0 : rxt(:ncol,:,usr_CO07_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2075 : end if
2076 0 : if( usr_CO08_OH_ndx > 0 ) then
2077 0 : rxt(:ncol,:,usr_CO08_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2078 : end if
2079 0 : if( usr_CO09_OH_ndx > 0 ) then
2080 0 : rxt(:ncol,:,usr_CO09_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2081 : end if
2082 0 : if( usr_CO10_OH_ndx > 0 ) then
2083 0 : rxt(:ncol,:,usr_CO10_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2084 : end if
2085 0 : if( usr_CO11_OH_ndx > 0 ) then
2086 0 : rxt(:ncol,:,usr_CO11_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2087 : end if
2088 0 : if( usr_CO12_OH_ndx > 0 ) then
2089 0 : rxt(:ncol,:,usr_CO12_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2090 : end if
2091 0 : if( usr_CO13_OH_ndx > 0 ) then
2092 0 : rxt(:ncol,:,usr_CO13_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2093 : end if
2094 0 : if( usr_CO14_OH_ndx > 0 ) then
2095 0 : rxt(:ncol,:,usr_CO14_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2096 : end if
2097 0 : if( usr_CO15_OH_ndx > 0 ) then
2098 0 : rxt(:ncol,:,usr_CO15_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2099 : end if
2100 0 : if( usr_CO16_OH_ndx > 0 ) then
2101 0 : rxt(:ncol,:,usr_CO16_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2102 : end if
2103 0 : if( usr_CO17_OH_ndx > 0 ) then
2104 0 : rxt(:ncol,:,usr_CO17_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2105 : end if
2106 0 : if( usr_CO18_OH_ndx > 0 ) then
2107 0 : rxt(:ncol,:,usr_CO18_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2108 : end if
2109 0 : if( usr_CO19_OH_ndx > 0 ) then
2110 0 : rxt(:ncol,:,usr_CO19_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2111 : end if
2112 0 : if( usr_CO20_OH_ndx > 0 ) then
2113 0 : rxt(:ncol,:,usr_CO20_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2114 : end if
2115 0 : if( usr_CO21_OH_ndx > 0 ) then
2116 0 : rxt(:ncol,:,usr_CO21_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2117 : end if
2118 0 : if( usr_CO22_OH_ndx > 0 ) then
2119 0 : rxt(:ncol,:,usr_CO22_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2120 : end if
2121 0 : if( usr_CO23_OH_ndx > 0 ) then
2122 0 : rxt(:ncol,:,usr_CO23_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2123 : end if
2124 0 : if( usr_CO24_OH_ndx > 0 ) then
2125 0 : rxt(:ncol,:,usr_CO24_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2126 : end if
2127 0 : if( usr_CO25_OH_ndx > 0 ) then
2128 0 : rxt(:ncol,:,usr_CO25_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2129 : end if
2130 0 : if( usr_CO26_OH_ndx > 0 ) then
2131 0 : rxt(:ncol,:,usr_CO26_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2132 : end if
2133 0 : if( usr_CO27_OH_ndx > 0 ) then
2134 0 : rxt(:ncol,:,usr_CO27_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2135 : end if
2136 0 : if( usr_CO28_OH_ndx > 0 ) then
2137 0 : rxt(:ncol,:,usr_CO28_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2138 : end if
2139 0 : if( usr_CO29_OH_ndx > 0 ) then
2140 0 : rxt(:ncol,:,usr_CO29_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2141 : end if
2142 0 : if( usr_CO30_OH_ndx > 0 ) then
2143 0 : rxt(:ncol,:,usr_CO30_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2144 : end if
2145 0 : if( usr_CO31_OH_ndx > 0 ) then
2146 0 : rxt(:ncol,:,usr_CO31_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2147 : end if
2148 0 : if( usr_CO32_OH_ndx > 0 ) then
2149 0 : rxt(:ncol,:,usr_CO32_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2150 : end if
2151 0 : if( usr_CO33_OH_ndx > 0 ) then
2152 0 : rxt(:ncol,:,usr_CO33_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2153 : end if
2154 0 : if( usr_CO34_OH_ndx > 0 ) then
2155 0 : rxt(:ncol,:,usr_CO34_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2156 : end if
2157 0 : if( usr_CO35_OH_ndx > 0 ) then
2158 0 : rxt(:ncol,:,usr_CO35_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2159 : end if
2160 0 : if( usr_CO36_OH_ndx > 0 ) then
2161 0 : rxt(:ncol,:,usr_CO36_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2162 : end if
2163 0 : if( usr_CO37_OH_ndx > 0 ) then
2164 0 : rxt(:ncol,:,usr_CO37_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2165 : end if
2166 0 : if( usr_CO38_OH_ndx > 0 ) then
2167 0 : rxt(:ncol,:,usr_CO38_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2168 : end if
2169 0 : if( usr_CO39_OH_ndx > 0 ) then
2170 0 : rxt(:ncol,:,usr_CO39_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2171 : end if
2172 0 : if( usr_CO40_OH_ndx > 0 ) then
2173 0 : rxt(:ncol,:,usr_CO40_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2174 : end if
2175 0 : if( usr_CO41_OH_ndx > 0 ) then
2176 0 : rxt(:ncol,:,usr_CO41_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2177 : end if
2178 0 : if( usr_CO42_OH_ndx > 0 ) then
2179 0 : rxt(:ncol,:,usr_CO42_OH_ndx) = rxt(:ncol,:,usr_CO_OH_ndx)
2180 : end if
2181 : end if
2182 : !lke--
2183 : !
2184 : ! jfl : additional BAM removal reactions. Zero out below the tropopause
2185 : !
2186 1352952 : do l=1,num_strat_tau
2187 : !
2188 1352952 : if ( usr_strat_tau_ndx(l) > 0 ) then
2189 0 : do i=1,ncol
2190 0 : rxt(i,tropchemlev(i)+1:pver,usr_strat_tau_ndx(l)) = 0._r8
2191 : end do
2192 : end if
2193 : !
2194 : end do
2195 : !
2196 :
2197 58824 : deallocate( sfc_array, dm_array )
2198 :
2199 117648 : end subroutine usrrxt
2200 :
2201 0 : subroutine usrrxt_hrates( rxt, tempn, tempi, tempe, &
2202 0 : h2ovmr, m, ncol, kbot )
2203 : !-----------------------------------------------------------------
2204 : ! ... set the user specified reaction rates for heating
2205 : !-----------------------------------------------------------------
2206 :
2207 58824 : use shr_kind_mod, only : r8 => shr_kind_r8
2208 : use chem_mods, only : rxntot
2209 : use ppgrid, only : pver, pcols
2210 :
2211 : implicit none
2212 :
2213 : !-----------------------------------------------------------------
2214 : ! ... dummy arguments
2215 : !-----------------------------------------------------------------
2216 : integer, intent(in) :: ncol ! number columns in chunk
2217 : integer, intent(in) :: kbot ! heating levels
2218 : real(r8), intent(in) :: tempn(pcols,pver) ! neutral temperature (K)
2219 : real(r8), intent(in) :: tempi(pcols,pver) ! ion temperature (K)
2220 : real(r8), intent(in) :: tempe(pcols,pver) ! electron temperature (K)
2221 : real(r8), intent(in) :: m(ncol,pver) ! total atm density (1/cm^3)
2222 : real(r8), intent(in) :: h2ovmr(ncol,pver) ! water vapor (vmr)
2223 : real(r8), intent(inout) :: rxt(ncol,pver,rxntot) ! gas phase rates
2224 :
2225 : !-----------------------------------------------------------------
2226 : ! ... local variables
2227 : !-----------------------------------------------------------------
2228 :
2229 : integer :: k
2230 : real(r8), dimension(ncol) :: &
2231 0 : tp, &
2232 0 : tinv, &
2233 0 : ko, &
2234 0 : kinf, &
2235 0 : fc
2236 :
2237 : !-----------------------------------------------------------------
2238 : ! ... o + o2 + m --> o3 + m
2239 : !-----------------------------------------------------------------
2240 0 : do k = 1,kbot
2241 0 : tinv(:ncol) = 1._r8 / tempn(:ncol,k)
2242 0 : tp(:) = 300._r8 * tinv(:)
2243 0 : rxt(:,k,usr_O_O2_ndx) = 6.e-34_r8 * tp(:)**2.4_r8
2244 :
2245 : !-----------------------------------------------------------------
2246 : ! ... o + o + m -> o2 + m
2247 : !-----------------------------------------------------------------
2248 0 : rxt(:,k,usr_O_O_ndx) = 2.76e-34_r8 * exp( 720.0_r8*tinv(:) )
2249 :
2250 : !-----------------------------------------------------------------
2251 : ! ... ho2 + ho2 --> h2o2
2252 : ! Note: this rate involves the water vapor number density
2253 : !-----------------------------------------------------------------
2254 0 : ko(:) = 3.0e-13_r8 * exp( 460._r8*tinv(:) )
2255 0 : kinf(:) = 2.1e-33_r8 * m(:,k) * exp( 920._r8*tinv(:) )
2256 0 : fc(:) = 1._r8 + 1.4e-21_r8 * m(:,k) * h2ovmr(:,k) * exp( 2200._r8*tinv(:) )
2257 0 : rxt(:,k,usr_HO2_HO2_ndx) = (ko(:) + kinf(:)) * fc(:)
2258 :
2259 : end do
2260 :
2261 : !-----------------------------------------------------------------
2262 : ! ... the ionic rates
2263 : !-----------------------------------------------------------------
2264 0 : if ( has_ion_rxts ) then
2265 0 : level_loop2 : do k = 1,kbot
2266 0 : tp(:ncol) = (2._r8*tempi(:ncol,k) + tempn(:ncol,k)) / ( 3._r8 * t0 )
2267 0 : tp(:) = max( min( tp(:),20._r8 ),1._r8 )
2268 0 : rxt(:,k,ion1_ndx) = 2.82e-11_r8 + tp(:)*(-7.74e-12_r8 + tp(:)*(1.073e-12_r8 &
2269 0 : + tp(:)*(-5.17e-14_r8 + 9.65e-16_r8*tp(:))))
2270 0 : tp(:ncol) = (.6363_r8*tempi(:ncol,k) + .3637_r8*tempn(:ncol,k)) / t0
2271 0 : tp(:) = max( min( tp(:),trlim2 ),1._r8 )
2272 0 : rxt(:,k,ion2_ndx) = 1.533e-12_r8 + tp(:)*(-5.92e-13_r8 + tp(:)*8.6e-14_r8)
2273 0 : tp(:ncol) = 2._r8 * t0 /(tempi(:ncol,k) + tempn(:ncol,k))
2274 0 : where( tp(:ncol) < trlim3 )
2275 : rxt(:,k,ion3_ndx) = 1.4e-10_r8 * tp(:)**.44_r8
2276 : elsewhere
2277 0 : rxt(:,k,ion3_ndx) = 5.2e-11_r8 / tp(:)**.2_r8
2278 : endwhere
2279 0 : tp(:ncol) = t0 / tempe(:ncol,k)
2280 0 : rxt(:,k,elec1_ndx) = 4.e-7_r8 * tp(:)**.85_r8
2281 0 : rxt(:,k,elec3_ndx) = 1.8e-7_r8 * tp(:)**.39_r8
2282 0 : where( tp(:ncol) < 4._r8 )
2283 : rxt(:,k,elec2_ndx) = 2.7e-7_r8 * tp(:)**.7_r8
2284 : elsewhere
2285 0 : rxt(:,k,elec2_ndx) = 1.6e-7_r8 * tp(:)**.55_r8
2286 : endwhere
2287 : end do level_loop2
2288 : endif
2289 0 : end subroutine usrrxt_hrates
2290 :
2291 : !-------------------------------------------------------------------------
2292 : !-------------------------------------------------------------------------
2293 32823792 : subroutine comp_exp( x, y, n )
2294 :
2295 : implicit none
2296 :
2297 : real(r8), intent(out) :: x(:)
2298 : real(r8), intent(in) :: y(:)
2299 : integer, intent(in) :: n
2300 :
2301 : #ifdef IBM
2302 : call vexp( x, y, n )
2303 : #else
2304 548080992 : x(:n) = exp( y(:n) )
2305 : #endif
2306 :
2307 32823792 : end subroutine comp_exp
2308 :
2309 : !-------------------------------------------------------------------------
2310 : ! Heterogeneous reaction rates for uptake of a gas on an aerosol:
2311 : !-------------------------------------------------------------------------
2312 0 : function hetrxtrate( sfc, dm_aer, dg_gas, c_gas, gamma_gas ) result(rate)
2313 :
2314 : real(r8), intent(in) :: sfc(:)
2315 : real(r8), intent(in) :: dm_aer(:)
2316 : real(r8), intent(in) :: dg_gas
2317 : real(r8), intent(in) :: c_gas
2318 : real(r8), intent(in) :: gamma_gas
2319 : real(r8) :: rate
2320 :
2321 0 : real(r8),allocatable :: rxt(:)
2322 : integer :: n, i
2323 :
2324 0 : n = size(sfc)
2325 :
2326 0 : allocate(rxt(n))
2327 0 : do i=1,n
2328 0 : rxt(i) = sfc(i) / (0.5_r8*dm_aer(i)/dg_gas + (4._r8/(c_gas*gamma_gas)))
2329 : enddo
2330 :
2331 0 : rate = sum(rxt)
2332 :
2333 0 : deallocate(rxt)
2334 :
2335 0 : endfunction hetrxtrate
2336 :
2337 : !-------------------------------------------------------------------------
2338 : ! Heterogeneous reaction rates for uptake of a glyoxal gas on an aerosol:
2339 : !-------------------------------------------------------------------------
2340 0 : function hetrxtrate_gly( sfc, c_gas, gamma_gas ) result(rate)
2341 :
2342 : real(r8), intent(in) :: sfc(:)
2343 : real(r8), intent(in) :: c_gas
2344 : real(r8), intent(in) :: gamma_gas
2345 : real(r8) :: rate
2346 :
2347 0 : real(r8),allocatable :: rxt(:)
2348 : integer :: n, i
2349 :
2350 0 : n = size(sfc)
2351 :
2352 0 : allocate(rxt(n))
2353 0 : do i=1,n
2354 0 : rxt(i) = 0.25_r8 * c_gas * sfc(i) * gamma_gas
2355 : enddo
2356 :
2357 0 : rate = sum(rxt)
2358 :
2359 0 : deallocate(rxt)
2360 :
2361 0 : endfunction hetrxtrate_gly
2362 :
2363 :
2364 : end module mo_usrrxt
|