module_cbmz.F
References to this file elsewhere.
1 !**********************************************************************************
2 ! This computer software was prepared by Battelle Memorial Institute, hereinafter
3 ! the Contractor, under Contract No. DE-AC05-76RL0 1830 with the Department of
4 ! Energy (DOE). NEITHER THE GOVERNMENT NOR THE CONTRACTOR MAKES ANY WARRANTY,
5 ! EXPRESS OR IMPLIED, OR ASSUMES ANY LIABILITY FOR THE USE OF THIS SOFTWARE.
6 !
7 ! Chemistry Option: CBMZ (Carbon Bond Mechanism IV - Zaveri)
8 ! * Primary investigator: Rahul A. Zaveri
9 ! * Co-investigator: Richard C. Easter, William I. Gustafson Jr.
10 ! Last update: September 2005
11 !
12 ! Contacts:
13 ! Rahul A. Zaveri, PhD Jerome D. Fast, PhD
14 ! Senior Research Scientist Staff Scientist
15 ! Pacific Northwest National Laboratory Pacific Northwest National Laboratory
16 ! P.O. Box 999, MSIN K9-30 P.O. Box 999, MSIN K9-30
17 ! Richland, WA 99352 Richland, WA, 99352
18 ! Phone: (509) 372-6159 Phone: (509) 372-6116
19 ! Email: Rahul.Zaveri@pnl.gov Email: Jerome.Fast@pnl.gov
20 !
21 ! Please report any bugs or problems to Rahul Zaveri, the primary author of the
22 ! code, or Jerome Fast, the WRF-chem implementation team leader
23 !
24 !Terms of Use:
25 ! 1) CBMZ and its sub-modules may not be included in any commerical package,
26 ! or used for any commercial applications without the primary author's
27 ! prior consent.
28 ! 2) The CBMZ source code is provided to the WRF modeling community; however,
29 ! no portion of CBMZ can be used separately or in another code without the
30 ! primary author's prior consent.
31 ! 3) The CBMZ source code may be used for research, educational, and non-profit
32 ! purposes only. Any other usage must be first approved by the primary author.
33 ! 4) Publications resulting from the usage of CBMZ must use one or more of the
34 ! references below (depending on the application) for proper acknowledgment.
35 !
36 ! References:
37 ! 1) Zaveri R.A., and L.K. Peters (1999), A new lumped structure photochemical
38 ! mechanism for large-scale applications, J. Geophys. Res., 104, 30387-30415.
39 ! 2) Fast, J.D., W.I. Gustafson Jr., R.C. Easter, R.A. Zaveri, J.C. Barnard, E.G.
40 ! Chapman, G.A. Grell, and S.E. Peckham (2005), Evolution of ozone, particulates,
41 ! and aerosol direct radiative forcing in the vicinity of Houston using a fully-
42 ! coupled meteorology-chemistry-aerosol model, Submitted to J. Geophys. Res.
43 !
44 ! Contact Jerome Fast for updates on the status of manuscripts under review.
45 !
46 ! Additional information:
47 ! 1) www.pnl.gov/atmos_sciences/raz
48 ! 2) www.pnl.gov/atmos_sciences/Jdf/wrfchem.html
49 !
50 ! Support:
51 ! Funding for developing and evaluating CBMZ was provided by the U.S. Department
52 ! of Energy under the auspices of Atmospheric Science Program of the Office of
53 ! Biological and Environmental Research the the PNNL Laboratory Research and
54 ! Directed Research and Development program.
55 !**********************************************************************************
56 module module_cbmz
57
58
59
60 use module_peg_util
61
62 contains
63
64
65 !***********************************************************************
66 ! < 1.> subr cbmz_driver
67 !
68 ! purpose: serves as an interface between subr. gas_chemistry and
69 ! the actual solver subr such as lsodes, rodas, etc.
70 !
71 ! grid : fixed i,j,k (box-model)
72 !
73 ! author : Rahul A. Zaveri
74 ! date : November 1998
75 !
76 !-----------------------------------------------------------------------
77
78 subroutine cbmz_driver( &
79 id, ktau, dtstep, ktauc, dtstepc, config_flags, &
80 gmt, julday, t_phy, moist, p8w, t8w, &
81 p_phy, chem, rho_phy, dz8w, z, z_at_w, vdrog3, &
82 ph_o31d, ph_o33p, ph_no2, ph_no3o2, ph_no3o, ph_hno2, &
83 ph_hno3, ph_hno4, ph_h2o2, ph_ch2or, ph_ch2om, &
84 ph_ch3o2h, ph_n2o5, &
85 ids,ide, jds,jde, kds,kde, &
86 ims,ime, jms,jme, kms,kme, &
87 its,ite, jts,jte, kts,kte )
88
89 USE module_configure, only: grid_config_rec_type, num_moist, num_chem, &
90 p_qv, p_so2, p_ho2, p_so4aj, p_corn, p_hcl, p_mtf
91 USE module_data_sorgam, only: ldrog
92 USE module_data_cbmz
93 IMPLICIT NONE
94
95
96 !-----------------------------------------------------------------------
97 ! subr arguments
98
99 INTEGER, INTENT(IN ) :: id, julday, &
100 ids,ide, jds,jde, kds,kde, &
101 ims,ime, jms,jme, kms,kme, &
102 its,ite, jts,jte, kts,kte
103
104 INTEGER, INTENT(IN ) :: ktau, ktauc
105
106 REAL, INTENT(IN ) :: dtstep, dtstepc, gmt
107 !
108 ! advected moisture variables
109 !
110 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_moist ), &
111 INTENT(IN ) :: moist
112 !
113 ! advected chemical tracers
114 !
115 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
116 INTENT(INOUT ) :: chem
117 !
118 ! arrays that hold photolysis rates
119 !
120 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
121 INTENT(INOUT ) :: &
122 ph_o31d, ph_o33p, ph_no2, ph_no3o2, ph_no3o, ph_hno2, &
123 ph_hno3, ph_hno4, ph_h2o2, ph_ch2or, ph_ch2om, &
124 ph_ch3o2h, ph_n2o5
125 !
126 ! on input from met model
127 !
128 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , &
129 INTENT(IN ) :: &
130 t_phy, & ! temperature
131 rho_phy, & ! air density (kg/m3)
132 p_phy, & ! NOT USED
133 z, z_at_w, & ! NOT USED
134 dz8w, & ! NOT USED
135 t8w, p8w ! NOT USED
136 !
137 ! for interaction with aerosols (really is output)
138 !
139 REAL, DIMENSION( ims:ime , kms:kme-0 , jms:jme , ldrog ) , &
140 INTENT(INOUT ) :: &
141 vdrog3 ! NOT USED
142
143 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
144
145
146 !-----------------------------------------------------------------------
147
148
149 ! local variables
150 integer :: idum, iok
151 integer :: iregime
152 integer :: igaschem_allowed_m1, igaschem_allowed_m2, &
153 igaschem_allowed_m3, igaschem_allowed_m4
154 integer :: igas_solver, iregime_forced
155 integer :: i_boxtest_units_convert
156 integer :: i_print_gasode_stats
157 integer :: i_force_dump, mode_force_dump
158 integer :: it, jt, kt
159 integer :: jsolver
160 integer :: lunerr, lunout, levdbg_err, levdbg_info
161 integer :: mgaschem
162
163
164 real :: abs_error, rel_error, trun
165 real :: tchem, dtchem
166 real :: tstart, tstop
167 real :: airdenbox, pressbox, tempbox
168 real :: cair_mlc
169 real :: h2o, ch4, oxygen, nitrogen, hydrogen
170 real :: cboxnew(ngas_z), cboxold(ngas_z)
171 real :: Aperox(nperox,nperox), Bperox(nperox,nperox)
172 real :: rk_param(nperox), rk_photo(nphoto)
173 real :: rk_m1(nrxn_m1), rk_m2(nrxn_m2), rk_m3(nrxn_m3), rk_m4(nrxn_m4)
174
175 integer, dimension(2,6), save :: inforodas=0
176 integer, dimension(6), save :: iodestatus_count=0, ioderegime_count=0
177
178 #ifdef CHEM_DBG_I
179 !rcetestb diagnostics --------------------------------------------------
180 print 93010, ' '
181 print 93010, 'rcetestb diagnostics from cbmz_driver'
182 print 93010, 'id, chem_opt, ktau, ktauc, julday ', &
183 id, config_flags%chem_opt, ktau, ktauc, julday
184 print 93020, 'dtstep, dtstepc, gmt ', &
185 dtstep, dtstepc, gmt
186 print 93010, 'ids/e, j, k', ids, ide, jds, jde, kds, kde
187 print 93010, 'ims/e, j, k', ims, ime, jms, jme, kms, kme
188 print 93010, 'its/e, j, k', its, ite, jts, jte, kts, kte
189 print 93010, 'num_moist, p_qv ', num_moist, p_qv
190 print 93010, 'num_chem, p_so2, p_ho2 ', num_chem, p_so2, p_ho2
191 print 93010, 'p_so4aj, p_corn, p_hcl, p_mtf', p_so4aj, p_corn, p_hcl, p_mtf
192 93010 format( a, 8(1x,i6) )
193 93020 format( a, 8(1p,e14.6) )
194 !rcetestb diagnostics --------------------------------------------------
195 #endif
196
197
198 ! set some control variables to their "standard for wrf-chem" values
199 igas_solver = 1
200 iregime_forced = -1
201 mgaschem = +1
202 i_boxtest_units_convert = 0
203
204 i_print_gasode_stats = 1
205 mode_force_dump = 0
206 lunerr = -1
207 lunout = -1
208 levdbg_err = 0
209 levdbg_info = 15
210
211 abs_error = 1.0e1 ! solver absolute tolerance (molecules/cm3)
212 rel_error = 1.0e-3 ! solver relative tolerance
213
214 ! set some control variables to non-standard values for testing
215 ! force dumps for center column, every 3rd level
216 ! mode_force_dump = +77
217 ! force dumps for center column, 1st level
218 ! mode_force_dump = +7
219 ! do levdbg_info output always
220 levdbg_info = 0
221
222
223 ! following call is for boxwrf testing only
224 ! it must be commented out for actual wrf applications
225 ! call boxtest_get_extra_args( &
226 ! igas_solver, iregime_forced, &
227 ! i_boxtest_units_convert, lunerr, lunout, &
228 ! abs_error, rel_error, trun )
229
230
231 ! currently nothing is done with vdrog3
232 ! vdrog3(its:ite,kts:kte,jts:jte,:) = 0.0 !This is already set to zero in chem_driver.
233
234
235 ! determine which regimes are allowed
236 ! based on which gas species are "active"
237 call set_gaschem_allowed_regimes( lunerr, &
238 igaschem_allowed_m1, igaschem_allowed_m2, &
239 igaschem_allowed_m3, igaschem_allowed_m4 )
240
241 !
242 ! main loop -- do gas chemistry at each i,j,k
243 !
244 do 2900 jt = jts, jte
245 do 2900 kt = kts, kte-1
246 do 2900 it = its, ite
247
248 trun = dtstep*(ktau-1) ! run time in s
249 tchem = gmt*3600.0 + dtstep*(ktau-1)
250 tchem = mod( tchem, 86400.0 ) ! time from 00 UTC in s
251 dtchem = dtstepc
252 tstart = tchem ! s
253 tstop = tstart + dtchem ! s
254
255 ! skip integration for very small dtchem
256 if ((tstop-tstart) .le. 1.0e-5) goto 2900
257
258
259 ! initial species mapping from host array
260 call mapgas_tofrom_host( 0, &
261 i_boxtest_units_convert, &
262 it,jt,kt, ims,ime, jms,jme, kms,kme, &
263 num_moist, num_chem, moist, chem, &
264 t_phy, p_phy, rho_phy, &
265 cboxold, tempbox, pressbox, airdenbox, &
266 cair_mlc, &
267 h2o, ch4, oxygen, nitrogen, hydrogen )
268 cboxnew(:) = cboxold(:)
269
270 ! determine regime
271 call selectgasregime( iregime, iregime_forced, cboxold, &
272 igaschem_allowed_m1, igaschem_allowed_m2, &
273 igaschem_allowed_m3, igaschem_allowed_m4 )
274 idum = iregime
275 if ((idum .lt. 1) .or. (idum .ge. 6)) idum = 6
276 ioderegime_count(idum) = ioderegime_count(idum) + 1
277 iodestatus_count(6) = iodestatus_count(6) + 1
278
279 ! compute rate constants
280 ! transfer/map incoming photolysis rate contants to local array
281 call gasphotoconstants( rk_photo, &
282 i_boxtest_units_convert, &
283 it,jt,kt, ims,ime, jms,jme, kms,kme, &
284 ph_o31d, ph_o33p, ph_no2, ph_no3o2, ph_no3o, ph_hno2, &
285 ph_hno3, ph_hno4, ph_h2o2, ph_ch2or, ph_ch2om, &
286 ph_ch3o2h, ph_n2o5 )
287 ! loads Aperox and Bperox
288 call loadperoxyparameters( Aperox, Bperox )
289 ! calculate parameterized rate constants
290 call peroxyrateconstants( tempbox, cboxold, &
291 Aperox, Bperox, rk_param )
292 ! calculate thermal rate constants
293 call gasrateconstants( iregime, tempbox, cair_mlc, &
294 rk_photo, rk_param, rk_m1, rk_m2, rk_m3, rk_m4 )
295
296 ! mode_force_dump selects a detailed dump of gaschem at either
297 ! first ijk grid, first ij column, all ijk, or no ijk
298 i_force_dump = 0
299 if (mode_force_dump .eq. 1) then
300 if ((it.eq.its) .and. (jt.eq.jts) &
301 .and. (kt.eq.kts)) i_force_dump = 1
302 else if (mode_force_dump .eq. 10) then
303 if ((it.eq.its) .and. (jt.eq.jts)) i_force_dump = 1
304 else if (mode_force_dump .eq. 100) then
305 i_force_dump = 1
306 else if (mode_force_dump .eq. 7) then
307 if ( (it .eq. (its+ite)/2) .and. &
308 (jt .eq. (jts+jte)/2) .and. &
309 (kt .eq. kts) ) i_force_dump = 1
310 else if (mode_force_dump .eq. 77) then
311 if ( (it .eq. (its+ite)/2) .and. &
312 (jt .eq. (jts+jte)/2) .and. &
313 (mod(kt-kts,3) .eq. 0) ) i_force_dump = 1
314 end if
315
316
317 ! rodas
318 iok = 0
319 jsolver = 0
320 if (igas_solver .eq. 1) then
321 jsolver = 1
322 call gasodesolver_rodas( tstart, tstop, iok, &
323 it, jt, kt, iregime, &
324 mgaschem, lunerr, lunout, levdbg_err, levdbg_info, &
325 i_force_dump, inforodas, iodestatus_count, &
326 abs_error, rel_error, trun, &
327 tempbox, pressbox, airdenbox, cboxnew, cboxold, &
328 rk_m1, rk_m2, rk_m3, rk_m4 )
329 endif
330
331 ! lsodes
332 if (igas_solver.eq.2 .or. iok.le.0) then
333 jsolver = 2
334 call gasodesolver_lsodes( tstart, tstop, iok, &
335 it, jt, kt, iregime, &
336 mgaschem, lunerr, lunout, levdbg_err, levdbg_info, &
337 i_force_dump, iodestatus_count, &
338 abs_error, rel_error, trun, &
339 tempbox, pressbox, airdenbox, cboxnew, cboxold, &
340 rk_m1, rk_m2, rk_m3, rk_m4 )
341 endif
342
343 ! final species mapping back to host array -- only when iok > 0
344 if (iok .gt. 0) then
345 call mapgas_tofrom_host( 1, &
346 i_boxtest_units_convert, &
347 it,jt,kt, ims,ime, jms,jme, kms,kme, &
348 num_moist, num_chem, moist, chem, &
349 t_phy, p_phy, rho_phy, &
350 cboxnew, tempbox, pressbox, airdenbox, &
351 cair_mlc, &
352 h2o, ch4, oxygen, nitrogen, hydrogen )
353 end if
354
355 ! following call is for boxwrf testing only
356 ! it must be commented out for actual wrf applications
357 ! call boxtest_set_extra_args( iregime, it, jt, kt )
358
359 2900 continue
360
361 if (i_print_gasode_stats .gt. 0) &
362 call print_gasode_stats( lunout, levdbg_info, &
363 inforodas, iodestatus_count, ioderegime_count )
364 return
365 end subroutine cbmz_driver
366
367
368
369 !***********************************************************************
370 ! < xx.> subr print_gasode_stats
371 !
372 ! purpose: writes some statistics on ode solver performance to unit lunout
373 !
374 !-----------------------------------------------------------------------
375
376 subroutine print_gasode_stats( lunout, levdbg, &
377 inforodas, iodestatus_count, ioderegime_count )
378
379 implicit none
380
381 ! subr arguments
382 integer lunout, levdbg
383 integer inforodas(2,6), iodestatus_count(6), ioderegime_count(6)
384
385 ! local variables
386 integer i, j
387 character*80 msg
388
389
390 msg = ' '
391 call peg_debugmsg( lunout, levdbg, msg )
392 msg = 'output from dump_cbmz_gasodeinfo'
393 call peg_debugmsg( lunout, levdbg, msg )
394 write(msg,9100) 'oderegime(1-6)', (ioderegime_count(i), i=1,6)
395 call peg_debugmsg( lunout, levdbg, msg )
396 write(msg,9100) 'odestatus(1-6)', (iodestatus_count(i), i=1,6)
397 call peg_debugmsg( lunout, levdbg, msg )
398
399 write(msg,9200) &
400 'inforodas(1-3)', ((inforodas(j,i), j=1,2), i=1,3)
401 call peg_debugmsg( lunout, levdbg, msg )
402 write(msg,9200) &
403 'inforodas(4-6)', ((inforodas(j,i), j=1,2), i=4,6)
404 call peg_debugmsg( lunout, levdbg, msg )
405
406 9100 format( a, 6i11 )
407 9200 format( a, 3( i11, '--', i9.9 ) )
408
409 return
410 end subroutine print_gasode_stats
411
412
413
414 !***********************************************************************
415 ! < xx.> subr gasodesolver_rodas
416 !
417 ! purpose: interfaces to rodas ode solver
418 !
419 !-----------------------------------------------------------------------
420
421 subroutine gasodesolver_rodas( tstart, tstop, iok, &
422 isvode, jsvode, ksvode, iregime, &
423 mgaschem, lunerr, lunout, levdbg_err, levdbg_info, &
424 i_force_dump, inforodas, iodestatus_count, &
425 abs_error, rel_error, trun, &
426 tempbox, pressbox, airdenbox, cboxnew, cboxold, &
427 rk_m1, rk_m2, rk_m3, rk_m4 )
428
429 use module_data_cbmz
430 use module_cbmz_rodas_prep, only: &
431 cbmz_v02r01_mapconcs, cbmz_v02r01_maprates, cbmz_v02r01_torodas, &
432 cbmz_v02r02_mapconcs, cbmz_v02r02_maprates, cbmz_v02r02_torodas, &
433 cbmz_v02r03_mapconcs, cbmz_v02r03_maprates, cbmz_v02r03_torodas, &
434 cbmz_v02r04_mapconcs, cbmz_v02r04_maprates, cbmz_v02r04_torodas, &
435 cbmz_v02r05_mapconcs, cbmz_v02r05_maprates, cbmz_v02r05_torodas, &
436 cbmz_v02r06_mapconcs, cbmz_v02r06_maprates, cbmz_v02r06_torodas
437
438 implicit none
439
440 ! subr arguments
441 integer iok, isvode, jsvode, ksvode, i_force_dump, iregime, &
442 mgaschem, lunerr, lunout, levdbg_err, levdbg_info
443 integer inforodas(2,6), iodestatus_count(6)
444 real tstart, tstop, abs_error, rel_error, trun
445 real tempbox, pressbox, airdenbox
446 real cboxnew(ngas_z), cboxold(ngas_z)
447 real rk_m1(nrxn_m1), rk_m2(nrxn_m2), rk_m3(nrxn_m3), rk_m4(nrxn_m4)
448
449 ! local variables
450 integer :: ia, idum, idydt_sngldble, ig, l, ntot
451 integer, save :: nrodas_failures = 0
452 integer, dimension(6) :: inforodas_cur
453
454 real hmin, hstart, taa, tzz
455 real atolvec(ngas_z), rtolvec(ngas_z), &
456 stot(ngas_z), &
457 yposlimit(ngas_z), yneglimit(ngas_z)
458 real sfixedkpp(nfixed_kppmax), rconstkpp(nreact_kppmax)
459
460 character*80 msg
461
462 ! map reaction rate constants (pegasus --> kpp)
463 ! map concentrations (cboxold --> stot)
464 ! dump rates (for debugging)
465 if (iregime .eq. 1) then
466 call cbmz_v02r01_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
467 rconstkpp )
468 call cbmz_v02r01_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
469
470 else if (iregime .eq. 2) then
471 call cbmz_v02r02_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
472 rconstkpp )
473 call cbmz_v02r02_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
474
475 else if (iregime .eq. 3) then
476 call cbmz_v02r03_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
477 rconstkpp )
478 call cbmz_v02r03_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
479
480 else if (iregime .eq. 4) then
481 call cbmz_v02r04_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
482 rconstkpp )
483 call cbmz_v02r04_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
484
485 else if (iregime .eq. 5) then
486 call cbmz_v02r05_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
487 rconstkpp )
488 call cbmz_v02r05_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
489
490 else
491 call cbmz_v02r06_maprates( rk_m1, rk_m2, rk_m3, rk_m4, &
492 rconstkpp )
493 call cbmz_v02r06_mapconcs( 0, ntot, stot, sfixedkpp, cboxold )
494 end if
495
496 ! set parameters for rodas call
497 do l = 1, ntot
498 atolvec(l) = abs_error
499 rtolvec(l) = rel_error
500 yposlimit(l) = 1.0e20
501 yneglimit(l) = -1.0e8
502 end do
503
504 taa = tstart
505 tzz = tstop
506 hmin = 1.0e-5
507 hstart = 60.0
508 idydt_sngldble = 1
509
510 ! call rodas integrator
511 ! subr cbmz_v02r06_torodas(
512 ! + ngas, taa, tzz,
513 ! + stot, atol, rtol, yposlimit, yneglimit,
514 ! + hmin, hstart,
515 ! + inforodas_cur, iok, lunerr, idydt_sngldble )
516
517 if (iregime .eq. 1) then
518 call cbmz_v02r01_torodas( &
519 ntot, taa, tzz, &
520 stot, atolvec, rtolvec, yposlimit, yneglimit, &
521 sfixedkpp, rconstkpp, &
522 hmin, hstart, &
523 inforodas_cur, iok, lunerr, idydt_sngldble )
524
525 else if (iregime .eq. 2) then
526 call cbmz_v02r02_torodas( &
527 ntot, taa, tzz, &
528 stot, atolvec, rtolvec, yposlimit, yneglimit, &
529 sfixedkpp, rconstkpp, &
530 hmin, hstart, &
531 inforodas_cur, iok, lunerr, idydt_sngldble )
532
533 else if (iregime .eq. 3) then
534 call cbmz_v02r03_torodas( &
535 ntot, taa, tzz, &
536 stot, atolvec, rtolvec, yposlimit, yneglimit, &
537 sfixedkpp, rconstkpp, &
538 hmin, hstart, &
539 inforodas_cur, iok, lunerr, idydt_sngldble )
540
541 else if (iregime .eq. 4) then
542 call cbmz_v02r04_torodas( &
543 ntot, taa, tzz, &
544 stot, atolvec, rtolvec, yposlimit, yneglimit, &
545 sfixedkpp, rconstkpp, &
546 hmin, hstart, &
547 inforodas_cur, iok, lunerr, idydt_sngldble )
548
549 else if (iregime .eq. 5) then
550 call cbmz_v02r05_torodas( &
551 ntot, taa, tzz, &
552 stot, atolvec, rtolvec, yposlimit, yneglimit, &
553 sfixedkpp, rconstkpp, &
554 hmin, hstart, &
555 inforodas_cur, iok, lunerr, idydt_sngldble )
556
557 else
558 call cbmz_v02r06_torodas( &
559 ntot, taa, tzz, &
560 stot, atolvec, rtolvec, yposlimit, yneglimit, &
561 sfixedkpp, rconstkpp, &
562 hmin, hstart, &
563 inforodas_cur, iok, lunerr, idydt_sngldble )
564 end if
565
566
567 ! increment odeinfo counters
568 if (iok .gt. 0) then
569 if (inforodas_cur(6) .le. 0) then
570 ia = 1
571 else
572 ia = 2
573 end if
574 else
575 ia = 3
576 end if
577 iodestatus_count(ia) = iodestatus_count(ia) + 1
578 ! do following to avoid overflow of the "inforodas" numbers
579 ! inforodas(2,i) contains rightmost 9 digits of each inforodas number
580 ! inforodas(1,i) contains any higher digits of each inforodas number
581 do ia = 1, 6
582 idum = inforodas(2,ia) + inforodas_cur(ia)
583 inforodas(1,ia) = inforodas(1,ia) + (idum/1000000000)
584 inforodas(2,ia) = mod(idum, 1000000000)
585 end do
586
587
588 ! map concentrations (stot --> cboxnew)
589 if (iregime .eq. 1) then
590 call cbmz_v02r01_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
591 else if (iregime .eq. 2) then
592 call cbmz_v02r02_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
593 else if (iregime .eq. 3) then
594 call cbmz_v02r03_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
595 else if (iregime .eq. 4) then
596 call cbmz_v02r04_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
597 else if (iregime .eq. 5) then
598 call cbmz_v02r05_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
599 else
600 call cbmz_v02r06_mapconcs( 1, ntot, stot, sfixedkpp, cboxnew )
601 end if
602
603
604 ! diagnostic output if integration fails OR if i_force_dump > 0
605 if (iok .gt. 0) then
606 if (i_force_dump .le. 0) goto 20000
607 else
608 nrodas_failures = nrodas_failures + 1
609 if (nrodas_failures .gt. 100) goto 20000
610 end if
611
612 msg = ' '
613 call peg_debugmsg( lunout, levdbg_err, msg )
614 if (iok .gt. 0) then
615 msg = '*** gasodesolver_rodas forced dump'
616 else
617 write(msg,*) '*** gasodesolver_rodas failure no.', &
618 nrodas_failures
619 end if
620 call peg_debugmsg( lunout, levdbg_err, msg )
621 msg = 'iregime, iok, i, j, k / t'
622 call peg_debugmsg( lunout, levdbg_err, msg )
623 write(msg,97010) iregime, iok, isvode, jsvode, ksvode
624 call peg_debugmsg( lunout, levdbg_err, msg )
625 write(msg,97020) trun
626 call peg_debugmsg( lunout, levdbg_err, msg )
627 msg = 'inforodas_cur(1-6) ='
628 call peg_debugmsg( lunout, levdbg_err, msg )
629 write(msg,97010) inforodas_cur
630 call peg_debugmsg( lunout, levdbg_err, msg )
631 msg = &
632 'tstart, tstop, abs_error, rel_error / temp, press, cair, cos_sza ='
633 call peg_debugmsg( lunout, levdbg_err, msg )
634 write(msg,97020) tstart, tstop, abs_error, rel_error
635 call peg_debugmsg( lunout, levdbg_err, msg )
636 write(msg,97020) tempbox, pressbox, airdenbox, -99.0
637 call peg_debugmsg( lunout, levdbg_err, msg )
638
639 idum = 0
640 do ig = nreact_kppmax, 1, -1
641 if ((idum .eq. 0) .and. (rconstkpp(ig) .ne. 0.0)) idum = ig
642 end do
643 msg = 'ngas_z, nrconst_nonzero ='
644 call peg_debugmsg( lunout, levdbg_err, msg )
645 write(msg,97010) ngas_z, idum
646 call peg_debugmsg( lunout, levdbg_err, msg )
647 msg = 'l, name, cboxold, cboxnew for l=1,ngas_z'
648 call peg_debugmsg( lunout, levdbg_err, msg )
649 do l = 1, ngas_z
650 write(msg,97030) l, name_z(l), cboxold(l), cboxnew(l)
651 call peg_debugmsg( lunout, levdbg_err, msg )
652 end do
653 msg = 'rconst for i=1,nrconst_nonzero'
654 call peg_debugmsg( lunout, levdbg_err, msg )
655 do ia = 1, idum, 4
656 write(msg,97020) ( rconstkpp(ig), ig = ia, min(ia+3,idum) )
657 call peg_debugmsg( lunout, levdbg_err, msg )
658 end do
659
660 97010 format( 6i12 )
661 97020 format( 4(1pe18.10) )
662 97030 format(( i3, 1x, a, 2(1pe18.10) ))
663
664
665 ! force non-negative values
666 20000 do l = 1, ngas_z
667 cboxnew(l) = max( cboxnew(l), 0.0 )
668 end do
669
670 return
671 end subroutine gasodesolver_rodas
672
673
674
675 !***********************************************************************
676 ! < xx.> subr gasodesolver_lsodes
677 !
678 ! purpose: interface to lsodes ode solver
679 !
680 ! author : Rahul A. Zaveri
681 ! date : May, 2000
682 !
683 !-----------------------------------------------------------------------
684
685 subroutine gasodesolver_lsodes( tstart, tstop, iok, &
686 isvode, jsvode, ksvode, iregime, &
687 mgaschem, lunerr, lunout, levdbg_err, levdbg_info, &
688 i_force_dump, iodestatus_count, &
689 abs_error, rel_error, trun, &
690 tempbox, pressbox, airdenbox, cboxnew, cboxold, &
691 rk_m1, rk_m2, rk_m3, rk_m4 )
692
693 use module_data_cbmz
694 use module_cbmz_lsodes_solver, only: lsodes_solver, xsetf, &
695 set_lsodes_common_vars
696 implicit none
697
698 ! subr arguments
699 integer i, iok, isvode, jsvode, ksvode, i_force_dump, iregime, &
700 mgaschem, lunerr, lunout, levdbg_err, levdbg_info
701 integer iodestatus_count(6)
702 real tstart, tstop, abs_error, rel_error, trun
703 real tempbox, pressbox, airdenbox
704 real cboxnew(ngas_z), cboxold(ngas_z)
705 real rk_m1(nrxn_m1), rk_m2(nrxn_m2), rk_m3(nrxn_m3), rk_m4(nrxn_m4)
706
707 ! lsodes parameters and local variables
708 integer itoler, itask, iopt, mf, lwm, nrdim, nidim
709 integer nruserpar, niuserpar
710 parameter( itoler = 1, itask = 1, iopt = 1, mf= 222 )
711 parameter( lwm = 3*ngas_tot*ngas_tot + 12*ngas_tot )
712 parameter( nrdim = 20 + 9*ngas_tot + lwm )
713 parameter( nidim = 31 + ngas_tot + ngas_tot*ngas_tot )
714 parameter( nruserpar = 5 + nrxn_m1 + nrxn_m2 + nrxn_m3 + nrxn_m4)
715 parameter( niuserpar = ngas_z + 1 )
716
717 integer ia, idum, ig, ioffset, istate, iwork(nidim), l
718 integer ntotvec(1), iuserpar(niuserpar)
719 integer indx(ngas_z)
720 integer, save :: iflagout = 0
721 integer, save :: nlsodes_failures = 0
722
723 real dtchem, rwork(nrdim), stot(ngas_tot)
724 real atolvec(1), rtolvec(1), ruserpar(nruserpar)
725
726 character*80 msg
727
728
729
730 iok = 1 ! reset
731
732 call set_lsodes_common_vars()
733
734 ! sets gas species indices for iregime
735 call setgasindices( iregime, indx )
736
737 ! map cboxold --> stot
738 call mapgasspecies( cboxold, stot, 0, iregime, indx )
739
740 !----------------------------------------------------------------------
741 ! set number of species (ntot) for the selected regime for LSODES
742 if (iregime .eq. 1) then
743 ntotvec(1) = ngas_r1
744 else if (iregime .eq. 2) then
745 ntotvec(1) = ngas_r2
746 else if (iregime .eq. 3) then
747 ntotvec(1) = ngas_r3
748 else if (iregime .eq. 4) then
749 ntotvec(1) = ngas_r4
750 else if (iregime .eq. 5) then
751 ntotvec(1) = ngas_r5
752 else
753 ntotvec(1) = ngas_r6
754 end if
755
756 100 continue
757
758 ! set other LSODES parameters...
759 iwork(6) = 1000 ! max iterations for a time step
760 iwork(7) = 1
761 istate = 1
762 rwork(6) = dtchem
763 if(iflagout.eq.0)then
764 call xsetf(iflagout)
765 endif
766
767 atolvec(1) = abs_error
768 rtolvec(1) = rel_error
769
770 do ig = 1, 5
771 ruserpar(ig) = ig*7.0
772 end do
773 ruserpar(1) = cboxold(ih2o_z)
774 ruserpar(2) = cboxold(ich4_z)
775 ruserpar(3) = cboxold(io2_z)
776 ruserpar(4) = cboxold(in2_z)
777 ruserpar(5) = cboxold(ih2_z)
778 ioffset = 5
779 do ig = 1, nrxn_m1
780 ruserpar(ioffset+ig) = rk_m1(ig)
781 end do
782 ioffset = ioffset + nrxn_m1
783 do ig = 1, nrxn_m2
784 ruserpar(ioffset+ig) = rk_m2(ig)
785 end do
786 ioffset = ioffset + nrxn_m2
787 do ig = 1, nrxn_m3
788 ruserpar(ioffset+ig) = rk_m3(ig)
789 end do
790 ioffset = ioffset + nrxn_m3
791 do ig = 1, nrxn_m4
792 ruserpar(ioffset+ig) = rk_m4(ig)
793 end do
794
795 iuserpar(1) = iregime
796 do ig = 1, ngas_z
797 iuserpar(1+ig) = indx(ig)
798 end do
799
800 call lsodes_solver( &
801 gasode_cbmz, ntotvec, stot, tstart, tstop, &
802 itoler, rtolvec, atolvec, itask, istate, iopt, &
803 rwork, nrdim, iwork, nidim, jcs, mf, &
804 ruserpar, nruserpar, iuserpar, niuserpar )
805
806 if (istate .le. 0) iok = -1
807
808
809 ! increment odeinfo counters
810 if (iok .gt. 0) then
811 ia = 4
812 else
813 ia = 5
814 end if
815 iodestatus_count(ia) = iodestatus_count(ia) + 1
816
817
818 ! map stot --> cboxnew
819 call mapgasspecies( cboxnew, stot, 1, iregime, indx )
820
821
822 ! do diagnostic output if integration fails OR if i_force_dump > 0
823 if (iok .gt. 0) then
824 if (i_force_dump .le. 0) goto 20000
825 else
826 nlsodes_failures = nlsodes_failures + 1
827 end if
828
829 msg = ' '
830 call peg_debugmsg( lunout, levdbg_err, msg )
831 if (iok .gt. 0) then
832 msg = '*** gasodesolver_lsodes forced dump'
833 else
834 write(msg,*) '*** gasodesolver_lsodes failure no.', &
835 nlsodes_failures
836 end if
837 call peg_debugmsg( lunout, levdbg_err, msg )
838 msg = 'iregime, iok, i, j, k / t'
839 call peg_debugmsg( lunout, levdbg_err, msg )
840 write(msg,97010) iregime, iok, isvode, jsvode, ksvode
841 call peg_debugmsg( lunout, levdbg_err, msg )
842 write(msg,97020) trun
843 call peg_debugmsg( lunout, levdbg_err, msg )
844 if (nlsodes_failures .gt. 1000) then
845 write(msg,*) '*** exceeded lsodes failure limit =', 1000
846 call peg_debugmsg( lunout, levdbg_err, msg )
847 call peg_error_fatal( lunerr, msg )
848 end if
849 if (nlsodes_failures .gt. 100) goto 20000
850
851 write(msg,*) 'istate -', istate
852 call peg_debugmsg( lunout, levdbg_err, msg )
853 msg = &
854 'tstart, tstop, abs_error, rel_error / temp, press, cair, cos_sza ='
855 call peg_debugmsg( lunout, levdbg_err, msg )
856 write(msg,97020) tstart, tstop, abs_error, rel_error
857 call peg_debugmsg( lunout, levdbg_err, msg )
858 write(msg,97020) tempbox, pressbox, airdenbox, -99.0
859 call peg_debugmsg( lunout, levdbg_err, msg )
860
861 idum = nrxn_m1 + nrxn_m2 + nrxn_m3 + nrxn_m4
862 msg = 'ngas_z, nrconst_m1+m2+m3+m4 ='
863 call peg_debugmsg( lunout, levdbg_err, msg )
864 write(msg,97010) ngas_z, idum
865 call peg_debugmsg( lunout, levdbg_err, msg )
866 msg = 'l, name, cboxold, cboxnew for l=1,ngas_z'
867 call peg_debugmsg( lunout, levdbg_err, msg )
868 do l = 1, ngas_z
869 write(msg,97030) l, name_z(l), cboxold(l), cboxnew(l)
870 call peg_debugmsg( lunout, levdbg_err, msg )
871 end do
872 msg = 'rconst for i=1,nrconst_nonzero'
873 call peg_debugmsg( lunout, levdbg_err, msg )
874 do ia = 1, nrxn_m1, 4
875 write(msg,97020) ( rk_m1(ig), ig = ia, min(ia+3,nrxn_m1) )
876 call peg_debugmsg( lunout, levdbg_err, msg )
877 end do
878 do ia = 1, nrxn_m2, 4
879 write(msg,97020) ( rk_m2(ig), ig = ia, min(ia+3,nrxn_m2) )
880 call peg_debugmsg( lunout, levdbg_err, msg )
881 end do
882 do ia = 1, nrxn_m3, 4
883 write(msg,97020) ( rk_m3(ig), ig = ia, min(ia+3,nrxn_m3) )
884 call peg_debugmsg( lunout, levdbg_err, msg )
885 end do
886 do ia = 1, nrxn_m4, 4
887 write(msg,97020) ( rk_m4(ig), ig = ia, min(ia+3,nrxn_m4) )
888 call peg_debugmsg( lunout, levdbg_err, msg )
889 end do
890
891 97010 format( 6i12 )
892 97020 format( 4(1pe18.10) )
893 97030 format(( i3, 1x, a, 2(1pe18.10) ))
894
895
896 ! force non-negative values
897 20000 do l = 1, ngas_z
898 cboxnew(l) = max( cboxnew(l), 0.0 )
899 end do
900
901 return
902 end subroutine gasodesolver_lsodes
903
904
905
906 !***********************************************************************
907 ! < 2.> subr selectgasregime
908 !
909 ! purpose: selects an optimum combination of gas-phase
910 ! mechanisms based on sensitivity of [OH]
911 ! concentrations to some lumped structure
912 ! hydrocarbon groups concentrations and [DMS]
913 ! concentration.
914 !
915 ! input : cbox = full species concentrations array (mol/cc)
916 !
917 ! output: iregime = 1 : com
918 ! = 2 : com + urb
919 ! = 3 : com + urb + bio
920 ! = 4 : com + mar
921 ! = 5 : com + urb + mar
922 ! = 6 : com + urb + bio + mar
923 ! ntot = number of gas-phase species in the selected mechanism
924 !
925 ! author: Rahul A. Zaveri
926 ! date : April 2000
927 !
928 !---------------------------------------------------------------------
929
930 subroutine selectgasregime( iregime, iregime_forced, cbox, &
931 igaschem_allowed_m1, igaschem_allowed_m2, &
932 igaschem_allowed_m3, igaschem_allowed_m4 )
933
934 use module_data_cbmz
935 implicit none
936
937 ! subr arguments
938 integer iregime, iregime_forced
939 integer igaschem_allowed_m1, igaschem_allowed_m2, &
940 igaschem_allowed_m3, igaschem_allowed_m4
941 real cbox(ngas_z)
942
943 ! local variables
944 integer iwork(6)
945 integer m_m1, m_m2, m_m3, m_m4
946 real cut_molecpcc
947
948
949 cut_molecpcc = 5.e+6 ! molecules/cc
950
951 ! initialize mechanism flags
952 m_m1 = 1 ! 1 (always)
953 m_m2 = 0 ! 0 or 1
954 m_m3 = 0 ! 0 or 2
955 m_m4 = 0 ! 0 or 3
956
957 if (igaschem_allowed_m2 .gt. 0) then
958 if (cbox(ipar_z ) .gt. cut_molecpcc) m_m2 = 1
959 if (cbox(iaone_z ) .gt. cut_molecpcc) m_m2 = 1
960 if (cbox(imgly_z ) .gt. cut_molecpcc) m_m2 = 1
961 if (cbox(ieth_z ) .gt. cut_molecpcc) m_m2 = 1
962 if (cbox(iolet_z ) .gt. cut_molecpcc) m_m2 = 1
963 if (cbox(iolei_z ) .gt. cut_molecpcc) m_m2 = 1
964 if (cbox(ixyl_z ) .gt. cut_molecpcc) m_m2 = 1
965 if (cbox(icres_z ) .gt. cut_molecpcc) m_m2 = 1
966 if (cbox(ito2_z ) .gt. cut_molecpcc) m_m2 = 1
967 if (cbox(icro_z ) .gt. cut_molecpcc) m_m2 = 1
968 if (cbox(iopen_z ) .gt. cut_molecpcc) m_m2 = 1
969 if (cbox(ionit_z ) .gt. cut_molecpcc) m_m2 = 1
970 if (cbox(irooh_z ) .gt. cut_molecpcc) m_m2 = 1
971 if (cbox(iro2_z ) .gt. cut_molecpcc) m_m2 = 1
972 if (cbox(iano2_z ) .gt. cut_molecpcc) m_m2 = 1
973 if (cbox(inap_z ) .gt. cut_molecpcc) m_m2 = 1
974 if (cbox(ixo2_z ) .gt. cut_molecpcc) m_m2 = 1
975 if (cbox(ixpar_z ) .gt. cut_molecpcc) m_m2 = 1
976 end if
977
978 if (igaschem_allowed_m3 .gt. 0) then
979 if (cbox(iisop_z ) .gt. cut_molecpcc) m_m3 = 2
980 if (cbox(iisoprd_z ) .gt. cut_molecpcc) m_m3 = 2
981 if (cbox(iisopp_z ) .gt. cut_molecpcc) m_m3 = 2
982 if (cbox(iisopn_z ) .gt. cut_molecpcc) m_m3 = 2
983 if (cbox(iisopo2_z ) .gt. cut_molecpcc) m_m3 = 2
984 end if
985
986 if (igaschem_allowed_m4 .gt. 0) then
987 if (cbox(idms_z ) .gt. cut_molecpcc) m_m4 = 3
988 if (cbox(imsa_z ) .gt. cut_molecpcc) m_m4 = 3
989 if (cbox(idmso_z ) .gt. cut_molecpcc) m_m4 = 3
990 if (cbox(idmso2_z ) .gt. cut_molecpcc) m_m4 = 3
991 if (cbox(ich3so2h_z ) .gt. cut_molecpcc) m_m4 = 3
992 if (cbox(ich3sch2oo_z ) .gt. cut_molecpcc) m_m4 = 3
993 if (cbox(ich3so2_z ) .gt. cut_molecpcc) m_m4 = 3
994 if (cbox(ich3so3_z ) .gt. cut_molecpcc) m_m4 = 3
995 if (cbox(ich3so2oo_z ) .gt. cut_molecpcc) m_m4 = 3
996 if (cbox(ich3so2ch2oo_z) .gt. cut_molecpcc) m_m4 = 3
997 if (cbox(imtf_z ) .gt. cut_molecpcc) m_m4 = 3
998 end if
999
1000 iregime = m_m1 + m_m2*((2-m_m3)/2) + m_m3 + m_m4
1001
1002 ! force iregime = iregime_forced
1003 if ((iregime_forced .ge. 1) .and. (iregime_forced .le. 6)) &
1004 iregime = iregime_forced
1005
1006 return
1007 end subroutine selectgasregime
1008
1009
1010
1011 !***********************************************************************
1012 ! < 3.> subr setgasindices
1013 !
1014 ! purpose: sets gas species indices
1015 !
1016 ! author : Rahul A. Zaveri
1017 ! date : May, 2000
1018 !
1019 !-----------------------------------------------------------------------
1020
1021 subroutine setgasindices( iregime, indx )
1022
1023 use module_data_cbmz
1024 implicit none
1025
1026 ! subr arguments
1027 integer iregime, indx(ngas_z)
1028
1029 ! local variables
1030 integer ilast
1031
1032 ilast = 0
1033 indx(:) = -999888777
1034
1035 goto (1,2,3,4,5,6), iregime
1036
1037 1 call setgasindex_m1( ilast, indx ) ! regime 1
1038 return
1039
1040 2 call setgasindex_m1( ilast, indx ) ! regime 2
1041 call setgasindex_m2( ilast, indx )
1042 return
1043
1044 3 call setgasindex_m1( ilast, indx ) ! regime 3
1045 call setgasindex_m2( ilast, indx )
1046 call setgasindex_m3( ilast, indx )
1047 return
1048
1049 4 call setgasindex_m1( ilast, indx ) ! regime 4
1050 call setgasindex_m4( ilast, indx )
1051 return
1052
1053 5 call setgasindex_m1( ilast, indx ) ! regime 5
1054 call setgasindex_m2( ilast, indx )
1055 call setgasindex_m4( ilast, indx )
1056 return
1057
1058 6 call setgasindex_m1( ilast, indx ) ! regime 6
1059 call setgasindex_m2( ilast, indx )
1060 call setgasindex_m3( ilast, indx )
1061 call setgasindex_m4( ilast, indx )
1062 return
1063
1064 end subroutine setgasindices
1065
1066
1067
1068
1069 !***********************************************************************
1070 ! < 4.> subr gasrateconstants
1071 !
1072 ! purpose: calls regime-dependent subrs for calculating
1073 ! gas-phase thermal reaction rate constants
1074 !
1075 ! author : Rahul A. Zaveri
1076 ! date : May, 2000
1077 !
1078 !-----------------------------------------------------------------------
1079
1080 subroutine gasrateconstants( iregime, tempbox, cair_mlc, &
1081 rk_photo, rk_param, rk_m1, rk_m2, rk_m3, rk_m4 )
1082
1083 use module_data_cbmz
1084 implicit none
1085
1086 ! subr arguments
1087 integer iregime
1088 real tempbox, cair_mlc
1089 real rk_photo(nphoto), rk_param(nperox)
1090 real rk_m1(nrxn_m1), rk_m2(nrxn_m2), rk_m3(nrxn_m3), rk_m4(nrxn_m4)
1091
1092
1093 ! iregime=1 --> do m1
1094 ! iregime=2 --> do m1, m2
1095 ! iregime=3 --> do m1, m2, m3
1096 ! iregime=4 --> do m1, --, --, m4
1097 ! iregime=5 --> do m1, m2, --, m4
1098 ! iregime=6 --> do m1, m2, m3, m4
1099
1100 call gasthermrk_m1( tempbox, cair_mlc, &
1101 rk_photo, rk_param, rk_m1, rk_m2 )
1102
1103 if ((iregime .eq. 2) .or. &
1104 (iregime .eq. 3) .or. &
1105 (iregime .eq. 5) .or. &
1106 (iregime .eq. 6)) &
1107 call gasthermrk_m2( tempbox, cair_mlc, &
1108 rk_photo, rk_param, rk_m2 )
1109
1110 if ((iregime .eq. 3) .or. &
1111 (iregime .eq. 6)) &
1112 call gasthermrk_m3( tempbox, cair_mlc, &
1113 rk_photo, rk_param, rk_m3 )
1114
1115 if ((iregime .eq. 4) .or. &
1116 (iregime .eq. 5) .or. &
1117 (iregime .eq. 6)) &
1118 call gasthermrk_m4( tempbox, cair_mlc, &
1119 rk_photo, rk_param, rk_m4 )
1120
1121 return
1122 end subroutine gasrateconstants
1123
1124
1125
1126 !***********************************************************************
1127 ! < 3.> subr setgasindex_m1
1128 !
1129 ! purpose: defines gas species indices for regime 1
1130 !
1131 ! author : Rahul A. Zaveri
1132 ! date : December, 1998
1133 !
1134 !-----------------------------------------------------------------------
1135
1136 subroutine setgasindex_m1( ilast, indx )
1137
1138 use module_data_cbmz
1139 implicit none
1140
1141 ! subr arguments
1142 integer ilast, indx(ngas_z)
1143
1144
1145 indx(ino_z) = 1
1146 indx(ino2_z) = 2
1147 indx(ino3_z) = 3
1148 indx(in2o5_z) = 4
1149 indx(ihono_z) = 5
1150 indx(ihno3_z) = 6
1151 indx(ihno4_z) = 7
1152 indx(io3_z) = 8
1153 indx(io1d_z) = 9
1154 indx(io3p_z) = 10
1155 indx(ioh_z) = 11
1156 indx(iho2_z) = 12
1157 indx(ih2o2_z) = 13
1158 indx(ico_z) = 14
1159 indx(iso2_z) = 15
1160 indx(ih2so4_z) = 16
1161 indx(inh3_z) = 17
1162 indx(ihcl_z) = 18
1163 indx(ic2h6_z) = 19
1164 indx(ich3o2_z) = 20
1165 indx(iethp_z) = 21
1166 indx(ihcho_z) = 22
1167 indx(ich3oh_z) = 23
1168 indx(ic2h5oh_z) = 24
1169 indx(ich3ooh_z) = 25
1170 indx(iethooh_z) = 26
1171 indx(iald2_z) = 27
1172 indx(ihcooh_z) = 28
1173 indx(ircooh_z) = 29
1174 indx(ic2o3_z) = 30
1175 indx(ipan_z) = 31
1176
1177 ilast = indx(ipan_z)
1178
1179 return
1180 end subroutine setgasindex_m1
1181
1182
1183
1184 !***********************************************************************
1185 ! < 4.> subr setgasindex_m2
1186 !
1187 ! purpose: defines gas species indices for regime 2
1188 !
1189 ! author : Rahul A. Zaveri
1190 ! date : December, 1998
1191 !
1192 !-----------------------------------------------------------------------
1193 subroutine setgasindex_m2( ilast, indx )
1194
1195 use module_data_cbmz
1196 implicit none
1197
1198 ! subr arguments
1199 integer ilast, indx(ngas_z)
1200
1201
1202 indx(ipar_z) = ilast + 1
1203 indx(iaone_z) = ilast + 2
1204 indx(imgly_z) = ilast + 3
1205 indx(ieth_z) = ilast + 4
1206 indx(iolet_z) = ilast + 5
1207 indx(iolei_z) = ilast + 6
1208 indx(itol_z) = ilast + 7
1209 indx(ixyl_z) = ilast + 8
1210 indx(icres_z) = ilast + 9
1211 indx(ito2_z) = ilast + 10
1212 indx(icro_z) = ilast + 11
1213 indx(iopen_z) = ilast + 12
1214 indx(ionit_z) = ilast + 13
1215 ! indx(ipan_z) = ilast + 14
1216 ! indx(ircooh_z) = ilast + 15
1217 indx(irooh_z) = ilast + 14
1218 ! indx(ic2o3_z) = ilast + 17
1219 indx(iro2_z) = ilast + 15
1220 indx(iano2_z) = ilast + 16
1221 indx(inap_z) = ilast + 17
1222 indx(ixo2_z) = ilast + 18
1223 indx(ixpar_z) = ilast + 19
1224
1225 ilast = indx(ixpar_z)
1226
1227 return
1228 end subroutine setgasindex_m2
1229
1230
1231
1232 !***********************************************************************
1233 ! < 5.> subr setgasindex_m3
1234 !
1235 ! purpose: defines gas species indices for regime 3
1236 !
1237 ! author : Rahul A. Zaveri
1238 ! date : December, 1998
1239 !
1240 !-----------------------------------------------------------------------
1241 subroutine setgasindex_m3( ilast, indx )
1242
1243 use module_data_cbmz
1244 implicit none
1245
1246 ! subr arguments
1247 integer ilast, indx(ngas_z)
1248
1249
1250 indx(iisop_z) = ilast + 1
1251 indx(iisoprd_z) = ilast + 2
1252 indx(iisopp_z) = ilast + 3
1253 indx(iisopn_z) = ilast + 4
1254 indx(iisopo2_z) = ilast + 5
1255
1256 ilast = indx(iisopo2_z)
1257
1258 return
1259 end subroutine setgasindex_m3
1260
1261
1262
1263 !***********************************************************************
1264 ! < 6.> subr setgasindex_m4
1265 !
1266 ! purpose: defines gas species indices for regime 4
1267 !
1268 ! author : Rahul A. Zaveri
1269 ! date : December, 1998
1270 !
1271 !-----------------------------------------------------------------------
1272
1273 subroutine setgasindex_m4( ilast, indx )
1274
1275 use module_data_cbmz
1276 implicit none
1277
1278 ! subr arguments
1279 integer ilast, indx(ngas_z)
1280
1281
1282 !
1283 indx(idms_z) = ilast + 1
1284 indx(imsa_z) = ilast + 2
1285 indx(idmso_z) = ilast + 3
1286 indx(idmso2_z) = ilast + 4
1287 indx(ich3so2h_z) = ilast + 5
1288 indx(ich3sch2oo_z) = ilast + 6
1289 indx(ich3so2_z) = ilast + 7
1290 indx(ich3so3_z) = ilast + 8
1291 indx(ich3so2oo_z) = ilast + 9
1292 indx(ich3so2ch2oo_z) = ilast + 10
1293 indx(imtf_z) = ilast + 11
1294
1295 ilast = indx(imtf_z)
1296
1297 return
1298 end subroutine setgasindex_m4
1299
1300
1301
1302 !***********************************************************************
1303 ! < 9.> subr mapgasspecies
1304 !
1305 ! purpose: map gas species between stot and cbox arrays
1306 !
1307 ! author : Rahul A. Zaveri
1308 ! date : December, 1998
1309 !
1310 ! ----------------------------------------------------------------------
1311
1312 subroutine mapgasspecies( cbox, stot, imap, iregime, indx )
1313
1314 use module_data_cbmz
1315 implicit none
1316
1317 ! subr arguments
1318 integer imap, iregime, indx(ngas_z)
1319 real cbox(ngas_z)
1320 real stot(ngas_tot)
1321 !
1322 !
1323 goto (1,2,3,4,5,6), iregime
1324 !
1325 !
1326 1 call mapgas_m1( cbox, stot, imap, indx )
1327 return
1328 !
1329 !
1330 2 call mapgas_m1( cbox, stot, imap, indx )
1331 call mapgas_m2( cbox, stot, imap, indx )
1332 return
1333 !
1334 !
1335 3 call mapgas_m1( cbox, stot, imap, indx )
1336 call mapgas_m2( cbox, stot, imap, indx )
1337 call mapgas_m3( cbox, stot, imap, indx )
1338 return
1339 !
1340 !
1341 4 call mapgas_m1( cbox, stot, imap, indx )
1342 call mapgas_m4( cbox, stot, imap, indx )
1343 return
1344 !
1345 !
1346 5 call mapgas_m1( cbox, stot, imap, indx )
1347 call mapgas_m2( cbox, stot, imap, indx )
1348 call mapgas_m4( cbox, stot, imap, indx )
1349 return
1350 !
1351 !
1352 6 call mapgas_m1( cbox, stot, imap, indx )
1353 call mapgas_m2( cbox, stot, imap, indx )
1354 call mapgas_m3( cbox, stot, imap, indx )
1355 call mapgas_m4( cbox, stot, imap, indx )
1356 return
1357
1358 end subroutine mapgasspecies
1359
1360
1361
1362 !***********************************************************************
1363 ! <10.> subr mapgas_m1
1364 !
1365 ! purpose: maps gas species between stot and cbox arrays for mechanism 1
1366 !
1367 ! author : Rahul A. Zaveri
1368 ! date : December, 1998
1369 !
1370 ! ----------------------------------------------------------------------
1371
1372 subroutine mapgas_m1( cbox, stot, imap, indx )
1373
1374 use module_data_cbmz
1375 implicit none
1376
1377 ! subr arguments
1378 integer imap, indx(ngas_z)
1379 real cbox(ngas_z)
1380 real stot(ngas_tot)
1381
1382 if(imap.eq.0)then ! map cbox --> stot (both molec/cc)
1383 stot(indx(ino_z)) = cbox(ino_z)
1384 stot(indx(ino2_z)) = cbox(ino2_z)
1385 stot(indx(ino3_z)) = cbox(ino3_z)
1386 stot(indx(in2o5_z)) = cbox(in2o5_z)
1387 stot(indx(ihono_z)) = cbox(ihono_z)
1388 stot(indx(ihno3_z)) = cbox(ihno3_z)
1389 stot(indx(ihno4_z)) = cbox(ihno4_z)
1390 stot(indx(io3_z)) = cbox(io3_z)
1391 stot(indx(io1d_z)) = cbox(io1d_z)
1392 stot(indx(io3p_z)) = cbox(io3p_z)
1393 stot(indx(ioh_z)) = cbox(ioh_z)
1394 stot(indx(iho2_z)) = cbox(iho2_z)
1395 stot(indx(ih2o2_z)) = cbox(ih2o2_z)
1396 stot(indx(ico_z)) = cbox(ico_z)
1397 stot(indx(iso2_z)) = cbox(iso2_z)
1398 stot(indx(ih2so4_z)) = cbox(ih2so4_z)
1399 stot(indx(inh3_z)) = cbox(inh3_z)
1400 stot(indx(ihcl_z)) = cbox(ihcl_z)
1401 stot(indx(ic2h6_z)) = cbox(ic2h6_z)
1402 stot(indx(ich3o2_z)) = cbox(ich3o2_z)
1403 stot(indx(iethp_z)) = cbox(iethp_z)
1404 stot(indx(ihcho_z)) = cbox(ihcho_z)
1405 stot(indx(ich3oh_z)) = cbox(ich3oh_z)
1406 stot(indx(ic2h5oh_z)) = cbox(ic2h5oh_z)
1407 stot(indx(ich3ooh_z)) = cbox(ich3ooh_z)
1408 stot(indx(iethooh_z)) = cbox(iethooh_z)
1409 stot(indx(iald2_z)) = cbox(iald2_z)
1410 stot(indx(ihcooh_z)) = cbox(ihcooh_z)
1411 stot(indx(ircooh_z)) = cbox(ircooh_z)
1412 stot(indx(ic2o3_z)) = cbox(ic2o3_z)
1413 stot(indx(ipan_z)) = cbox(ipan_z)
1414 !
1415 else ! map stot --> cbox (both molec/cc)
1416 cbox(ino_z) = stot(indx(ino_z))
1417 cbox(ino2_z) = stot(indx(ino2_z))
1418 cbox(ino3_z) = stot(indx(ino3_z))
1419 cbox(in2o5_z) = stot(indx(in2o5_z))
1420 cbox(ihono_z) = stot(indx(ihono_z))
1421 cbox(ihno3_z) = stot(indx(ihno3_z))
1422 cbox(ihno4_z) = stot(indx(ihno4_z))
1423 cbox(io3_z) = stot(indx(io3_z))
1424 cbox(io1d_z) = stot(indx(io1d_z))
1425 cbox(io3p_z) = stot(indx(io3p_z))
1426 cbox(ioh_z) = stot(indx(ioh_z))
1427 cbox(iho2_z) = stot(indx(iho2_z))
1428 cbox(ih2o2_z) = stot(indx(ih2o2_z))
1429 cbox(ico_z) = stot(indx(ico_z))
1430 cbox(iso2_z) = stot(indx(iso2_z))
1431 cbox(ih2so4_z) = stot(indx(ih2so4_z))
1432 cbox(inh3_z) = stot(indx(inh3_z))
1433 cbox(ihcl_z) = stot(indx(ihcl_z))
1434 cbox(ic2h6_z) = stot(indx(ic2h6_z))
1435 cbox(ich3o2_z) = stot(indx(ich3o2_z))
1436 cbox(iethp_z) = stot(indx(iethp_z))
1437 cbox(ihcho_z) = stot(indx(ihcho_z))
1438 cbox(ich3oh_z) = stot(indx(ich3oh_z))
1439 cbox(ic2h5oh_z) = stot(indx(ic2h5oh_z))
1440 cbox(ich3ooh_z) = stot(indx(ich3ooh_z))
1441 cbox(iethooh_z) = stot(indx(iethooh_z))
1442 cbox(iald2_z) = stot(indx(iald2_z))
1443 cbox(ihcooh_z) = stot(indx(ihcooh_z))
1444 cbox(ircooh_z) = stot(indx(ircooh_z))
1445 cbox(ic2o3_z) = stot(indx(ic2o3_z))
1446 cbox(ipan_z) = stot(indx(ipan_z))
1447 endif
1448
1449 return
1450 end subroutine mapgas_m1
1451
1452
1453
1454 !***********************************************************************
1455 ! <11.> subr mapgas_m2
1456 !
1457 ! purpose: maps gas species between stot and cbox arrays for mechanism 2
1458 !
1459 ! author : Rahul A. Zaveri
1460 ! date : December, 1998
1461 !
1462 ! ----------------------------------------------------------------------
1463
1464 subroutine mapgas_m2( cbox, stot, imap, indx )
1465
1466 use module_data_cbmz
1467 implicit none
1468
1469 ! subr arguments
1470 integer imap, indx(ngas_z)
1471 real cbox(ngas_z)
1472 real stot(ngas_tot)
1473
1474 if(imap.eq.0)then ! map cbox --> stot (both molec/cc)
1475 stot(indx(ipar_z)) = cbox(ipar_z)
1476 stot(indx(iaone_z)) = cbox(iaone_z)
1477 stot(indx(imgly_z)) = cbox(imgly_z)
1478 stot(indx(ieth_z)) = cbox(ieth_z)
1479 stot(indx(iolet_z)) = cbox(iolet_z)
1480 stot(indx(iolei_z)) = cbox(iolei_z)
1481 stot(indx(itol_z)) = cbox(itol_z)
1482 stot(indx(ixyl_z)) = cbox(ixyl_z)
1483 stot(indx(icres_z)) = cbox(icres_z)
1484 stot(indx(ito2_z)) = cbox(ito2_z)
1485 stot(indx(icro_z)) = cbox(icro_z)
1486 stot(indx(iopen_z)) = cbox(iopen_z)
1487 stot(indx(ionit_z)) = cbox(ionit_z)
1488 ! stot(indx(ipan_z)) = cbox(ipan_z)
1489 ! stot(indx(ircooh_z)) = cbox(ircooh_z)
1490 stot(indx(irooh_z)) = cbox(irooh_z)
1491 ! stot(indx(ic2o3_z)) = cbox(ic2o3_z)
1492 stot(indx(iro2_z)) = cbox(iro2_z)
1493 stot(indx(iano2_z)) = cbox(iano2_z)
1494 stot(indx(inap_z)) = cbox(inap_z)
1495 stot(indx(ixo2_z)) = cbox(ixo2_z)
1496 stot(indx(ixpar_z)) = cbox(ixpar_z)
1497 !
1498 else ! map stot --> cbox (both molec/cc)
1499 cbox(ipar_z) = stot(indx(ipar_z))
1500 cbox(iaone_z) = stot(indx(iaone_z))
1501 cbox(imgly_z) = stot(indx(imgly_z))
1502 cbox(ieth_z) = stot(indx(ieth_z))
1503 cbox(iolet_z) = stot(indx(iolet_z))
1504 cbox(iolei_z) = stot(indx(iolei_z))
1505 cbox(itol_z) = stot(indx(itol_z))
1506 cbox(ixyl_z) = stot(indx(ixyl_z))
1507 cbox(icres_z) = stot(indx(icres_z))
1508 cbox(ito2_z) = stot(indx(ito2_z))
1509 cbox(icro_z) = stot(indx(icro_z))
1510 cbox(iopen_z) = stot(indx(iopen_z))
1511 cbox(ionit_z) = stot(indx(ionit_z))
1512 ! cbox(ipan_z) = stot(indx(ipan_z))
1513 ! cbox(ircooh_z) = stot(indx(ircooh_z))
1514 cbox(irooh_z) = stot(indx(irooh_z))
1515 ! cbox(ic2o3_z) = stot(indx(ic2o3_z))
1516 cbox(iro2_z) = stot(indx(iro2_z))
1517 cbox(iano2_z) = stot(indx(iano2_z))
1518 cbox(inap_z) = stot(indx(inap_z))
1519 cbox(ixo2_z) = stot(indx(ixo2_z))
1520 cbox(ixpar_z) = stot(indx(ixpar_z))
1521 endif
1522
1523 return
1524 end subroutine mapgas_m2
1525
1526
1527
1528 !***********************************************************************
1529 ! <12.> subr mapgas_m3
1530 !
1531 ! purpose: maps gas species between stot and cbox arrays for mechanism 3
1532 !
1533 ! author : Rahul A. Zaveri
1534 ! date : December, 1998
1535 !
1536 ! ----------------------------------------------------------------------
1537
1538 subroutine mapgas_m3( cbox, stot, imap, indx )
1539
1540 use module_data_cbmz
1541 implicit none
1542
1543 ! subr arguments
1544 integer imap, indx(ngas_z)
1545 real cbox(ngas_z)
1546 real stot(ngas_tot)
1547
1548 if(imap.eq.0)then ! map cbox --> stot (both molec/cc)
1549 stot(indx(iisop_z)) = cbox(iisop_z)
1550 stot(indx(iisoprd_z)) = cbox(iisoprd_z)
1551 stot(indx(iisopp_z)) = cbox(iisopp_z)
1552 stot(indx(iisopn_z)) = cbox(iisopn_z)
1553 stot(indx(iisopo2_z)) = cbox(iisopo2_z)
1554 !
1555 else ! map stot --> cbox (both molec/cc)
1556 cbox(iisop_z) = stot(indx(iisop_z))
1557 cbox(iisoprd_z) = stot(indx(iisoprd_z))
1558 cbox(iisopp_z) = stot(indx(iisopp_z))
1559 cbox(iisopn_z) = stot(indx(iisopn_z))
1560 cbox(iisopo2_z) = stot(indx(iisopo2_z))
1561 endif
1562
1563 return
1564 end subroutine mapgas_m3
1565
1566
1567
1568 !***********************************************************************
1569 ! <13.> subr mapgas_m4
1570 !
1571 ! purpose: maps gas species between stot and cbox arrays for mechanism 4
1572 !
1573 ! author : Rahul A. Zaveri
1574 ! date : December, 1998
1575 !
1576 ! ----------------------------------------------------------------------
1577
1578 subroutine mapgas_m4( cbox, stot, imap, indx )
1579
1580 use module_data_cbmz
1581 implicit none
1582
1583 ! subr arguments
1584 integer imap, indx(ngas_z)
1585 real cbox(ngas_z)
1586 real stot(ngas_tot)
1587
1588 if(imap.eq.0)then ! map cbox --> stot (both molec/cc)
1589 stot(indx(idms_z)) = cbox(idms_z)
1590 stot(indx(imsa_z)) = cbox(imsa_z)
1591 stot(indx(idmso_z)) = cbox(idmso_z)
1592 stot(indx(idmso2_z)) = cbox(idmso2_z)
1593 stot(indx(ich3so2h_z)) = cbox(ich3so2h_z)
1594 stot(indx(ich3sch2oo_z)) = cbox(ich3sch2oo_z)
1595 stot(indx(ich3so2_z)) = cbox(ich3so2_z)
1596 stot(indx(ich3so3_z)) = cbox(ich3so3_z)
1597 stot(indx(ich3so2oo_z)) = cbox(ich3so2oo_z)
1598 stot(indx(ich3so2ch2oo_z))= cbox(ich3so2ch2oo_z)
1599 stot(indx(imtf_z)) = cbox(imtf_z)
1600 !
1601 else ! map stot --> cbox (both molec/cc)
1602 cbox(idms_z) = stot(indx(idms_z))
1603 cbox(imsa_z) = stot(indx(imsa_z))
1604 cbox(idmso_z) = stot(indx(idmso_z))
1605 cbox(idmso2_z) = stot(indx(idmso2_z))
1606 cbox(ich3so2h_z) = stot(indx(ich3so2h_z))
1607 cbox(ich3sch2oo_z) = stot(indx(ich3sch2oo_z))
1608 cbox(ich3so2_z) = stot(indx(ich3so2_z))
1609 cbox(ich3so3_z) = stot(indx(ich3so3_z))
1610 cbox(ich3so2oo_z) = stot(indx(ich3so2oo_z))
1611 cbox(ich3so2ch2oo_z)= stot(indx(ich3so2ch2oo_z))
1612 cbox(imtf_z) = stot(indx(imtf_z))
1613 endif
1614
1615 return
1616 end subroutine mapgas_m4
1617
1618
1619
1620 !***********************************************************************
1621 ! <xx.> subr check_userpar
1622 !
1623 ! purpose: called by lsodes (external)
1624 ! computes time derivatives of species concentrations ds/dt.
1625 ! by calling subr. gasrate and ode
1626 !
1627 ! author : Rahul A. Zaveri
1628 ! date : December 1998
1629 !
1630 !---------------------------------------------------------------------------
1631
1632 subroutine check_userpar( ruserpar, nruserpar, iuserpar, niuserpar )
1633
1634 use module_data_cbmz
1635 implicit none
1636
1637 ! subr arguments
1638 integer nruserpar, niuserpar
1639 integer iuserpar(niuserpar)
1640 real ruserpar(nruserpar)
1641
1642 ! local variables
1643 integer i
1644 real dum
1645 character*80 msg
1646
1647 if (nruserpar .ne. (5 + nrxn_m1 + nrxn_m2 + nrxn_m3 + nrxn_m4)) then
1648 write(msg,9010) 'nruserpar', -1, nruserpar
1649 call wrf_error_fatal( msg )
1650 end if
1651
1652 if (niuserpar .ne. (ngas_z + 1)) then
1653 write(msg,9010) 'niuserpar', -1, niuserpar
1654 call wrf_error_fatal( msg )
1655 end if
1656
1657 9010 format( '*** check_userpar error -- ', a, 1x, i8, 1x, i8 )
1658
1659
1660 return
1661 end subroutine check_userpar
1662
1663
1664
1665 !***********************************************************************
1666 ! <14.> subr gasode_cbmz
1667 !
1668 ! purpose: called by lsodes (external)
1669 ! computes time derivatives of species concentrations ds/dt.
1670 ! by calling subr. gasrate and ode
1671 !
1672 ! author : Rahul A. Zaveri
1673 ! date : December 1998
1674 !
1675 !---------------------------------------------------------------------------
1676
1677 subroutine gasode_cbmz( ntot, tt, s, sdot, &
1678 ruserpar, nruserpar, iuserpar, niuserpar )
1679
1680 use module_data_cbmz
1681 implicit none
1682
1683 ! subr arguments
1684 integer ntot, nruserpar, niuserpar
1685 integer iuserpar(niuserpar)
1686 real tt
1687 real s(ngas_tot), sdot(ngas_tot), ruserpar(nruserpar)
1688
1689 ! local variables
1690 integer ig, ioffset, iregime, irxn
1691 integer indx(ngas_z)
1692 real h2o, ch4, oxygen, nitrogen, hydrogen
1693 real rk_m1(nrxn_m1), r_m1(nrxn_m1)
1694 real rk_m2(nrxn_m2), r_m2(nrxn_m2)
1695 real rk_m3(nrxn_m3), r_m3(nrxn_m3)
1696 real rk_m4(nrxn_m4), r_m4(nrxn_m4)
1697 real p_m1(ngas_tot), d_m1(ngas_tot)
1698 real p_m2(ngas_tot), d_m2(ngas_tot)
1699 real p_m3(ngas_tot), d_m3(ngas_tot)
1700 real p_m4(ngas_tot), d_m4(ngas_tot)
1701
1702
1703 ! test on userpar
1704 call check_userpar( ruserpar, nruserpar, iuserpar, niuserpar )
1705
1706 iregime = iuserpar(1)
1707 do ig = 1, ngas_z
1708 indx(ig) = iuserpar(ig+1)
1709 end do
1710
1711 h2o = ruserpar(1)
1712 ch4 = ruserpar(2)
1713 oxygen = ruserpar(3)
1714 nitrogen = ruserpar(4)
1715 hydrogen = ruserpar(5)
1716 ioffset = 5
1717 do ig = 1, nrxn_m1
1718 rk_m1(ig) = ruserpar(ioffset+ig)
1719 end do
1720 ioffset = ioffset + nrxn_m1
1721 do ig = 1, nrxn_m2
1722 rk_m2(ig) = ruserpar(ioffset+ig)
1723 end do
1724 ioffset = ioffset + nrxn_m2
1725 do ig = 1, nrxn_m3
1726 rk_m3(ig) = ruserpar(ioffset+ig)
1727 end do
1728 ioffset = ioffset + nrxn_m3
1729 do ig = 1, nrxn_m4
1730 rk_m4(ig) = ruserpar(ioffset+ig)
1731 end do
1732
1733
1734 ! initialize to zero
1735 do irxn=1,nrxn_m1
1736 r_m1(irxn) = 0.
1737 enddo
1738
1739 do irxn=1,nrxn_m2
1740 r_m2(irxn) = 0.
1741 enddo
1742
1743 do irxn=1,nrxn_m3
1744 r_m3(irxn) = 0.
1745 enddo
1746
1747 do irxn=1,nrxn_m4
1748 r_m4(irxn) = 0.
1749 enddo
1750 !
1751 !
1752 ! initialize to zero
1753 do ig=1,ngas_tot
1754 p_m1(ig) = 0.
1755 p_m2(ig) = 0.
1756 p_m3(ig) = 0.
1757 p_m4(ig) = 0.
1758 !
1759 d_m1(ig) = 0.
1760 d_m2(ig) = 0.
1761 d_m3(ig) = 0.
1762 d_m4(ig) = 0.
1763 enddo
1764
1765
1766 ! iregime=1 --> do m1
1767 ! iregime=2 --> do m1, m2
1768 ! iregime=3 --> do m1, m2, m3
1769 ! iregime=4 --> do m1, --, --, m4
1770 ! iregime=5 --> do m1, m2, --, m4
1771 ! iregime=6 --> do m1, m2, m3, m4
1772
1773 call gasrate_m1( indx, s, r_m1, r_m2, rk_m1, rk_m2, &
1774 h2o, ch4, oxygen, nitrogen, hydrogen )
1775
1776 if ((iregime .eq. 2) .or. &
1777 (iregime .eq. 3) .or. &
1778 (iregime .eq. 5) .or. &
1779 (iregime .eq. 6)) &
1780 call gasrate_m2( indx, s, r_m2, rk_m2 )
1781
1782 if ((iregime .eq. 3) .or. &
1783 (iregime .eq. 6)) &
1784 call gasrate_m3( indx, s, r_m3, rk_m3 )
1785
1786 if ((iregime .eq. 4) .or. &
1787 (iregime .eq. 5) .or. &
1788 (iregime .eq. 6)) &
1789 call gasrate_m4( indx, s, r_m4, rk_m4, oxygen )
1790
1791 call gasode_m1( indx, r_m1, p_m1, d_m1, r_m2 )
1792
1793 if ((iregime .eq. 2) .or. &
1794 (iregime .eq. 3) .or. &
1795 (iregime .eq. 5) .or. &
1796 (iregime .eq. 6)) &
1797 call gasode_m2( indx, r_m2, p_m2, d_m2 )
1798
1799 if ((iregime .eq. 3) .or. &
1800 (iregime .eq. 6)) &
1801 call gasode_m3( indx, r_m3, p_m3, d_m3 )
1802
1803 if ((iregime .eq. 4) .or. &
1804 (iregime .eq. 5) .or. &
1805 (iregime .eq. 6)) &
1806 call gasode_m4( indx, r_m4, p_m4, d_m4 )
1807
1808
1809 if (iregime .eq. 1) then
1810 ! regime = 1
1811 do ig = 1, ngas_r1
1812 sdot(ig) = real( dble(p_m1(ig)) - &
1813 dble(d_m1(ig)) )
1814 end do
1815
1816 else if (iregime .eq. 2) then
1817 ! regime = 2
1818 do ig = 1, ngas_r2
1819 sdot(ig) = real( dble(p_m1(ig)+p_m2(ig)) - &
1820 dble(d_m1(ig)+d_m2(ig)) )
1821 end do
1822
1823 else if (iregime .eq. 3) then
1824 ! regime = 3
1825 do ig = 1, ngas_r3
1826 sdot(ig) = real( dble(p_m1(ig)+p_m2(ig)+p_m3(ig)) - &
1827 dble(d_m1(ig)+d_m2(ig)+d_m3(ig)) )
1828 end do
1829
1830 else if (iregime .eq. 4) then
1831 ! regime = 4
1832 do ig = 1, ngas_r4
1833 sdot(ig) = real( dble(p_m1(ig)+p_m4(ig)) - &
1834 dble(d_m1(ig)+d_m4(ig)) )
1835 end do
1836
1837 else if (iregime .eq. 5) then
1838 ! regime = 5
1839 do ig = 1, ngas_r5
1840 sdot(ig) = real( dble(p_m1(ig)+p_m2(ig)+p_m4(ig)) - &
1841 dble(d_m1(ig)+d_m2(ig)+d_m4(ig)) )
1842 end do
1843
1844 else if (iregime .eq. 6) then
1845 ! regime = 6
1846 do ig = 1, ngas_r6
1847 sdot(ig) = real( dble(p_m1(ig)+p_m2(ig)+p_m3(ig)+p_m4(ig)) - &
1848 dble(d_m1(ig)+d_m2(ig)+d_m3(ig)+d_m4(ig)) )
1849 end do
1850
1851 end if
1852
1853 return
1854 end subroutine gasode_cbmz
1855
1856
1857
1858 !***********************************************************************
1859 ! <15.> subr jcs
1860 !
1861 ! purpose: external dummy jacobian evaluation for LSODES (when mf=222)
1862 !
1863 !-----------------------------------------------------------------------
1864
1865 subroutine jcs( ngas, tt, s, j, ian, jan, pdj, &
1866 ruserpar, nruserpar, iuserpar, niuserpar )
1867
1868 implicit none
1869 integer ngas, j, ian(*), jan(*), nruserpar, niuserpar
1870 integer iuserpar(niuserpar)
1871 real tt, s(*), pdj(*)
1872 real ruserpar(nruserpar)
1873
1874 ! test on userpar
1875 call check_userpar( ruserpar, nruserpar, iuserpar, niuserpar )
1876
1877 return
1878 end subroutine jcs
1879
1880
1881
1882 !***********************************************************************
1883 ! <16.> subr gasode_m1
1884 !
1885 ! purpose: updates production and destruction rates for mechanism 1
1886 ! background troposphere
1887 !
1888 ! author : Rahul A. Zaveri
1889 ! date : December, 1998
1890 !
1891 ! ----------------------------------------------------------------------
1892
1893 subroutine gasode_m1( indx, r_m1, p_m1, d_m1, r_m2 )
1894
1895 use module_data_cbmz
1896 implicit none
1897
1898 ! subr arguments
1899 integer indx(ngas_z)
1900 real r_m1(nrxn_m1), p_m1(ngas_tot), d_m1(ngas_tot)
1901 real r_m2(nrxn_m2)
1902
1903
1904 p_m1(indx(ino_z))= r_m1(1)+0.11*r_m1(2) &
1905 +r_m1(3)+r_m1(15)+r_m1(38)
1906 d_m1(indx(ino_z))= r_m1(17)+r_m1(18)+r_m1(23) &
1907 +r_m1(33)+r_m1(37)+r_m1(57) &
1908 +r_m1(58) &
1909 +r_m2(34)
1910 !
1911 p_m1(indx(ino2_z))= 0.89*r_m1(2)+r_m1(4) &
1912 +r_m1(5)+r_m1(6)+r_m1(17) &
1913 +r_m1(18)+r_m1(25) &
1914 +r_m1(26)+r_m1(28) &
1915 +r_m1(33)+r_m1(36) &
1916 +r_m1(37)+r_m1(37) &
1917 +r_m1(38)+r_m1(40) &
1918 +r_m1(40)+.7*r_m1(41) &
1919 +r_m1(43)+r_m1(57) &
1920 +r_m1(58)+r_m1(59) &
1921 +r_m1(60) &
1922 +r_m2(32)+r_m2(34)+r_m2(39)
1923 d_m1(indx(ino2_z))= r_m1(1)+r_m1(15)+r_m1(16) &
1924 +r_m1(19)+r_m1(24) &
1925 +r_m1(34)+r_m1(35) &
1926 +r_m1(38)+r_m1(39) &
1927 +r_m2(31)
1928 !
1929 p_m1(indx(ino3_z))= r_m1(6)+r_m1(16)+r_m1(19) &
1930 +r_m1(27)+r_m1(43)
1931 d_m1(indx(ino3_z))= r_m1(2)+r_m1(25)+r_m1(37) &
1932 +r_m1(38)+r_m1(39) &
1933 +r_m1(40)+r_m1(40) &
1934 +r_m1(41)+r_m1(52) &
1935 +r_m1(59)+r_m1(60) &
1936 +r_m2(4)+r_m2(39)
1937 !
1938 p_m1(indx(in2o5_z))= r_m1(39)
1939 d_m1(indx(in2o5_z))= r_m1(6)+r_m1(42) &
1940 +r_m1(43)
1941 !
1942 p_m1(indx(ihono_z))= r_m1(23)+r_m1(35)
1943 d_m1(indx(ihono_z))= r_m1(3)+r_m1(26)
1944 !
1945 p_m1(indx(ihno3_z))= r_m1(24)+.3*r_m1(41) &
1946 +r_m1(42)+r_m1(42) &
1947 +r_m1(52) &
1948 +r_m2(4)
1949 d_m1(indx(ihno3_z))= r_m1(4)+r_m1(27)
1950 !
1951 p_m1(indx(ihno4_z))= r_m1(34)
1952 d_m1(indx(ihno4_z))= r_m1(5)+r_m1(28) &
1953 +r_m1(36)
1954 !
1955 p_m1(indx(io3_z))= r_m1(13) &
1956 +.4*r_m2(44)
1957 d_m1(indx(io3_z))= r_m1(7)+r_m1(8)+r_m1(14) &
1958 +r_m1(18)+r_m1(19)+r_m1(20) &
1959 +r_m1(21)
1960 !
1961 p_m1(indx(io1d_z))= r_m1(8)
1962 d_m1(indx(io1d_z))= r_m1(10)+r_m1(11) &
1963 +r_m1(12)
1964 !
1965 p_m1(indx(io3p_z))= r_m1(1)+0.89*r_m1(2) &
1966 +r_m1(7)+r_m1(10)+r_m1(11)
1967 d_m1(indx(io3p_z))= r_m1(13)+r_m1(14) &
1968 +r_m1(15)+r_m1(16) &
1969 +r_m1(17)
1970 !
1971 p_m1(indx(ioh_z))= r_m1(3)+r_m1(4)+2*r_m1(9) &
1972 +2*r_m1(12)+r_m1(21) &
1973 +r_m1(33)+.7*r_m1(41) &
1974 +r_m1(53)+r_m1(54)+.3*r_m1(55) &
1975 +.5*r_m1(56)
1976 d_m1(indx(ioh_z))= r_m1(20)+r_m1(22)+r_m1(23) &
1977 +r_m1(24)+r_m1(25)+r_m1(26) &
1978 +r_m1(27)+r_m1(28)+r_m1(29) &
1979 +r_m1(30)+r_m1(44)+r_m1(45) &
1980 +r_m1(46)+r_m1(47)+r_m1(48) &
1981 +r_m1(51)+r_m1(55)+r_m1(56) &
1982 +r_m1(65) &
1983 +r_m2(3)
1984 !
1985 p_m1(indx(iho2_z))= r_m1(5)+r_m1(20)+r_m1(22) &
1986 +r_m1(25)+r_m1(30) &
1987 +r_m1(36)+r_m1(44) &
1988 +r_m1(45)+r_m1(48) &
1989 +2*r_m1(49)+r_m1(51) &
1990 +r_m1(52)+r_m1(53) &
1991 +r_m1(54)+r_m1(57) &
1992 +r_m1(58)+r_m1(59) &
1993 +r_m1(60)+.32*r_m1(63) &
1994 +.6*r_m1(64)+r_m1(65) &
1995 +r_m2(2)
1996 d_m1(indx(iho2_z))= r_m1(21)+r_m1(29) &
1997 +r_m1(31)+r_m1(31) &
1998 +r_m1(32)+r_m1(32) &
1999 +r_m1(33)+r_m1(34) &
2000 +r_m1(35)+r_m1(41) &
2001 +r_m1(61)+r_m1(62) &
2002 +r_m2(44)
2003 !
2004 p_m1(indx(ih2o2_z))= r_m1(31)+r_m1(32)
2005 d_m1(indx(ih2o2_z))= r_m1(9)+r_m1(30)
2006 !
2007 p_m1(indx(ico_z))= r_m1(49)+r_m1(50)+r_m1(51) &
2008 +r_m1(52) &
2009 +r_m2(2)
2010 d_m1(indx(ico_z))= r_m1(44)
2011 !
2012 p_m1(indx(iso2_z))= 0.0
2013 d_m1(indx(iso2_z))= r_m1(45)
2014 !
2015 p_m1(indx(ih2so4_z))= r_m1(45)
2016 d_m1(indx(ih2so4_z))= 0.0
2017 !
2018 p_m1(indx(inh3_z))= 0.0
2019 d_m1(indx(inh3_z))= 0.0
2020 !
2021 p_m1(indx(ihcl_z))= 0.0
2022 d_m1(indx(ihcl_z))= 0.0
2023 !
2024 p_m1(indx(ic2h6_z))= .2*r_m1(64)
2025 d_m1(indx(ic2h6_z))= r_m1(47)
2026 !
2027 p_m1(indx(ich3o2_z))= r_m1(46)+.7*r_m1(55) &
2028 +r_m2(2)+r_m2(34)+r_m2(39)+r_m2(49)
2029 d_m1(indx(ich3o2_z))= r_m1(57)+r_m1(59) &
2030 +r_m1(61)+r_m1(63)
2031 !
2032 p_m1(indx(iethp_z))= r_m1(47)+.5*r_m1(56)
2033 d_m1(indx(iethp_z))= r_m1(58)+r_m1(60) &
2034 +r_m1(62)+r_m1(64)
2035 !
2036 p_m1(indx(ihcho_z))= r_m1(48)+r_m1(53) &
2037 +.3*r_m1(55)+r_m1(57) &
2038 +r_m1(59)+.66*r_m1(63)
2039 d_m1(indx(ihcho_z))= r_m1(49)+r_m1(50) &
2040 +r_m1(51)+r_m1(52)
2041 !
2042 p_m1(indx(ich3oh_z))= .34*r_m1(63)
2043 d_m1(indx(ich3oh_z))= r_m1(48)
2044 !
2045 p_m1(indx(ic2h5oh_z))= 0.0
2046 d_m1(indx(ic2h5oh_z))= r_m1(65)
2047 !
2048 p_m1(indx(ich3ooh_z))= r_m1(61)
2049 d_m1(indx(ich3ooh_z))= r_m1(53)+r_m1(55)
2050 !
2051 p_m1(indx(iethooh_z))= r_m1(62)
2052 d_m1(indx(iethooh_z))= r_m1(54)+r_m1(56)
2053 !
2054 p_m1(indx(iald2_z))= r_m1(54)+.5*r_m1(56) &
2055 +r_m1(58)+r_m1(60) &
2056 +.8*r_m1(64)+r_m1(65)
2057 d_m1(indx(iald2_z))= r_m2(2)+r_m2(3)+r_m2(4)
2058 !
2059 p_m1(indx(ihcooh_z))= 0.0
2060 d_m1(indx(ihcooh_z))= 0.0
2061 !
2062 p_m1(indx(ircooh_z))= .4*r_m2(44)
2063 d_m1(indx(ircooh_z))= 0.0
2064 !
2065 p_m1(indx(ic2o3_z))= r_m2(3)+r_m2(4)+r_m2(32)
2066 d_m1(indx(ic2o3_z))= r_m2(31)+r_m2(34)+r_m2(39)+r_m2(44)+r_m2(49)
2067 !
2068 p_m1(indx(ipan_z))= r_m2(31)
2069 d_m1(indx(ipan_z))= r_m2(32)
2070
2071 return
2072 end subroutine gasode_m1
2073
2074
2075
2076 !***********************************************************************
2077 ! <17.> subr gasode_m2
2078 !
2079 ! purpose: updates production and destruction rates for mechanism 2
2080 ! anthropogenic hydrocarbons
2081 !
2082 ! author : Rahul A. Zaveri
2083 ! date : December, 1998
2084 !
2085 ! ----------------------------------------------------------------------
2086
2087 subroutine gasode_m2( indx, r_m2, p_m2, d_m2 )
2088
2089 use module_data_cbmz
2090 implicit none
2091
2092 ! subr arguments
2093 integer indx(ngas_z)
2094 real r_m2(nrxn_m2), p_m2(ngas_tot), d_m2(ngas_tot)
2095
2096
2097 p_m2(indx(ino_z))= 0.0
2098 d_m2(indx(ino_z))= r_m2(20)+r_m2(33) +r_m2(35) &
2099 +r_m2(36)+r_m2(37)
2100 !
2101 p_m2(indx(ino2_z))= .95*r_m2(20)+r_m2(30) +.84*r_m2(33) &
2102 +r_m2(35)+1.5*r_m2(36)+r_m2(37) &
2103 +r_m2(38) +r_m2(40)+1.5*r_m2(41)+r_m2(42) &
2104 +.5*r_m2(51)
2105 d_m2(indx(ino2_z))= r_m2(23)
2106 !
2107 p_m2(indx(ino3_z))= 0.0
2108 d_m2(indx(ino3_z))= +r_m2(9)+r_m2(16)+r_m2(17)+r_m2(22) &
2109 +r_m2(38) +r_m2(40)+r_m2(41)+r_m2(42)
2110 !
2111 p_m2(indx(ihno3_z))= +r_m2(9)+r_m2(22)
2112 d_m2(indx(ihno3_z))= 0.0
2113 !
2114 p_m2(indx(io3_z))= 0.0
2115 d_m2(indx(io3_z))= r_m2(10)+r_m2(12)+r_m2(13)+r_m2(26)
2116 !
2117 p_m2(indx(ioh_z))= .12*r_m2(10)+.33*r_m2(12)+.60*r_m2(13) &
2118 +.08*r_m2(26)+r_m2(27)+.23*r_m2(28)
2119 d_m2(indx(ioh_z))= r_m2(1) +r_m2(6)+r_m2(8)+r_m2(11) &
2120 +r_m2(14)+r_m2(15)+r_m2(18)+r_m2(19)+r_m2(21) &
2121 +r_m2(24)+r_m2(28)+r_m2(29)
2122 !
2123 p_m2(indx(iho2_z))= +r_m2(7)+.22*r_m2(10)+r_m2(11) &
2124 +.26*r_m2(12)+.22*r_m2(13)+r_m2(14)+r_m2(15) &
2125 +.2*r_m2(18)+.55*r_m2(19)+.95*r_m2(20) &
2126 +.6*r_m2(21)+2*r_m2(24)+r_m2(25)+.76*r_m2(26) &
2127 +.9*r_m2(27)+.9*r_m2(30)+.76*r_m2(33)+.5*r_m2(36) &
2128 +.9*r_m2(38)+.5*r_m2(41)+.54*r_m2(48)
2129 d_m2(indx(iho2_z))= r_m2(43) +r_m2(45)+r_m2(46)+r_m2(47)
2130 !
2131 p_m2(indx(ico_z))= +r_m2(7)+r_m2(9)+.24*r_m2(10) &
2132 +.31*r_m2(12)+.30*r_m2(13)+2*r_m2(24)+r_m2(25) &
2133 +.69*r_m2(26)
2134 d_m2(indx(ico_z))= 0.0
2135 !
2136 p_m2(indx(ipar_z))= 1.1*r_m2(19)
2137 d_m2(indx(ipar_z))= r_m2(1) + r_m2(53)
2138 !
2139 p_m2(indx(ich3oh_z))= .03*r_m2(12)+.04*r_m2(13)
2140 d_m2(indx(ich3oh_z))= 0.0
2141 !
2142 p_m2(indx(ihcho_z))= r_m2(10)+1.56*r_m2(11)+.57*r_m2(12)+r_m2(14) &
2143 +r_m2(24)+.7*r_m2(26)+r_m2(35)+.5*r_m2(36) &
2144 +r_m2(40)+.5*r_m2(41)+.7*r_m2(50)+.5*r_m2(51)
2145 d_m2(indx(ihcho_z))= 0.0
2146 !
2147 p_m2(indx(iald2_z))= .22*r_m2(11)+.47*r_m2(12)+1.03*r_m2(13) &
2148 +r_m2(14)+1.77*r_m2(15)+.03*r_m2(26)+.3*r_m2(27) &
2149 +.04*r_m2(28)+.3*r_m2(30)+.25*r_m2(33)+.5*r_m2(36) &
2150 +.3*r_m2(38)+.5*r_m2(41)+.21*r_m2(48)+.5*r_m2(51)
2151 d_m2(indx(iald2_z))= 0.0
2152 !
2153 p_m2(indx(ihcooh_z))= .52*r_m2(10)+.22*r_m2(12)
2154 d_m2(indx(ihcooh_z))= 0.0
2155 !
2156 p_m2(indx(iaone_z))= .07*r_m2(13)+.23*r_m2(15)+.74*r_m2(27) &
2157 +.74*r_m2(30)+.62*r_m2(33)+.74*r_m2(38) &
2158 +.57*r_m2(48)+.15*r_m2(50)
2159 d_m2(indx(iaone_z))= r_m2(5)+r_m2(6)
2160 !
2161 p_m2(indx(imgly_z))= .04*r_m2(12)+.07*r_m2(13)+.8*r_m2(19) &
2162 +.2*r_m2(26)+.19*r_m2(28)+.15*r_m2(50)
2163 d_m2(indx(imgly_z))= r_m2(7)+r_m2(8)+r_m2(9)
2164 !
2165 p_m2(indx(ieth_z))= 0.0
2166 d_m2(indx(ieth_z))= r_m2(10)+r_m2(11)
2167 !
2168 p_m2(indx(iolet_z))= 0.0
2169 d_m2(indx(iolet_z))= r_m2(12)+r_m2(14)+r_m2(16)
2170 !
2171 p_m2(indx(iolei_z))= 0.0
2172 d_m2(indx(iolei_z))= r_m2(13)+r_m2(15)+r_m2(17)
2173 !
2174 p_m2(indx(itol_z))= 0.0
2175 d_m2(indx(itol_z))= r_m2(18)
2176 !
2177 p_m2(indx(ixyl_z))= 0.0
2178 d_m2(indx(ixyl_z))= r_m2(19)
2179 !
2180 p_m2(indx(icres_z))= .12*r_m2(18)+.05*r_m2(19)
2181 d_m2(indx(icres_z))= r_m2(21)+r_m2(22)
2182 !
2183 p_m2(indx(ito2_z))= .8*r_m2(18)+.45*r_m2(19)
2184 d_m2(indx(ito2_z))= r_m2(20)
2185 !
2186 p_m2(indx(icro_z))= .4*r_m2(21)+r_m2(22)
2187 d_m2(indx(icro_z))= r_m2(23)
2188 !
2189 p_m2(indx(iopen_z))= .95*r_m2(20)+.3*r_m2(21)
2190 d_m2(indx(iopen_z))= r_m2(24)+r_m2(25)+r_m2(26)
2191 !
2192 p_m2(indx(ionit_z))= .05*r_m2(20)+r_m2(23)+.16*r_m2(33) &
2193 +.5*r_m2(36)+.5*r_m2(41)+r_m2(46)+.5*r_m2(51)
2194 d_m2(indx(ionit_z))= r_m2(29)+r_m2(30)
2195 !
2196 p_m2(indx(ipan_z))= 0.0
2197 d_m2(indx(ipan_z))= 0.0
2198 !
2199 p_m2(indx(ircooh_z))= .09*r_m2(12)+.16*r_m2(13)
2200 d_m2(indx(ircooh_z))= 0.0
2201 !
2202 p_m2(indx(irooh_z))= r_m2(43)+r_m2(45)
2203 d_m2(indx(irooh_z))= r_m2(27)+r_m2(28)
2204 !
2205 p_m2(indx(ich3o2_z))= +r_m2(5)+.07*r_m2(12)+.10*r_m2(13)
2206 d_m2(indx(ich3o2_z))= 0.0
2207 !
2208 p_m2(indx(iethp_z))= .06*r_m2(12)+.05*r_m2(13)+.1*r_m2(27) &
2209 +.1*r_m2(30)+.08*r_m2(33)+.1*r_m2(38)+.06*r_m2(48)
2210 d_m2(indx(iethp_z))= 0.0
2211 !
2212 p_m2(indx(ic2o3_z))= +r_m2(5)+r_m2(7)+r_m2(8) &
2213 +r_m2(9)+.13*r_m2(12)+.19*r_m2(13)+r_m2(24) &
2214 +r_m2(25)+.62*r_m2(26) +r_m2(35) &
2215 +r_m2(40)+.7*r_m2(50)
2216 d_m2(indx(ic2o3_z))= 0.0
2217 !
2218 p_m2(indx(iro2_z))= r_m2(1)+.03*r_m2(12)+.09*r_m2(13)+.77*r_m2(28)
2219 d_m2(indx(iro2_z))= r_m2(33)+r_m2(38)+r_m2(43)+r_m2(48)
2220 !
2221 p_m2(indx(iano2_z))= r_m2(6)+.11*r_m2(13)
2222 d_m2(indx(iano2_z))= r_m2(35)+r_m2(40)+r_m2(45)+r_m2(50)
2223 !
2224 p_m2(indx(inap_z))= r_m2(16)+r_m2(17)+r_m2(29)
2225 d_m2(indx(inap_z))= r_m2(36)+r_m2(41)+r_m2(46)+r_m2(51)
2226 !
2227 p_m2(indx(ixo2_z))= r_m2(8)+r_m2(11)+r_m2(14)+r_m2(15) &
2228 +.08*r_m2(18)+.5*r_m2(19)+.6*r_m2(21) &
2229 +r_m2(24)+.03*r_m2(26)+.4*r_m2(27)+.4*r_m2(30) &
2230 +.34*r_m2(33)+.4*r_m2(38)+.24*r_m2(48)
2231 d_m2(indx(ixo2_z))= r_m2(37)+r_m2(42)+r_m2(47)+r_m2(52)
2232
2233 p_m2(indx(ixpar_z))= 1.06*r_m2(12)+2.26*r_m2(13)+r_m2(14) &
2234 +2.23*r_m2(15)+1.98*r_m2(27)+.42*r_m2(28) &
2235 +1.98*r_m2(30)+1.68*r_m2(33)+r_m2(36) &
2236 +1.98*r_m2(38)+r_m2(41)+1.25*r_m2(48)+r_m2(51)
2237 d_m2(indx(ixpar_z))= r_m2(53)
2238 !
2239 return
2240 end subroutine gasode_m2
2241
2242
2243
2244 !***********************************************************************
2245 ! <18.> subr gasode_m3
2246 !
2247 ! purpose: updates production and destruction rates for mechanism 3
2248 ! isoprene
2249 !
2250 ! author : Rahul A. Zaveri
2251 ! date : December, 1998
2252 !
2253 ! ----------------------------------------------------------------------
2254
2255 subroutine gasode_m3( indx, r_m3, p_m3, d_m3 )
2256
2257 use module_data_cbmz
2258 implicit none
2259
2260 ! subr arguments
2261 integer indx(ngas_z)
2262 real r_m3(nrxn_m3), p_m3(ngas_tot), d_m3(ngas_tot)
2263
2264
2265 p_m3(indx(ino_z))= 0.0
2266 d_m3(indx(ino_z))= r_m3(8)+r_m3(9)+r_m3(10)
2267 !
2268 p_m3(indx(ino2_z))= .91*r_m3(8)+1.2*r_m3(9)+r_m3(10)
2269 d_m3(indx(ino2_z))= 0.0
2270 !
2271 p_m3(indx(ino3_z))= 0.0
2272 d_m3(indx(ino3_z))= r_m3(3)+r_m3(7)
2273 !
2274 p_m3(indx(ihno3_z))= .07*r_m3(7)
2275 d_m3(indx(ihno3_z))= 0.0
2276 !
2277 p_m3(indx(io3_z))= 0.0
2278 d_m3(indx(io3_z))= r_m3(2)+r_m3(6)
2279 !
2280 p_m3(indx(ioh_z))= .27*r_m3(2)+.27*r_m3(6)
2281 d_m3(indx(ioh_z))= r_m3(1)+r_m3(5)
2282 !
2283 p_m3(indx(iho2_z))= .07*r_m3(2)+.33*r_m3(4)+.1*r_m3(6)+.93*r_m3(7) &
2284 +.91*r_m3(8)+.8*r_m3(9)+r_m3(10)
2285 d_m3(indx(iho2_z))= r_m3(11)+r_m3(12)+r_m3(13)
2286 !
2287 p_m3(indx(ico_z))= .07*r_m3(2)+.33*r_m3(4)+.16*r_m3(6)+.64*r_m3(7) &
2288 +.59*r_m3(10)
2289 d_m3(indx(ico_z))= 0.0
2290 !
2291 p_m3(indx(ipar_z))= 1.86*r_m3(7)+0.18*r_m3(8)+1.6*r_m3(9)+2*r_m3(12) &
2292 +2*r_m3(15)
2293 d_m3(indx(ipar_z))= 0.0
2294 !
2295 p_m3(indx(ihcho_z))= .6*r_m3(2)+.2*r_m3(4)+.15*r_m3(6)+.28*r_m3(7) &
2296 +.63*r_m3(8)+.25*r_m3(10)
2297 d_m3(indx(ihcho_z))= 0.0
2298 !
2299 p_m3(indx(iald2_z))= .15*r_m3(2)+.07*r_m3(4)+.02*r_m3(6)+.28*r_m3(7) &
2300 +.8*r_m3(9)+.55*r_m3(10)+r_m3(15)+.5*r_m3(16)
2301 d_m3(indx(iald2_z))= 0.0
2302 !
2303 p_m3(indx(iaone_z))= .03*r_m3(4)+.09*r_m3(6)+.63*r_m3(10)+.5*r_m3(16)
2304 d_m3(indx(iaone_z))= 0.0
2305 !
2306 p_m3(indx(imgly_z))= .85*r_m3(6)+.34*r_m3(10)
2307 d_m3(indx(imgly_z))= 0.0
2308 !
2309 p_m3(indx(ionit_z))= .93*r_m3(7)+.09*r_m3(8)+.8*r_m3(9)+r_m3(12) &
2310 +r_m3(15)
2311 d_m3(indx(ionit_z))= 0.0
2312 !
2313 p_m3(indx(ihcooh_z))= .39*r_m3(2)+.46*r_m3(6)
2314 d_m3(indx(ihcooh_z))= 0.0
2315 !
2316 p_m3(indx(irooh_z))= r_m3(11)+r_m3(13)
2317 d_m3(indx(irooh_z))= 0.0
2318 !
2319 p_m3(indx(ich3o2_z))= .7*r_m3(4)+.05*r_m3(6)
2320 d_m3(indx(ich3o2_z))= 0.0
2321 !
2322 p_m3(indx(ic2o3_z))= .2*r_m3(2)+.97*r_m3(4)+.5*r_m3(5)+.11*r_m3(6) &
2323 +.07*r_m3(7)
2324 d_m3(indx(ic2o3_z))= 0.0
2325 !
2326 p_m3(indx(ixo2_z))= .08*r_m3(1)+.2*r_m3(2)+.2*r_m3(5)+.07*r_m3(6) &
2327 +.93*r_m3(7)
2328 d_m3(indx(ixo2_z))= 0.0
2329 !
2330 p_m3(indx(iisop_z))= 0.0
2331 d_m3(indx(iisop_z))= r_m3(1)+r_m3(2)+r_m3(3)
2332 !
2333 p_m3(indx(iisoprd_z))= .65*r_m3(2)+.91*r_m3(8)+.2*r_m3(9)+r_m3(14)
2334 d_m3(indx(iisoprd_z))= r_m3(4)+r_m3(5)+r_m3(6)+r_m3(7)
2335 !
2336 p_m3(indx(iisopp_z))= r_m3(1)
2337 d_m3(indx(iisopp_z))= r_m3(8)+r_m3(11)+r_m3(14)
2338 !
2339 p_m3(indx(iisopn_z))= r_m3(3)
2340 d_m3(indx(iisopn_z))= r_m3(9)+r_m3(12)+r_m3(15)
2341 !
2342 p_m3(indx(iisopo2_z))= .5*r_m3(5)
2343 d_m3(indx(iisopo2_z))= r_m3(10)+r_m3(13)+r_m3(16)
2344 !
2345 return
2346 end subroutine gasode_m3
2347
2348
2349
2350 !***********************************************************************
2351 ! <19.> subr gasode_m4
2352 !
2353 ! purpose: updates production and destruction rates for mechanism 4
2354 ! dimethylsulfide
2355 !
2356 ! author : Rahul A. Zaveri
2357 ! date : December, 1998
2358 !
2359 ! ----------------------------------------------------------------------
2360
2361 subroutine gasode_m4( indx, r_m4, p_m4, d_m4 )
2362
2363 use module_data_cbmz
2364 implicit none
2365
2366 ! subr arguments
2367 integer indx(ngas_z)
2368 real r_m4(nrxn_m4), p_m4(ngas_tot), d_m4(ngas_tot)
2369
2370
2371 p_m4(indx(ino_z))= r_m4(19)
2372 d_m4(indx(ino_z))= r_m4(5)+r_m4(11)+r_m4(26)+r_m4(30)
2373 !
2374 p_m4(indx(ino2_z))= r_m4(5)+r_m4(11)+r_m4(26)
2375 d_m4(indx(ino2_z))= r_m4(19)+r_m4(29)
2376 !
2377 p_m4(indx(ino3_z))= 0.0
2378 d_m4(indx(ino3_z))= r_m4(2)+r_m4(14)
2379 !
2380 p_m4(indx(ihono_z))= r_m4(30)
2381 d_m4(indx(ihono_z))= 0.0
2382 !
2383 p_m4(indx(ihno3_z))= r_m4(2)+r_m4(14)+r_m4(29)
2384 d_m4(indx(ihno3_z))= 0.0
2385 !
2386 p_m4(indx(io3_z))= 0.0
2387 d_m4(indx(io3_z))= r_m4(20)
2388 !
2389 p_m4(indx(io3p_z))= 0.0
2390 d_m4(indx(io3p_z))= r_m4(3)
2391 !
2392 p_m4(indx(ioh_z))= r_m4(21)
2393 d_m4(indx(ioh_z))= r_m4(1)+r_m4(4)+r_m4(9)+r_m4(10)+r_m4(16)+r_m4(23)
2394 !
2395 p_m4(indx(iho2_z))= .965*r_m4(4)+r_m4(6)+.27*r_m4(9)+r_m4(12)+r_m4(22) &
2396 +r_m4(27)+r_m4(32)
2397 d_m4(indx(iho2_z))= r_m4(13)+r_m4(21)+r_m4(31)
2398 !
2399 p_m4(indx(ih2o2_z))= r_m4(13)
2400 d_m4(indx(ih2o2_z))= 0.0
2401 !
2402 p_m4(indx(ico_z))= r_m4(32)
2403 d_m4(indx(ico_z))= 0.0
2404 !
2405 p_m4(indx(iso2_z))= r_m4(18)
2406 d_m4(indx(iso2_z))= 0.0
2407 !
2408 p_m4(indx(ih2so4_z))= r_m4(28)
2409 d_m4(indx(ih2so4_z))= 0.0
2410 !
2411 p_m4(indx(ihcho_z))= r_m4(5)+2*r_m4(6)+r_m4(7)+r_m4(11)+2*r_m4(12) &
2412 +r_m4(22)+r_m4(27)
2413 d_m4(indx(ihcho_z))= r_m4(32)
2414 !
2415 p_m4(indx(ich3o2_z))= r_m4(3)+.035*r_m4(4)+.73*r_m4(9)+r_m4(18)+r_m4(28)
2416 d_m4(indx(ich3o2_z))= r_m4(6)+r_m4(12)+r_m4(15)+r_m4(22)+r_m4(27)
2417 !
2418 p_m4(indx(ich3ooh_z))= r_m4(15)
2419 d_m4(indx(ich3ooh_z))= 0.0
2420 !
2421 p_m4(indx(idms_z))= 0.0
2422 d_m4(indx(idms_z))= r_m4(1)+r_m4(2)+r_m4(3)+r_m4(4)
2423 !
2424 p_m4(indx(imsa_z))= r_m4(17)+r_m4(23)+r_m4(29)+r_m4(30)+r_m4(31)+r_m4(32)
2425 d_m4(indx(imsa_z))= 0.0
2426 !
2427 p_m4(indx(idmso_z))= .965*r_m4(4)
2428 d_m4(indx(idmso_z))= r_m4(9)
2429 !
2430 p_m4(indx(idmso2_z))= .27*r_m4(9)
2431 d_m4(indx(idmso2_z))= r_m4(10)
2432 !
2433 p_m4(indx(ich3so2h_z))= .73*r_m4(9)
2434 d_m4(indx(ich3so2h_z))= r_m4(13)+r_m4(14)+r_m4(15)+r_m4(16)+r_m4(17)
2435 !
2436 p_m4(indx(ich3sch2oo_z))= r_m4(1)+r_m4(2)
2437 d_m4(indx(ich3sch2oo_z))= r_m4(5)+r_m4(6)+r_m4(7)+r_m4(8)+r_m4(8)
2438 !
2439 p_m4(indx(ich3so2_z))= r_m4(3)+.035*r_m4(4)+r_m4(5)+r_m4(6)+r_m4(7) &
2440 +1.85*r_m4(8) &
2441 +r_m4(11)+r_m4(12)+r_m4(13)+r_m4(14)+r_m4(15) &
2442 +r_m4(16)+r_m4(17)+r_m4(25)
2443 d_m4(indx(ich3so2_z))= r_m4(7)+r_m4(18)+r_m4(19)+r_m4(20)+r_m4(21) &
2444 +r_m4(22)+r_m4(23)+r_m4(24)
2445 !
2446 p_m4(indx(ich3so3_z))= r_m4(7)+r_m4(19)+r_m4(20)+r_m4(21)+r_m4(22) &
2447 +r_m4(26)+r_m4(27)
2448 d_m4(indx(ich3so3_z))= r_m4(17)+r_m4(28)+r_m4(29)+r_m4(30)+r_m4(31) &
2449 +r_m4(32)
2450 !
2451 p_m4(indx(ich3so2oo_z))= r_m4(24)
2452 d_m4(indx(ich3so2oo_z))= r_m4(25)+r_m4(26)+r_m4(27)
2453 !
2454 p_m4(indx(ich3so2ch2oo_z))= r_m4(10)
2455 d_m4(indx(ich3so2ch2oo_z))= r_m4(11)+r_m4(12)
2456 !
2457 p_m4(indx(imtf_z))= .15*r_m4(8)
2458 d_m4(indx(imtf_z))= 0.0
2459 !
2460 return
2461 end subroutine gasode_m4
2462
2463
2464
2465 !***********************************************************************
2466 ! <20.> subr gasrate_m1
2467 !
2468 ! purpose: computes reaction rates for mechanism 1
2469 !
2470 ! author : Rahul A. Zaveri
2471 ! date : December 1998
2472 !
2473 !-------------------------------------------------------------------------
2474
2475 subroutine gasrate_m1( indx, s, r_m1, r_m2, rk_m1, rk_m2, &
2476 h2o, ch4, oxygen, nitrogen, hydrogen )
2477
2478 use module_data_cbmz
2479 implicit none
2480
2481 ! subr arguments
2482 integer indx(ngas_z)
2483 real s(ngas_tot), r_m1(nrxn_m1), r_m2(nrxn_m2)
2484 real rk_m1(nrxn_m1), rk_m2(nrxn_m2)
2485 real h2o, ch4, oxygen, nitrogen, hydrogen
2486
2487 r_m1(1) = rk_m1(1)*s(indx(ino2_z))
2488 r_m1(2) = rk_m1(2)*s(indx(ino3_z))
2489 r_m1(3) = rk_m1(3)*s(indx(ihono_z))
2490 r_m1(4) = rk_m1(4)*s(indx(ihno3_z))
2491 r_m1(5) = rk_m1(5)*s(indx(ihno4_z))
2492 r_m1(6) = rk_m1(6)*s(indx(in2o5_z))
2493 r_m1(7) = rk_m1(7)*s(indx(io3_z))
2494 r_m1(8) = rk_m1(8)*s(indx(io3_z))
2495 r_m1(9) = rk_m1(9)*s(indx(ih2o2_z))
2496 r_m1(10) = rk_m1(10)*s(indx(io1d_z))*oxygen
2497 r_m1(11) = rk_m1(11)*s(indx(io1d_z))*nitrogen
2498 r_m1(12) = rk_m1(12)*s(indx(io1d_z))*h2o
2499 r_m1(13) = rk_m1(13)*s(indx(io3p_z))*oxygen
2500 r_m1(14) = rk_m1(14)*s(indx(io3p_z))*s(indx(io3_z))
2501 r_m1(15) = rk_m1(15)*s(indx(io3p_z))*s(indx(ino2_z))
2502 r_m1(16) = rk_m1(16)*s(indx(io3p_z))*s(indx(ino2_z))
2503 r_m1(17) = rk_m1(17)*s(indx(io3p_z))*s(indx(ino_z))
2504 r_m1(18) = rk_m1(18)*s(indx(io3_z))*s(indx(ino_z))
2505 r_m1(19) = rk_m1(19)*s(indx(io3_z))*s(indx(ino2_z))
2506 r_m1(20) = rk_m1(20)*s(indx(io3_z))*s(indx(ioh_z))
2507 r_m1(21) = rk_m1(21)*s(indx(io3_z))*s(indx(iho2_z))
2508 r_m1(22) = rk_m1(22)*s(indx(ioh_z))*hydrogen
2509 r_m1(23) = rk_m1(23)*s(indx(ioh_z))*s(indx(ino_z))
2510 r_m1(24) = rk_m1(24)*s(indx(ioh_z))*s(indx(ino2_z))
2511 r_m1(25) = rk_m1(25)*s(indx(ioh_z))*s(indx(ino3_z))
2512 r_m1(26) = rk_m1(26)*s(indx(ioh_z))*s(indx(ihono_z))
2513 r_m1(27) = rk_m1(27)*s(indx(ioh_z))*s(indx(ihno3_z))
2514 r_m1(28) = rk_m1(28)*s(indx(ioh_z))*s(indx(ihno4_z))
2515 r_m1(29) = rk_m1(29)*s(indx(ioh_z))*s(indx(iho2_z))
2516 r_m1(30) = rk_m1(30)*s(indx(ioh_z))*s(indx(ih2o2_z))
2517 r_m1(31) = rk_m1(31)*s(indx(iho2_z))*s(indx(iho2_z))
2518 r_m1(32) = rk_m1(32)*s(indx(iho2_z))*s(indx(iho2_z))*h2o
2519 r_m1(33) = rk_m1(33)*s(indx(iho2_z))*s(indx(ino_z))
2520 r_m1(34) = rk_m1(34)*s(indx(iho2_z))*s(indx(ino2_z))
2521 r_m1(35) = rk_m1(35)*s(indx(iho2_z))*s(indx(ino2_z))
2522 r_m1(36) = rk_m1(36)*s(indx(ihno4_z))
2523 r_m1(37) = rk_m1(37)*s(indx(ino3_z))*s(indx(ino_z))
2524 r_m1(38) = rk_m1(38)*s(indx(ino3_z))*s(indx(ino2_z))
2525 r_m1(39) = rk_m1(39)*s(indx(ino3_z))*s(indx(ino2_z))
2526 r_m1(40) = rk_m1(40)*s(indx(ino3_z))*s(indx(ino3_z))
2527 r_m1(41) = rk_m1(41)*s(indx(ino3_z))*s(indx(iho2_z))
2528 r_m1(42) = rk_m1(42)*s(indx(in2o5_z))*h2o
2529 r_m1(43) = rk_m1(43)*s(indx(in2o5_z))
2530 r_m1(44) = rk_m1(44)*s(indx(ico_z))*s(indx(ioh_z))
2531 r_m1(45) = rk_m1(45)*s(indx(iso2_z))*s(indx(ioh_z))
2532 r_m1(46) = rk_m1(46)*s(indx(ioh_z))*ch4
2533 r_m1(47) = rk_m1(47)*s(indx(ic2h6_z))*s(indx(ioh_z))
2534 r_m1(48) = rk_m1(48)*s(indx(ich3oh_z))*s(indx(ioh_z))
2535 r_m1(49) = rk_m1(49)*s(indx(ihcho_z))
2536 r_m1(50) = rk_m1(50)*s(indx(ihcho_z))
2537 r_m1(51) = rk_m1(51)*s(indx(ihcho_z))*s(indx(ioh_z))
2538 r_m1(52) = rk_m1(52)*s(indx(ihcho_z))*s(indx(ino3_z))
2539 r_m1(53) = rk_m1(53)*s(indx(ich3ooh_z))
2540 r_m1(54) = rk_m1(54)*s(indx(iethooh_z))
2541 r_m1(55) = rk_m1(55)*s(indx(ich3ooh_z))*s(indx(ioh_z))
2542 r_m1(56) = rk_m1(56)*s(indx(iethooh_z))*s(indx(ioh_z))
2543 r_m1(57) = rk_m1(57)*s(indx(ich3o2_z))*s(indx(ino_z))
2544 r_m1(58) = rk_m1(58)*s(indx(iethp_z))*s(indx(ino_z))
2545 r_m1(59) = rk_m1(59)*s(indx(ich3o2_z))*s(indx(ino3_z))
2546 r_m1(60) = rk_m1(60)*s(indx(iethp_z))*s(indx(ino3_z))
2547 r_m1(61) = rk_m1(61)*s(indx(ich3o2_z))*s(indx(iho2_z))
2548 r_m1(62) = rk_m1(62)*s(indx(iethp_z))*s(indx(iho2_z))
2549 r_m1(63) = rk_m1(63)*s(indx(ich3o2_z))
2550 r_m1(64) = rk_m1(64)*s(indx(iethp_z))
2551 r_m1(65) = rk_m1(65)*s(indx(ic2h5oh_z))*s(indx(ioh_z))
2552
2553 r_m2(2) = rk_m2(2)*s(indx(iald2_z))
2554 r_m2(3) = rk_m2(3)*s(indx(iald2_z))*s(indx(ioh_z))
2555 r_m2(4) = rk_m2(4)*s(indx(iald2_z))*s(indx(ino3_z))
2556 r_m2(31) = rk_m2(31)*s(indx(ic2o3_z))*s(indx(ino2_z))
2557 r_m2(32) = rk_m2(32)*s(indx(ipan_z))
2558 r_m2(34) = rk_m2(34)*s(indx(ic2o3_z))*s(indx(ino_z))
2559 r_m2(39) = rk_m2(39)*s(indx(ic2o3_z))*s(indx(ino3_z))
2560 r_m2(44) = rk_m2(44)*s(indx(ic2o3_z))*s(indx(iho2_z))
2561 r_m2(49) = rk_m2(49)*s(indx(ic2o3_z))
2562
2563 return
2564 end subroutine gasrate_m1
2565
2566
2567
2568 !***********************************************************************
2569 ! <21.> subr gasrate_m2
2570 !
2571 ! purpose: computes reaction rates for mechanism 2
2572 !
2573 ! author : Rahul A. Zaveri
2574 ! date : December 1998
2575 !
2576 !-------------------------------------------------------------------------
2577 !
2578 subroutine gasrate_m2( indx, s, r_m2, rk_m2 )
2579
2580 use module_data_cbmz
2581 implicit none
2582
2583 ! subr arguments
2584 integer indx(ngas_z)
2585 real s(ngas_tot), r_m2(nrxn_m2), rk_m2(nrxn_m2)
2586
2587 r_m2(1) = rk_m2(1)*s(indx(ipar_z))*s(indx(ioh_z))
2588
2589 r_m2(5) = rk_m2(5)*s(indx(iaone_z))
2590 r_m2(6) = rk_m2(6)*s(indx(iaone_z))*s(indx(ioh_z))
2591 r_m2(7) = rk_m2(7)*s(indx(imgly_z))
2592 r_m2(8) = rk_m2(8)*s(indx(imgly_z))*s(indx(ioh_z))
2593 r_m2(9) = rk_m2(9)*s(indx(imgly_z))*s(indx(ino3_z))
2594 r_m2(10) = rk_m2(10)*s(indx(ieth_z))*s(indx(io3_z))
2595 r_m2(11) = rk_m2(11)*s(indx(ieth_z))*s(indx(ioh_z))
2596 r_m2(12) = rk_m2(12)*s(indx(iolet_z))*s(indx(io3_z))
2597 r_m2(13) = rk_m2(13)*s(indx(iolei_z))*s(indx(io3_z))
2598 r_m2(14) = rk_m2(14)*s(indx(iolet_z))*s(indx(ioh_z))
2599 r_m2(15) = rk_m2(15)*s(indx(iolei_z))*s(indx(ioh_z))
2600 r_m2(16) = rk_m2(16)*s(indx(iolet_z))*s(indx(ino3_z))
2601 r_m2(17) = rk_m2(17)*s(indx(iolei_z))*s(indx(ino3_z))
2602 r_m2(18) = rk_m2(18)*s(indx(itol_z))*s(indx(ioh_z))
2603 r_m2(19) = rk_m2(19)*s(indx(ixyl_z))*s(indx(ioh_z))
2604 r_m2(20) = rk_m2(20)*s(indx(ito2_z))*s(indx(ino_z))
2605 r_m2(21) = rk_m2(21)*s(indx(icres_z))*s(indx(ioh_z))
2606 r_m2(22) = rk_m2(22)*s(indx(icres_z))*s(indx(ino3_z))
2607 r_m2(23) = rk_m2(23)*s(indx(icro_z))*s(indx(ino2_z))
2608 r_m2(24) = rk_m2(24)*s(indx(iopen_z))*s(indx(ioh_z))
2609 r_m2(25) = rk_m2(25)*s(indx(iopen_z))
2610 r_m2(26) = rk_m2(26)*s(indx(iopen_z))*s(indx(io3_z))
2611 r_m2(27) = rk_m2(27)*s(indx(irooh_z))
2612 r_m2(28) = rk_m2(28)*s(indx(irooh_z))*s(indx(ioh_z))
2613 r_m2(29) = rk_m2(29)*s(indx(ionit_z))*s(indx(ioh_z))
2614 r_m2(30) = rk_m2(30)*s(indx(ionit_z))
2615
2616 r_m2(33) = rk_m2(33)*s(indx(iro2_z))*s(indx(ino_z))
2617
2618 r_m2(35) = rk_m2(35)*s(indx(iano2_z))*s(indx(ino_z))
2619 r_m2(36) = rk_m2(36)*s(indx(inap_z))*s(indx(ino_z))
2620 r_m2(37) = rk_m2(37)*s(indx(ixo2_z))*s(indx(ino_z))
2621 r_m2(38) = rk_m2(38)*s(indx(iro2_z))*s(indx(ino3_z))
2622
2623 r_m2(40) = rk_m2(40)*s(indx(iano2_z))*s(indx(ino3_z))
2624 r_m2(41) = rk_m2(41)*s(indx(inap_z))*s(indx(ino3_z))
2625 r_m2(42) = rk_m2(42)*s(indx(ixo2_z))*s(indx(ino3_z))
2626 r_m2(43) = rk_m2(43)*s(indx(iro2_z))*s(indx(iho2_z))
2627
2628 r_m2(45) = rk_m2(45)*s(indx(iano2_z))*s(indx(iho2_z))
2629 r_m2(46) = rk_m2(46)*s(indx(inap_z))*s(indx(iho2_z))
2630 r_m2(47) = rk_m2(47)*s(indx(ixo2_z))*s(indx(iho2_z))
2631 r_m2(48) = rk_m2(48)*s(indx(iro2_z))
2632
2633 r_m2(50) = rk_m2(50)*s(indx(iano2_z))
2634 r_m2(51) = rk_m2(51)*s(indx(inap_z))
2635 r_m2(52) = rk_m2(52)*s(indx(ixo2_z))
2636 r_m2(53) = rk_m2(53)*s(indx(ipar_z))*s(indx(ixpar_z))
2637 !
2638 return
2639 end subroutine gasrate_m2
2640
2641
2642
2643 !***********************************************************************
2644 ! <22.> subr gasrate_m3
2645 !
2646 ! purpose: computes reaction rates for mechanism 3
2647 !
2648 ! author : Rahul A. Zaveri
2649 ! date : December 1998
2650 !
2651 !-------------------------------------------------------------------------
2652 !
2653 subroutine gasrate_m3( indx, s, r_m3, rk_m3 )
2654
2655 use module_data_cbmz
2656 implicit none
2657
2658 ! subr arguments
2659 integer indx(ngas_z)
2660 real s(ngas_tot), r_m3(nrxn_m3), rk_m3(nrxn_m3)
2661
2662 r_m3(1) = rk_m3(1)*s(indx(iisop_z))*s(indx(ioh_z))
2663 r_m3(2) = rk_m3(2)*s(indx(iisop_z))*s(indx(io3_z))
2664 r_m3(3) = rk_m3(3)*s(indx(iisop_z))*s(indx(ino3_z))
2665 r_m3(4) = rk_m3(4)*s(indx(iisoprd_z))
2666 r_m3(5) = rk_m3(5)*s(indx(iisoprd_z))*s(indx(ioh_z))
2667 r_m3(6) = rk_m3(6)*s(indx(iisoprd_z))*s(indx(io3_z))
2668 r_m3(7) = rk_m3(7)*s(indx(iisoprd_z))*s(indx(ino3_z))
2669 r_m3(8) = rk_m3(8)*s(indx(iisopp_z))*s(indx(ino_z))
2670 r_m3(9) = rk_m3(9)*s(indx(iisopn_z))*s(indx(ino_z))
2671 r_m3(10) = rk_m3(10)*s(indx(iisopo2_z))*s(indx(ino_z))
2672 r_m3(11) = rk_m3(11)*s(indx(iisopp_z))*s(indx(iho2_z))
2673 r_m3(12) = rk_m3(12)*s(indx(iisopn_z))*s(indx(iho2_z))
2674 r_m3(13) = rk_m3(13)*s(indx(iisopo2_z))*s(indx(iho2_z))
2675 r_m3(14) = rk_m3(14)*s(indx(iisopp_z))
2676 r_m3(15) = rk_m3(15)*s(indx(iisopn_z))
2677 r_m3(16) = rk_m3(16)*s(indx(iisopo2_z))
2678
2679 return
2680 end subroutine gasrate_m3
2681
2682
2683
2684 !***********************************************************************
2685 ! <23.> subr gasrate_m4
2686 !
2687 ! purpose: computes reaction rates for mechanism 4
2688 !
2689 ! author : Rahul A. Zaveri
2690 ! date : December 1998
2691 !
2692 !-------------------------------------------------------------------------
2693 !
2694 subroutine gasrate_m4( indx, s, r_m4, rk_m4, oxygen )
2695
2696 use module_data_cbmz
2697 implicit none
2698
2699 ! subr arguments
2700 integer indx(ngas_z)
2701 real s(ngas_tot), r_m4(nrxn_m4), rk_m4(nrxn_m4)
2702 real oxygen
2703
2704 r_m4(1) = rk_m4(1)*s(indx(idms_z))*s(indx(ioh_z))
2705 r_m4(2) = rk_m4(2)*s(indx(idms_z))*s(indx(ino3_z))
2706 r_m4(3) = rk_m4(3)*s(indx(idms_z))*s(indx(io3p_z))
2707 r_m4(4) = rk_m4(4)*s(indx(idms_z))*s(indx(ioh_z))
2708 r_m4(5) = rk_m4(5)*s(indx(ich3sch2oo_z))*s(indx(ino_z))
2709 r_m4(6) = rk_m4(6)*s(indx(ich3sch2oo_z))*s(indx(ich3o2_z))
2710 r_m4(7) = rk_m4(7)*s(indx(ich3sch2oo_z))*s(indx(ich3so2_z))
2711 r_m4(8) = rk_m4(8)*s(indx(ich3sch2oo_z))*s(indx(ich3sch2oo_z))
2712 r_m4(9) = rk_m4(9)*s(indx(idmso_z))*s(indx(ioh_z))
2713 r_m4(10) = rk_m4(10)*s(indx(idmso2_z))*s(indx(ioh_z))
2714 r_m4(11) = rk_m4(11)*s(indx(ich3so2ch2oo_z))*s(indx(ino_z))
2715 r_m4(12) = rk_m4(12)*s(indx(ich3so2ch2oo_z))*s(indx(ich3o2_z))
2716 r_m4(13) = rk_m4(13)*s(indx(ich3so2h_z))*s(indx(iho2_z))
2717 r_m4(14) = rk_m4(14)*s(indx(ich3so2h_z))*s(indx(ino3_z))
2718 r_m4(15) = rk_m4(15)*s(indx(ich3so2h_z))*s(indx(ich3o2_z))
2719 r_m4(16) = rk_m4(16)*s(indx(ich3so2h_z))*s(indx(ioh_z))
2720 r_m4(17) = rk_m4(17)*s(indx(ich3so2h_z))*s(indx(ich3so3_z))
2721 r_m4(18) = rk_m4(18)*s(indx(ich3so2_z))
2722 r_m4(19) = rk_m4(19)*s(indx(ich3so2_z))*s(indx(ino2_z))
2723 r_m4(20) = rk_m4(20)*s(indx(ich3so2_z))*s(indx(io3_z))
2724 r_m4(21) = rk_m4(21)*s(indx(ich3so2_z))*s(indx(iho2_z))
2725 r_m4(22) = rk_m4(22)*s(indx(ich3so2_z))*s(indx(ich3o2_z))
2726 r_m4(23) = rk_m4(23)*s(indx(ich3so2_z))*s(indx(ioh_z))
2727 r_m4(24) = rk_m4(24)*s(indx(ich3so2_z))*oxygen
2728 r_m4(25) = rk_m4(25)*s(indx(ich3so2oo_z))
2729 r_m4(26) = rk_m4(26)*s(indx(ich3so2oo_z))*s(indx(ino_z))
2730 r_m4(27) = rk_m4(27)*s(indx(ich3so2oo_z))*s(indx(ich3o2_z))
2731 r_m4(28) = rk_m4(28)*s(indx(ich3so3_z))
2732 r_m4(29) = rk_m4(29)*s(indx(ich3so3_z))*s(indx(ino2_z))
2733 r_m4(30) = rk_m4(30)*s(indx(ich3so3_z))*s(indx(ino_z))
2734 r_m4(31) = rk_m4(31)*s(indx(ich3so3_z))*s(indx(iho2_z))
2735 r_m4(32) = rk_m4(32)*s(indx(ich3so3_z))*s(indx(ihcho_z))
2736
2737 return
2738 end subroutine gasrate_m4
2739
2740
2741
2742 !**************************************************************************
2743 ! <24.> subr loadperoxyparameters
2744 !
2745 ! purpose: loads thermal rate coefficients for peroxy-peroxy
2746 ! permutation reactions
2747 !
2748 ! author : Rahul A. Zaveri
2749 ! date : June 1998
2750 !
2751 ! nomenclature:
2752 ! Aperox = Pre-exponential factor (molec-cc-s)
2753 ! Bperox = activation energy (-E/R) (K)
2754 !
2755 !-------------------------------------------------------------------------
2756 subroutine loadperoxyparameters( Aperox, Bperox )
2757
2758 use module_data_cbmz
2759 implicit none
2760
2761 ! subr arguments
2762 real Aperox(nperox,nperox), Bperox(nperox,nperox)
2763
2764 ! local variables
2765 integer i, j
2766
2767 Aperox(jch3o2,jch3o2) = 2.5e-13
2768 Aperox(jethp,jethp) = 6.8e-14
2769 Aperox(jc2o3,jc2o3) = 2.9e-12
2770 Aperox(jano2,jano2) = 8.0e-12
2771 Aperox(jnap,jnap) = 1.0e-12
2772 Aperox(jro2,jro2) = 5.3e-16
2773 Aperox(jisopp,jisopp) = 3.1e-14
2774 Aperox(jisopn,jisopn) = 3.1e-14
2775 Aperox(jisopo2,jisopo2) = 3.1e-14
2776 Aperox(jxo2,jxo2) = 3.1e-14
2777
2778 Bperox(jch3o2,jch3o2) = 190.
2779 Bperox(jethp,jethp) = 0.0
2780 Bperox(jc2o3,jc2o3) = 500.
2781 Bperox(jano2,jano2) = 0.0
2782 Bperox(jnap,jnap) = 0.0
2783 Bperox(jro2,jro2) = 1980.
2784 Bperox(jisopp,jisopp) = 1000.
2785 Bperox(jisopn,jisopn) = 1000.
2786 Bperox(jisopo2,jisopo2) = 1000.
2787 Bperox(jxo2,jxo2) = 1000.
2788
2789 do i = 1, nperox
2790 do j = 1, nperox
2791 if(i.ne.j)then
2792 Aperox(i,j) = 2.0*sqrt(Aperox(i,i)*Aperox(j,j))
2793 Bperox(i,j) = 0.5*(Bperox(i,i) + Bperox(j,j))
2794 endif
2795 enddo
2796 enddo
2797
2798 ! except for
2799 Aperox(jc2o3,jch3o2) = 1.3e-12
2800 Aperox(jch3o2,jc2o3) = 1.3e-12
2801 Bperox(jc2o3,jch3o2) = 640.
2802 Bperox(jch3o2,jc2o3) = 640.
2803
2804 return
2805 end subroutine loadperoxyparameters
2806
2807
2808
2809
2810 !**************************************************************************
2811 ! <25.> subr peroxyrateconstants
2812 !
2813 ! purpose: computes parameterized thermal rate coefficients
2814 ! for the alkylperoxy radical permutation reactions
2815 ! for the entire mechanism.
2816 !
2817 ! author : Rahul A. Zaveri
2818 ! date : June 1998
2819 !
2820 ! nomenclature:
2821 ! rk_param = parameterized reaction rate constants (1/s)
2822 ! rk_perox = individual permutation reaction rate constants (molec-cc-s)
2823 ! te = ambient atmospheric temperature (K)
2824 !
2825 !-------------------------------------------------------------------------
2826 subroutine peroxyrateconstants( tempbox, cbox, &
2827 Aperox, Bperox, rk_param )
2828
2829 use module_data_cbmz
2830 implicit none
2831
2832 ! subr arguments
2833 real tempbox, cbox(ngas_z)
2834 real Aperox(nperox,nperox), Bperox(nperox,nperox), rk_param(nperox)
2835
2836 ! local variables
2837 integer i, j
2838 real te
2839 real sperox(nperox), rk_perox(nperox,nperox)
2840
2841
2842 te = tempbox
2843
2844 sperox(jch3o2) = cbox(ich3o2_z)
2845 sperox(jethp) = cbox(iethp_z)
2846 sperox(jro2) = cbox(iro2_z)
2847 sperox(jc2o3) = cbox(ic2o3_z)
2848 sperox(jano2) = cbox(iano2_z)
2849 sperox(jnap) = cbox(inap_z)
2850 sperox(jisopp) = cbox(iisopp_z)
2851 sperox(jisopn) = cbox(iisopn_z)
2852 sperox(jisopo2) = cbox(iisopo2_z)
2853 sperox(jxo2) = cbox(ixo2_z)
2854
2855 !
2856 ! initialize to zero
2857 do i = 1, nperox
2858 rk_param(i) = 0.0
2859 enddo
2860
2861 do i = 1, nperox
2862 do j = 1, nperox
2863 rk_perox(i,j) = arr( Aperox(i,j), Bperox(i,j), te )
2864 rk_param(i) = rk_param(i) + rk_perox(i,j)*sperox(j)
2865 enddo
2866 enddo
2867
2868 return
2869 end subroutine peroxyrateconstants
2870
2871
2872
2873 !***********************************************************************
2874 ! <26.> subr gasthermrk_m1
2875 !
2876 ! purpose: computes thermal reaction rate coefficients for
2877 ! mechanism 1
2878 !
2879 ! author : Rahul A. Zaveri
2880 ! date : December 1998
2881 !
2882 !-------------------------------------------------------------------------
2883
2884 subroutine gasthermrk_m1( tempbox, cair_mlc, &
2885 rk_photo, rk_param, rk_m1, rk_m2 )
2886
2887 use module_data_cbmz
2888 implicit none
2889
2890 ! subr arguments
2891 real tempbox, cair_mlc
2892 real rk_photo(nphoto), rk_param(nperox)
2893 real rk_m1(nrxn_m1), rk_m2(nrxn_m2)
2894 ! local variables
2895 integer i
2896 real rk0, rk2, rk3, rki, rko, rmm, rnn, te
2897 ! real arr, troe
2898
2899
2900 te = tempbox
2901
2902 rk_m1(1) = rk_photo(jphoto_no2)
2903 rk_m1(2) = rk_photo(jphoto_no3)
2904 rk_m1(3) = rk_photo(jphoto_hono)
2905 rk_m1(4) = rk_photo(jphoto_hno3)
2906 rk_m1(5) = rk_photo(jphoto_hno4)
2907 rk_m1(6) = rk_photo(jphoto_n2o5)
2908 rk_m1(7) = rk_photo(jphoto_o3a)
2909 rk_m1(8) = rk_photo(jphoto_o3b)
2910 rk_m1(9) = rk_photo(jphoto_h2o2)
2911 rk_m1(10) = arr(3.2e-11, 70., te)
2912 rk_m1(11) = arr(1.8e-11, 110., te)
2913 rk_m1(12) = 2.2e-10
2914 rk_m1(13) = cair_mlc*6.e-34*(te/300.)**(-2.3)
2915 rk_m1(14) = arr(8.0e-12, -2060., te)
2916 rk_m1(15) = arr(6.5e-12, -120., te)
2917 !
2918 rk0 = 9.0e-32
2919 rnn = 2.0
2920 rki = 2.2e-11
2921 rmm = 0.0
2922 rk_m1(16) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2923 !
2924 rk0 = 9.0e-32
2925 rnn = 1.5
2926 rki = 3.0e-11
2927 rmm = 0.0
2928 rk_m1(17) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2929 !
2930 rk_m1(18) = arr(2.0e-12, -1400., te)
2931 rk_m1(19) = arr(1.2e-13, -2450., te)
2932 rk_m1(20) = arr(1.6e-12, -940., te)
2933 rk_m1(21) = arr(1.1e-14, -500., te)
2934 rk_m1(22) = arr(5.5e-12, -2000., te)
2935 !
2936 rk0 = 7.0e-31
2937 rnn = 2.6
2938 rki = 3.6e-11
2939 rmm = 0.1
2940 rk_m1(23) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2941 !
2942 rk0 = 2.5e-30
2943 rnn = 4.4
2944 rki = 1.6e-11
2945 rmm = 1.7
2946 rk_m1(24) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2947 !
2948 rk_m1(25) = 2.2e-11
2949 rk_m1(26) = arr(1.8e-11, -390., te)
2950 rko = 7.2e-15 * exp(785./te)
2951 rk2 = 4.1e-16 * exp(1440./te)
2952 rk3 = 1.9e-33 * exp(725./te)*cair_mlc
2953 rk_m1(27) = rko + rk3/(1.+rk3/rk2)
2954 rk_m1(28) = arr(1.3e-12, 380., te)
2955 rk_m1(29) = arr(4.8e-11, 250., te)
2956 rk_m1(30) = arr(2.9e-12, -160., te)
2957 rk_m1(31) = 2.3e-13 * exp(600./te) + &
2958 1.7e-33 * exp(1000./te)*cair_mlc ! ho2 + ho2 --> h2o2
2959 rk_m1(32) = rk_m1(31)*1.4e-21*exp(2200./te) ! ho2 + ho2 + h2o --> h2o2
2960 rk_m1(33) = arr(3.5e-12, 250., te)
2961 !
2962 rk0 = 1.8e-31
2963 rnn = 3.2
2964 rki = 4.7e-12
2965 rmm = 1.4
2966 rk_m1(34) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2967 !
2968 rk_m1(35) = 5.0e-16
2969 rk_m1(36) = rk_m1(34)*arr(4.8e26, -10900., te)
2970 rk_m1(37) = arr(1.5e-11, 170., te)
2971 rk_m1(38) = arr(4.5e-14, -1260., te)
2972 !
2973 rk0 = 2.2e-30
2974 rnn = 3.9
2975 rki = 1.5e-12
2976 rmm = 0.7
2977 rk_m1(39) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2978 !
2979 rk_m1(40) = arr(8.5e-13, -2450., te)
2980 rk_m1(41) = 3.5e-12
2981 rk_m1(42) = 2.0e-21
2982 rk_m1(43) = rk_m1(39)*arr(3.7e26, -11000., te)
2983 rk_m1(44) = 1.5e-13 * (1.+8.18e-23*te*cair_mlc) ! co + oh --> ho2
2984
2985 rk0 = 3.0e-31
2986 rnn = 3.3
2987 rki = 1.5e-12
2988 rmm = 0.0
2989 rk_m1(45) = Troe(cair_mlc,te,rk0,rnn,rki,rmm)
2990
2991 rk_m1(46) = te**.667*arr(2.8e-14, -1575., te)
2992 rk_m1(47) = te**2*arr(1.5e-17, -492., te)
2993 rk_m1(48) = arr(6.7e-12, -600., te)
2994 rk_m1(49) = rk_photo(jphoto_hchoa) ! hcho + hv --> 2ho2 + co
2995 rk_m1(50) = rk_photo(jphoto_hchob) ! hcho + hv --> co
2996 rk_m1(51) = 1.0e-11
2997 rk_m1(52) = arr(3.4e-13, -1900., te)
2998 rk_m1(53) = rk_photo(jphoto_ch3ooh)
2999 rk_m1(54) = rk_photo(jphoto_ethooh)
3000 rk_m1(55) = arr(3.8e-12, 200., te)
3001 rk_m1(56) = arr(3.8e-12, 200., te)
3002 rk_m1(57) = arr(3.0e-12, 280., te)
3003 rk_m1(58) = arr(2.6e-12, 365., te)
3004 rk_m1(59) = 1.1e-12
3005 rk_m1(60) = 2.5e-12
3006 rk_m1(61) = arr(3.8e-13, 800., te)
3007 rk_m1(62) = arr(7.5e-13, 700., te)
3008 rk_m1(63) = rk_param(jch3o2)
3009 rk_m1(64) = rk_param(jethp)
3010 rk_m1(65) = arr(7.0e-12, -235.,te)
3011
3012 rk_m2(2) = rk_photo(jphoto_ald2)
3013 rk_m2(3) = arr(5.6e-12, 270., te)
3014 rk_m2(4) = arr(1.4e-12, -1900., te)
3015 !
3016 rk0 = 9.7e-29
3017 rnn = 5.6
3018 rki = 9.3e-12
3019 rmm = 1.5
3020 rk_m2(31) = troe(cair_mlc,te,rk0,rnn,rki,rmm)
3021 !
3022 rk_m2(32) = rk_m2(31)*arr(1.1e28, -14000., te)
3023 rk_m2(34) = arr(5.3e-12, 360., te)
3024 rk_m2(39) = 4.0e-12
3025 rk_m2(44) = arr(4.5e-13, 1000., te)
3026 rk_m2(49) = rk_param(jc2o3)
3027
3028 ! Heterogeneous reactions
3029 ! rk_m1(65) = rk_het(1) ! O3 -->
3030 ! rk_m1(66) = rk_het(2) ! HO2 --> 0.5H2O2
3031 ! rk_m1(67) = rk_het(3) ! NO2 --> 0.5HONO + 0.5HNO3
3032 ! rk_m1(68) = rk_het(4) ! N2O5 --> 2HNO3
3033 ! rk_m1(69) = rk_het(5) ! HNO3 --> NO2
3034 ! rk_m1(70) = rk_het(6) ! HNO3 --> NO
3035 ! rk_m1(71) = rk_het(7) ! NO3 --> NO + O2
3036
3037 ! all rate constants but be >= 0
3038 do i = 1, nrxn_m1
3039 rk_m1(i) = max( rk_m1(i), 0.0 )
3040 end do
3041 do i = 1, nrxn_m2
3042 rk_m2(i) = max( rk_m2(i), 0.0 )
3043 end do
3044
3045 return
3046 end subroutine gasthermrk_m1
3047
3048
3049
3050 !***********************************************************************
3051 ! <27.> subr gasthermrk_m2
3052 !
3053 ! purpose: computes thermal reaction rate coefficients for
3054 ! mechanism 2
3055 !
3056 ! author : Rahul A. Zaveri
3057 ! date : December 1998
3058 !
3059 !-------------------------------------------------------------------------
3060
3061 subroutine gasthermrk_m2( tempbox, cair_mlc, &
3062 rk_photo, rk_param, rk_m2 )
3063
3064 use module_data_cbmz
3065 implicit none
3066
3067 ! subr arguments
3068 real tempbox, cair_mlc
3069 real rk_photo(nphoto), rk_param(nperox), rk_m2(nrxn_m2)
3070 ! local variables
3071 integer i
3072 real rk0, rki, rmm, rnn, te
3073 ! real arr, troe
3074
3075
3076 te = tempbox
3077
3078 rk_m2(1) = 8.1e-13
3079
3080 rk_m2(5) = rk_photo(jphoto_aone)
3081 rk_m2(6) = te**2*arr(5.3e-18, -230., te)
3082 rk_m2(7) = rk_photo(jphoto_mgly)
3083 rk_m2(8) = 1.7e-11
3084 rk_m2(9) = arr(1.4e-12, -1900., te)
3085 rk_m2(10) = arr(1.2e-14, -2630., te)
3086 !
3087 rk0 = 1.0e-28
3088 rnn = 0.8
3089 rki = 8.8e-12
3090 rmm = 0.0
3091 rk_m2(11) = troe(cair_mlc,te,rk0,rnn,rki,rmm)
3092 !
3093 rk_m2(12) = arr(4.2e-15, -1800., te)
3094 rk_m2(13) = arr(8.9e-16, -392., te)
3095 rk_m2(14) = arr(5.8e-12, 478., te)
3096 rk_m2(15) = arr(2.9e-11, 255., te)
3097 rk_m2(16) = arr(3.1e-13, -1010., te)
3098 rk_m2(17) = 2.5e-12
3099 rk_m2(18) = arr(2.1e-12, 322., te)
3100 rk_m2(19) = arr(1.7e-11, 116., te)
3101 rk_m2(20) = 8.1e-12
3102 rk_m2(21) = 4.1e-11
3103 rk_m2(22) = 2.2e-11
3104 rk_m2(23) = 1.4e-11
3105 rk_m2(24) = 3.0e-11
3106 rk_m2(25) = rk_photo(jphoto_open)
3107 rk_m2(26) = arr(5.4e-17, -500., te)
3108 rk_m2(27) = rk_photo(jphoto_rooh)
3109 rk_m2(28) = arr(3.8e-12, 200., te)
3110 rk_m2(29) = arr(1.6e-11, -540., te)
3111 rk_m2(30) = rk_photo(jphoto_onit)
3112
3113 rk_m2(33) = 4.0e-12
3114
3115 rk_m2(35) = 4.0e-12
3116 rk_m2(36) = 4.0e-12
3117 rk_m2(37) = 4.0e-12
3118 rk_m2(38) = 2.5e-12
3119
3120 rk_m2(40) = 1.2e-12
3121 rk_m2(41) = 4.0e-12
3122 rk_m2(42) = 2.5e-12
3123 rk_m2(43) = arr(1.7e-13, 1300., te)
3124
3125 rk_m2(45) = arr(1.2e-13, 1300., te)
3126 rk_m2(46) = arr(1.7e-13, 1300., te)
3127 rk_m2(47) = arr(1.7e-13, 1300., te)
3128 rk_m2(48) = rk_param(jro2)
3129
3130 rk_m2(50) = rk_param(jano2)
3131 rk_m2(51) = rk_param(jnap)
3132 rk_m2(52) = rk_param(jxo2)
3133 rk_m2(53) = 1.0e-11 ! XPAR + PAR -->
3134
3135 ! all rate constants but be >= 0
3136 do i = 1, nrxn_m2
3137 rk_m2(i) = max( rk_m2(i), 0.0 )
3138 end do
3139
3140 return
3141 end subroutine gasthermrk_m2
3142
3143
3144
3145 !***********************************************************************
3146 ! <28.> subr gasthermrk_m3
3147 !
3148 ! purpose: computes thermal reaction rate coefficients for
3149 ! mechanism 3
3150 !
3151 ! author : Rahul A. Zaveri
3152 ! date : December 1998
3153 !
3154 !-------------------------------------------------------------------------
3155
3156 subroutine gasthermrk_m3( tempbox, cair_mlc, &
3157 rk_photo, rk_param, rk_m3 )
3158
3159 use module_data_cbmz
3160 implicit none
3161
3162 ! subr arguments
3163 real tempbox, cair_mlc
3164 real rk_photo(nphoto), rk_param(nperox), rk_m3(nrxn_m3)
3165 ! local variables
3166 integer i
3167 real te
3168 ! real arr
3169
3170
3171 te = tempbox
3172 !
3173 rk_m3(1) = arr(2.6e-11, 409., te)
3174 rk_m3(2) = arr(1.2e-14, -2013., te)
3175 rk_m3(3) = arr(3.0e-12, -446., te)
3176 rk_m3(4) = rk_photo(jphoto_isoprd)
3177 rk_m3(5) = 3.3e-11
3178 rk_m3(6) = 7.0e-18
3179 rk_m3(7) = 1.0e-15
3180 rk_m3(8) = 4.0e-12
3181 rk_m3(9) = 4.0e-12
3182 rk_m3(10) = 4.0e-12
3183 rk_m3(11) = arr(1.7e-13, 1300., te)
3184 rk_m3(12) = arr(1.7e-13, 1300., te)
3185 rk_m3(13) = arr(1.7e-13, 1300., te)
3186 rk_m3(14) = rk_param(jisopp)
3187 rk_m3(15) = rk_param(jisopn)
3188 rk_m3(16) = rk_param(jisopo2)
3189
3190 ! all rate constants but be >= 0
3191 do i = 1, nrxn_m3
3192 rk_m3(i) = max( rk_m3(i), 0.0 )
3193 end do
3194
3195 return
3196 end subroutine gasthermrk_m3
3197
3198
3199
3200 !***********************************************************************
3201 ! <29.> subr gasthermrk_m4
3202 !
3203 ! purpose: computes thermal reaction rate coefficients for
3204 ! mechanism 4
3205 !
3206 ! author : Rahul A. Zaveri
3207 ! date : December 1998
3208 !
3209 !-------------------------------------------------------------------------
3210
3211 subroutine gasthermrk_m4( tempbox, cair_mlc, &
3212 rk_photo, rk_param, rk_m4 )
3213
3214 use module_data_cbmz
3215 implicit none
3216
3217 ! subr arguments
3218 real tempbox, cair_mlc
3219 real rk_photo(nphoto), rk_param(nperox), rk_m4(nrxn_m4)
3220 ! local variables
3221 integer i
3222 real B_abs, B_add, rk_tot, rk_tot_den, rk_tot_num, te
3223 ! real arr
3224
3225
3226 te = tempbox
3227 !
3228 rk_m4(1) = arr(9.6e-12, -234., te) ! ch3sch3 + oh --> ch3sch2
3229 rk_m4(2) = arr(1.4e-13, 500., te)
3230 rk_m4(3) = arr(1.3e-11, 409., te)
3231
3232 ! Hynes et al. (1986)
3233 rk_tot_num = te * exp(-234./te) + &
3234 8.46e-10 * exp(7230./te) + &
3235 2.68e-10 * exp(7810./te)
3236 rk_tot_den = 1.04e+11 * te + 88.1 * exp(7460./te)
3237 rk_tot = rk_tot_num/rk_tot_den
3238 B_abs = rk_m4(1)/rk_tot
3239 B_add = 1. - B_abs
3240
3241 rk_m4(4) = B_add*rk_tot ! ch3sch3 + oh --> ch3s(oh)ch3
3242 rk_m4(5) = 8.0e-12
3243 rk_m4(6) = 1.8e-13
3244 rk_m4(7) = 2.5e-13
3245 rk_m4(8) = 8.6e-14
3246 rk_m4(9) = 5.8e-11
3247 rk_m4(10) = 1.0e-14
3248 rk_m4(11) = 5.0e-12
3249 rk_m4(12) = 1.8e-13
3250 rk_m4(13) = 1.0e-15
3251 rk_m4(14) = 1.0e-13
3252 rk_m4(15) = 1.0e-15
3253 rk_m4(16) = 1.6e-11
3254 rk_m4(17) = 1.0e-13
3255 rk_m4(18) = arr(2.5e-13, -8686., te)
3256 rk_m4(19) = 1.0e-14
3257 rk_m4(20) = 5.0e-15
3258 rk_m4(21) = 2.5e-13
3259 rk_m4(22) = 2.5e-13
3260 rk_m4(23) = 5.0e-11
3261 rk_m4(24) = 2.6e-18
3262 rk_m4(25) = 3.3
3263 rk_m4(26) = 1.0e-11
3264 rk_m4(27) = 5.5e-12
3265 rk_m4(28) = arr(2.0e17, -12626., te)
3266 rk_m4(29) = 3.0e-15
3267 rk_m4(30) = 3.0e-15
3268 rk_m4(31) = 5.0e-11
3269 rk_m4(32) = 1.6e-15
3270
3271 ! all rate constants but be >= 0
3272 do i = 1, nrxn_m4
3273 rk_m4(i) = max( rk_m4(i), 0.0 )
3274 end do
3275
3276 return
3277 end subroutine gasthermrk_m4
3278
3279
3280
3281
3282 !***********************************************************************
3283 ! <26.> subr hetrateconstants
3284 !
3285 ! purpose: computes heterogeneous reaction rate coefficients
3286 !
3287 ! author : Rahul A. Zaveri
3288 ! date : May 2000
3289 !
3290 !-------------------------------------------------------------------------
3291
3292 subroutine hetrateconstants
3293 implicit none
3294
3295 return
3296 end subroutine hetrateconstants
3297
3298
3299
3300 !***********************************************************************
3301 ! <31.> func troe
3302 !
3303 ! purpose: calculates Troe reaction rate coefficient
3304 !
3305 ! author : Rahul A. Zaveri
3306 ! date : December 1998
3307 !-----------------------------------------------------------------------
3308
3309 real function troe( cairmlc, te, rk0, rnn, rki, rmm )
3310 implicit none
3311 ! func parameters
3312 real cairmlc, te, rk0, rnn, rki, rmm
3313 ! local variables
3314 real expo
3315
3316 rk0 = rk0*cairmlc*(te/300.)**(-rnn)
3317 rki = rki*(te/300.)**(-rmm)
3318 expo= 1./(1. + (ALOG10(rk0/rki))**2)
3319 troe = (rk0*rki/(rk0+rki))*.6**expo
3320 return
3321 end function troe
3322
3323
3324
3325 !***********************************************************************
3326 ! <32.> func arr
3327 !
3328 ! purpose: calculates arrhenius rate coefficient
3329 !
3330 ! author : Rahul A. Zaveri
3331 ! date : December 1998
3332 !-----------------------------------------------------------------------
3333
3334 real function arr( aa, bb, te )
3335 implicit none
3336 ! func parameters
3337 real aa, bb, te
3338
3339 arr = aa*exp(bb/te)
3340 return
3341 end function arr
3342
3343
3344
3345 !***********************************************************************
3346 ! <xx.> subr mapgas_tofrom_host
3347 !
3348 ! purpose: maps gas species between cboxold/new and host arrays
3349 !
3350 ! author : R. C. Easter
3351 ! date : November, 2003
3352 !
3353 ! ----------------------------------------------------------------------
3354
3355 subroutine mapgas_tofrom_host( imap, &
3356 i_boxtest_units_convert, &
3357 it,jt,kt, ims,ime, jms,jme, kms,kme, &
3358 num_moist, num_chem, moist, chem, &
3359 t_phy, p_phy, rho_phy, &
3360 cbox, tempbox, pressbox, airdenbox, &
3361 cair_mlc, &
3362 h2o, ch4, oxygen, nitrogen, hydrogen )
3363
3364 use module_configure, only: &
3365 p_qv, &
3366 p_so2, p_sulf, p_no2, p_no, p_o3, &
3367 p_hno3, p_h2o2, p_ald, p_hcho, p_op1, &
3368 p_op2, p_paa, p_ora1, p_ora2, p_nh3, &
3369 p_n2o5, p_no3, p_pan, p_hc3, p_hc5, &
3370 p_hc8, p_eth, p_co, p_ol2, p_olt, &
3371 p_oli, p_tol, p_xyl, p_aco3, p_tpan, &
3372 p_hono, p_hno4, p_ket, p_gly, p_mgly, &
3373 p_dcb, p_onit, p_csl, p_iso, p_ho, &
3374 p_ho2, &
3375 p_hcl, p_ch3o2, p_ethp, p_ch3oh, p_c2h5oh, &
3376 p_par, p_to2, p_cro, p_open, p_op3, &
3377 p_c2o3, p_ro2, p_ano2, p_nap, p_xo2, &
3378 p_xpar, p_isoprd, p_isopp, p_isopn, p_isopo2, &
3379 p_dms, p_msa, p_dmso, p_dmso2, p_ch3so2h, &
3380 p_ch3sch2oo, p_ch3so2, p_ch3so3, p_ch3so2oo, p_ch3so2ch2oo, &
3381 p_mtf
3382 use module_data_cbmz
3383 implicit none
3384
3385 ! subr arguments
3386 INTEGER, INTENT(IN) :: imap, it,jt,kt, ims,ime, jms,jme, kms,kme, &
3387 num_moist, num_chem, i_boxtest_units_convert
3388 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_moist ), &
3389 INTENT(IN) :: moist
3390 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
3391 INTENT(INOUT) :: chem
3392 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , &
3393 INTENT(IN) :: t_phy, & ! temperature
3394 p_phy, & ! air pressure (Pa)
3395 rho_phy ! air density (kg/m3)
3396 REAL, INTENT(INOUT) :: cbox(ngas_z)
3397 REAL, INTENT(INOUT) :: tempbox, pressbox, airdenbox
3398 REAL, INTENT(INOUT) :: cair_mlc
3399 REAL, INTENT(INOUT) :: h2o, ch4, oxygen, nitrogen, hydrogen
3400
3401 ! local variables
3402 integer l
3403 real factoraa
3404 real, parameter :: eps=0.622
3405
3406
3407 tempbox = t_phy(it,kt,jt)
3408 ! p_phy = (Pa); pressbox = (dynes/cm2)
3409 pressbox = p_phy(it,kt,jt)*10.0
3410 ! rho_phy = (kg_air/m3); airdenbox = (mole_air/cm3)
3411 airdenbox = rho_phy(it,kt,jt)/28.966e3
3412 if (i_boxtest_units_convert .eq. 10) then
3413 airdenbox = rho_phy(it,kt,jt)
3414 end if
3415
3416 if (imap .gt. 0) goto 2000
3417
3418 !
3419 ! imap==0 -- initial species mapping from host array to cboxold
3420 ! chem --> czz --> cbox
3421 !
3422 ! note: do not map nh3, hcl
3423 !
3424 cbox(:) = 0.0
3425
3426 ! cair_mlc = (molecules_air/cm3)
3427 cair_mlc = airdenbox*avognumkpp
3428
3429 ! moist = (kg_h2o/kg_air); czz = (mole_h2o/cm3); h2o = (molecules_h2o/cm3)
3430 ! czz(ih2o_z) = (moist(it,kt,jt,p_qv)/eps)*airdenbox
3431 ! if (i_boxtest_units_convert .eq. 10) then
3432 ! czz(ih2o_z) = moist(it,kt,jt,p_qv)*airdenbox
3433 ! end if
3434 ! h2o = czz(ih2o_z)*avognumkpp
3435 h2o = (moist(it,kt,jt,p_qv)/eps)*airdenbox
3436 if (i_boxtest_units_convert .eq. 10) then
3437 h2o = moist(it,kt,jt,p_qv)*airdenbox
3438 end if
3439 h2o = h2o*avognumkpp
3440
3441 ! czz(ich4_z) = 1.7e-6*airdenbox ! ch4 conc. in mol/cc
3442 ! ch4 = czz(ich4_z)*avognumkpp ! ch4 conc. in molec/cc
3443 ch4 = 1.7e-6*airdenbox*avognumkpp ! ch4 conc. in molec/cc
3444
3445 oxygen = 0.21*cair_mlc ! o2 conc. in molec/cc
3446 nitrogen = 0.79*cair_mlc ! n2 conc. in molec/cc
3447 hydrogen = 0.58e-6*cair_mlc ! h2 conc. in molec/cc
3448
3449 ! chem units = (ppm); czz units = (mole/cm3); cbox units = (molecules/cm3)
3450 factoraa = airdenbox*1.0e-6
3451 if (i_boxtest_units_convert .eq. 10) factoraa = airdenbox
3452 factoraa = factoraa*avognumkpp
3453
3454 cbox(iso2_z) = chem(it,kt,jt,p_so2)*factoraa
3455 cbox(ih2so4_z) = chem(it,kt,jt,p_sulf)*factoraa
3456 cbox(ino2_z) = chem(it,kt,jt,p_no2)*factoraa
3457 cbox(ino_z) = chem(it,kt,jt,p_no)*factoraa
3458 cbox(io3_z) = chem(it,kt,jt,p_o3)*factoraa
3459 cbox(ihno3_z) = chem(it,kt,jt,p_hno3)*factoraa
3460 cbox(ih2o2_z) = chem(it,kt,jt,p_h2o2)*factoraa
3461 cbox(iald2_z) = chem(it,kt,jt,p_ald)*factoraa
3462 cbox(ihcho_z) = chem(it,kt,jt,p_hcho)*factoraa
3463 cbox(ich3ooh_z) = chem(it,kt,jt,p_op1)*factoraa
3464 cbox(iethooh_z) = chem(it,kt,jt,p_op2)*factoraa
3465 cbox(ihcooh_z) = chem(it,kt,jt,p_ora1)*factoraa
3466 cbox(ircooh_z) = chem(it,kt,jt,p_ora2)*factoraa
3467 cbox(inh3_z) = chem(it,kt,jt,p_nh3)*factoraa
3468 cbox(in2o5_z) = chem(it,kt,jt,p_n2o5)*factoraa
3469 cbox(ino3_z) = chem(it,kt,jt,p_no3)*factoraa
3470 cbox(ipan_z) = chem(it,kt,jt,p_pan)*factoraa
3471 cbox(ic2h6_z) = chem(it,kt,jt,p_eth)*factoraa
3472 cbox(ico_z) = chem(it,kt,jt,p_co)*factoraa
3473 cbox(ieth_z) = chem(it,kt,jt,p_ol2)*factoraa
3474 cbox(iolet_z) = chem(it,kt,jt,p_olt)*factoraa
3475 cbox(iolei_z) = chem(it,kt,jt,p_oli)*factoraa
3476 cbox(itol_z) = chem(it,kt,jt,p_tol)*factoraa
3477 cbox(ixyl_z) = chem(it,kt,jt,p_xyl)*factoraa
3478 cbox(ihono_z) = chem(it,kt,jt,p_hono)*factoraa
3479 cbox(ihno4_z) = chem(it,kt,jt,p_hno4)*factoraa
3480 cbox(iaone_z) = chem(it,kt,jt,p_ket)*factoraa
3481 cbox(imgly_z) = chem(it,kt,jt,p_mgly)*factoraa
3482 cbox(ionit_z) = chem(it,kt,jt,p_onit)*factoraa
3483 cbox(icres_z) = chem(it,kt,jt,p_csl)*factoraa
3484 cbox(iisop_z) = chem(it,kt,jt,p_iso)*factoraa
3485 cbox(ioh_z) = chem(it,kt,jt,p_ho)*factoraa
3486 cbox(iho2_z) = chem(it,kt,jt,p_ho2)*factoraa
3487
3488 cbox(ihcl_z) = chem(it,kt,jt,p_hcl)*factoraa
3489 cbox(ich3o2_z) = chem(it,kt,jt,p_ch3o2)*factoraa
3490 cbox(iethp_z) = chem(it,kt,jt,p_ethp)*factoraa
3491 cbox(ich3oh_z) = chem(it,kt,jt,p_ch3oh)*factoraa
3492 cbox(ic2h5oh_z) = chem(it,kt,jt,p_c2h5oh)*factoraa
3493 cbox(ipar_z) = chem(it,kt,jt,p_par)*factoraa
3494 cbox(ito2_z) = chem(it,kt,jt,p_to2)*factoraa
3495 cbox(icro_z) = chem(it,kt,jt,p_cro)*factoraa
3496 cbox(iopen_z) = chem(it,kt,jt,p_open)*factoraa
3497 cbox(irooh_z) = chem(it,kt,jt,p_op3)*factoraa
3498 cbox(ic2o3_z) = chem(it,kt,jt,p_c2o3)*factoraa
3499 cbox(iro2_z) = chem(it,kt,jt,p_ro2)*factoraa
3500 cbox(iano2_z) = chem(it,kt,jt,p_ano2)*factoraa
3501 cbox(inap_z) = chem(it,kt,jt,p_nap)*factoraa
3502 cbox(ixo2_z) = chem(it,kt,jt,p_xo2)*factoraa
3503 cbox(ixpar_z) = chem(it,kt,jt,p_xpar)*factoraa
3504 cbox(iisoprd_z) = chem(it,kt,jt,p_isoprd)*factoraa
3505 cbox(iisopp_z) = chem(it,kt,jt,p_isopp)*factoraa
3506 cbox(iisopn_z) = chem(it,kt,jt,p_isopn)*factoraa
3507 cbox(iisopo2_z) = chem(it,kt,jt,p_isopo2)*factoraa
3508 cbox(idms_z) = chem(it,kt,jt,p_dms)*factoraa
3509 cbox(imsa_z) = chem(it,kt,jt,p_msa)*factoraa
3510 cbox(idmso_z) = chem(it,kt,jt,p_dmso)*factoraa
3511 cbox(idmso2_z) = chem(it,kt,jt,p_dmso2)*factoraa
3512 cbox(ich3so2h_z) = chem(it,kt,jt,p_ch3so2h)*factoraa
3513 cbox(ich3sch2oo_z) = chem(it,kt,jt,p_ch3sch2oo)*factoraa
3514 cbox(ich3so2_z) = chem(it,kt,jt,p_ch3so2)*factoraa
3515 cbox(ich3so3_z) = chem(it,kt,jt,p_ch3so3)*factoraa
3516 cbox(ich3so2oo_z) = chem(it,kt,jt,p_ch3so2oo)*factoraa
3517 cbox(ich3so2ch2oo_z) = chem(it,kt,jt,p_ch3so2ch2oo)*factoraa
3518 cbox(imtf_z) = chem(it,kt,jt,p_mtf)*factoraa
3519
3520 cbox(ih2o_z) = h2o
3521 cbox(ich4_z) = ch4
3522 cbox(io2_z) = oxygen
3523 cbox(in2_z) = nitrogen
3524 cbox(ih2_z) = hydrogen
3525
3526 return
3527
3528 !
3529 ! imap==1 -- final species mapping from cbox back to host array
3530 ! cbox --> czz --> chem
3531 !
3532 ! note1: do not map nh3, hcl, ch4
3533 !
3534 2000 continue
3535 ! chem = (ppm); czz = (mole/cm3); cbox = (molecules/cm3)
3536 factoraa = airdenbox*1.0e-6
3537 if (i_boxtest_units_convert .eq. 10) factoraa = airdenbox
3538 factoraa = factoraa*avognumkpp
3539
3540 chem(it,kt,jt,p_so2) = cbox(iso2_z)/factoraa
3541 chem(it,kt,jt,p_sulf) = cbox(ih2so4_z)/factoraa
3542 chem(it,kt,jt,p_no2) = cbox(ino2_z)/factoraa
3543 chem(it,kt,jt,p_no) = cbox(ino_z)/factoraa
3544 chem(it,kt,jt,p_o3) = cbox(io3_z)/factoraa
3545 chem(it,kt,jt,p_hno3) = cbox(ihno3_z)/factoraa
3546 chem(it,kt,jt,p_h2o2) = cbox(ih2o2_z)/factoraa
3547 chem(it,kt,jt,p_ald) = cbox(iald2_z)/factoraa
3548 chem(it,kt,jt,p_hcho) = cbox(ihcho_z)/factoraa
3549 chem(it,kt,jt,p_op1) = cbox(ich3ooh_z)/factoraa
3550 chem(it,kt,jt,p_op2) = cbox(iethooh_z)/factoraa
3551 chem(it,kt,jt,p_ora1) = cbox(ihcooh_z)/factoraa
3552 chem(it,kt,jt,p_ora2) = cbox(ircooh_z)/factoraa
3553 chem(it,kt,jt,p_nh3) = cbox(inh3_z)/factoraa
3554 chem(it,kt,jt,p_n2o5) = cbox(in2o5_z)/factoraa
3555 chem(it,kt,jt,p_no3) = cbox(ino3_z)/factoraa
3556 chem(it,kt,jt,p_pan) = cbox(ipan_z)/factoraa
3557 chem(it,kt,jt,p_eth) = cbox(ic2h6_z)/factoraa
3558 chem(it,kt,jt,p_co) = cbox(ico_z)/factoraa
3559 chem(it,kt,jt,p_ol2) = cbox(ieth_z)/factoraa
3560 chem(it,kt,jt,p_olt) = cbox(iolet_z)/factoraa
3561 chem(it,kt,jt,p_oli) = cbox(iolei_z)/factoraa
3562 chem(it,kt,jt,p_tol) = cbox(itol_z)/factoraa
3563 chem(it,kt,jt,p_xyl) = cbox(ixyl_z)/factoraa
3564 chem(it,kt,jt,p_hono) = cbox(ihono_z)/factoraa
3565 chem(it,kt,jt,p_hno4) = cbox(ihno4_z)/factoraa
3566 chem(it,kt,jt,p_ket) = cbox(iaone_z)/factoraa
3567 chem(it,kt,jt,p_mgly) = cbox(imgly_z)/factoraa
3568 chem(it,kt,jt,p_onit) = cbox(ionit_z)/factoraa
3569 chem(it,kt,jt,p_csl) = cbox(icres_z)/factoraa
3570 chem(it,kt,jt,p_iso) = cbox(iisop_z)/factoraa
3571 chem(it,kt,jt,p_ho) = cbox(ioh_z)/factoraa
3572 chem(it,kt,jt,p_ho2) = cbox(iho2_z)/factoraa
3573
3574 chem(it,kt,jt,p_hcl) = cbox(ihcl_z)/factoraa
3575 chem(it,kt,jt,p_ch3o2) = cbox(ich3o2_z)/factoraa
3576 chem(it,kt,jt,p_ethp) = cbox(iethp_z)/factoraa
3577 chem(it,kt,jt,p_ch3oh) = cbox(ich3oh_z)/factoraa
3578 chem(it,kt,jt,p_c2h5oh) = cbox(ic2h5oh_z)/factoraa
3579 chem(it,kt,jt,p_par) = cbox(ipar_z)/factoraa
3580 chem(it,kt,jt,p_to2) = cbox(ito2_z)/factoraa
3581 chem(it,kt,jt,p_cro) = cbox(icro_z)/factoraa
3582 chem(it,kt,jt,p_open) = cbox(iopen_z)/factoraa
3583 chem(it,kt,jt,p_op3) = cbox(irooh_z)/factoraa
3584 chem(it,kt,jt,p_c2o3) = cbox(ic2o3_z)/factoraa
3585 chem(it,kt,jt,p_ro2) = cbox(iro2_z)/factoraa
3586 chem(it,kt,jt,p_ano2) = cbox(iano2_z)/factoraa
3587 chem(it,kt,jt,p_nap) = cbox(inap_z)/factoraa
3588 chem(it,kt,jt,p_xo2) = cbox(ixo2_z)/factoraa
3589 chem(it,kt,jt,p_xpar) = cbox(ixpar_z)/factoraa
3590 chem(it,kt,jt,p_isoprd) = cbox(iisoprd_z)/factoraa
3591 chem(it,kt,jt,p_isopp) = cbox(iisopp_z)/factoraa
3592 chem(it,kt,jt,p_isopn) = cbox(iisopn_z)/factoraa
3593 chem(it,kt,jt,p_isopo2) = cbox(iisopo2_z)/factoraa
3594 chem(it,kt,jt,p_dms) = cbox(idms_z)/factoraa
3595 chem(it,kt,jt,p_msa) = cbox(imsa_z)/factoraa
3596 chem(it,kt,jt,p_dmso) = cbox(idmso_z)/factoraa
3597 chem(it,kt,jt,p_dmso2) = cbox(idmso2_z)/factoraa
3598 chem(it,kt,jt,p_ch3so2h) = cbox(ich3so2h_z)/factoraa
3599 chem(it,kt,jt,p_ch3sch2oo) = cbox(ich3sch2oo_z)/factoraa
3600 chem(it,kt,jt,p_ch3so2) = cbox(ich3so2_z)/factoraa
3601 chem(it,kt,jt,p_ch3so3) = cbox(ich3so3_z)/factoraa
3602 chem(it,kt,jt,p_ch3so2oo) = cbox(ich3so2oo_z)/factoraa
3603 chem(it,kt,jt,p_ch3so2ch2oo) = cbox(ich3so2ch2oo_z)/factoraa
3604 chem(it,kt,jt,p_mtf) = cbox(imtf_z)/factoraa
3605
3606 return
3607 end subroutine mapgas_tofrom_host
3608
3609
3610
3611 !***********************************************************************
3612 ! <xx.> subr set_gaschem_allowed_regimes
3613 !
3614 ! purpose: determines which gas-phase chemistry regimes are allowed based
3615 ! on which species are active in the simulation
3616 !
3617 ! author :
3618 ! date :
3619 !
3620 ! ----------------------------------------------------------------------
3621
3622 subroutine set_gaschem_allowed_regimes( lunerr, &
3623 igaschem_allowed_m1, igaschem_allowed_m2, &
3624 igaschem_allowed_m3, igaschem_allowed_m4 )
3625 !
3626 ! determines which gas-phase chemistry regimes are allowed based
3627 ! on which species are active in the simulation
3628 !
3629 use module_configure, only: &
3630 p_qv, &
3631 p_so2, p_sulf, p_no2, p_no, p_o3, &
3632 p_hno3, p_h2o2, p_ald, p_hcho, p_op1, &
3633 p_op2, p_paa, p_ora1, p_ora2, p_nh3, &
3634 p_n2o5, p_no3, p_pan, p_hc3, p_hc5, &
3635 p_hc8, p_eth, p_co, p_ol2, p_olt, &
3636 p_oli, p_tol, p_xyl, p_aco3, p_tpan, &
3637 p_hono, p_hno4, p_ket, p_gly, p_mgly, &
3638 p_dcb, p_onit, p_csl, p_iso, p_ho, &
3639 p_ho2, &
3640 p_hcl, p_ch3o2, p_ethp, p_ch3oh, p_c2h5oh, &
3641 p_par, p_to2, p_cro, p_open, p_op3, &
3642 p_c2o3, p_ro2, p_ano2, p_nap, p_xo2, &
3643 p_xpar, p_isoprd, p_isopp, p_isopn, p_isopo2, &
3644 p_dms, p_msa, p_dmso, p_dmso2, p_ch3so2h, &
3645 p_ch3sch2oo, p_ch3so2, p_ch3so3, p_ch3so2oo, p_ch3so2ch2oo, &
3646 p_mtf
3647 use module_state_description, only: param_first_scalar
3648 use module_data_cbmz
3649 implicit none
3650
3651 ! subr arguments
3652 integer lunerr
3653 integer igaschem_allowed_m1, igaschem_allowed_m2, &
3654 igaschem_allowed_m3, igaschem_allowed_m4
3655
3656 ! local variables
3657 integer nactive, ndum, p1st
3658 character*80 msg
3659
3660
3661 ! index for first "active" scalar (= 2)
3662 p1st = param_first_scalar
3663
3664 ! determine if regime 1 is allowed
3665 ! (note: p_xxx>1 if xxx is active, p_xxx=1 if inactive)
3666 if (p_qv .lt. p1st) then
3667 msg = '*** subr set_gaschem_allowed_regimes'
3668 call peg_message( lunerr, msg )
3669 msg = '*** water vapor IS NOT ACTIVE'
3670 call peg_message( lunerr, msg )
3671 call peg_error_fatal( lunerr, msg )
3672 end if
3673
3674 ! determine if regime 1 is allowed
3675 ! (note: p_xxx>1 if xxx is active, p_xxx=1 if inactive)
3676 nactive = 0
3677 ndum = 27
3678 if (p_no .ge. p1st) nactive = nactive + 1
3679 if (p_no2 .ge. p1st) nactive = nactive + 1
3680 if (p_no3 .ge. p1st) nactive = nactive + 1
3681 if (p_n2o5 .ge. p1st) nactive = nactive + 1
3682 if (p_hono .ge. p1st) nactive = nactive + 1
3683 if (p_hno3 .ge. p1st) nactive = nactive + 1
3684 if (p_hno4 .ge. p1st) nactive = nactive + 1
3685 if (p_o3 .ge. p1st) nactive = nactive + 1
3686 ! if (p_o1d .ge. p1st) nactive = nactive + 1
3687 ! if (p_o3p .ge. p1st) nactive = nactive + 1
3688 if (p_ho .ge. p1st) nactive = nactive + 1
3689 if (p_ho2 .ge. p1st) nactive = nactive + 1
3690 if (p_h2o2 .ge. p1st) nactive = nactive + 1
3691 if (p_co .ge. p1st) nactive = nactive + 1
3692 if (p_so2 .ge. p1st) nactive = nactive + 1
3693 if (p_sulf .ge. p1st) nactive = nactive + 1
3694 ! if (p_nh3 .ge. p1st) nactive = nactive + 1
3695 ! if (p_hcl .ge. p1st) nactive = nactive + 1
3696 if (p_eth .ge. p1st) nactive = nactive + 1
3697 if (p_ch3o2 .ge. p1st) nactive = nactive + 1
3698 if (p_ethp .ge. p1st) nactive = nactive + 1
3699 if (p_hcho .ge. p1st) nactive = nactive + 1
3700 if (p_ch3oh .ge. p1st) nactive = nactive + 1
3701 if (p_c2h5oh .ge. p1st) nactive = nactive + 1
3702 if (p_op1 .ge. p1st) nactive = nactive + 1
3703 if (p_op2 .ge. p1st) nactive = nactive + 1
3704 if (p_ald .ge. p1st) nactive = nactive + 1
3705 if (p_ora1 .ge. p1st) nactive = nactive + 1
3706 if (p_pan .ge. p1st) nactive = nactive + 1
3707 if (p_ora2 .ge. p1st) nactive = nactive + 1
3708 if (p_c2o3 .ge. p1st) nactive = nactive + 1
3709
3710 if (nactive .le. 0) then
3711 igaschem_allowed_m1 = 0
3712 else if (nactive .eq. ndum) then
3713 igaschem_allowed_m1 = 1
3714 else
3715 msg = '*** subr set_gaschem_allowed_regimes'
3716 call peg_message( lunerr, msg )
3717 write(msg,90200) 1, nactive, ndum
3718 call peg_message( lunerr, msg )
3719 call peg_error_fatal( lunerr, msg )
3720 end if
3721 90200 format( ' error for regime ', i1, ', nactive, nexpected = ', 2i5 )
3722
3723 ! determine if regime 2 is allowed
3724 nactive = 0
3725 ndum = 19
3726 if (p_par .ge. p1st) nactive = nactive + 1
3727 if (p_ket .ge. p1st) nactive = nactive + 1
3728 if (p_mgly .ge. p1st) nactive = nactive + 1
3729 if (p_ol2 .ge. p1st) nactive = nactive + 1
3730 if (p_olt .ge. p1st) nactive = nactive + 1
3731 if (p_oli .ge. p1st) nactive = nactive + 1
3732 if (p_tol .ge. p1st) nactive = nactive + 1
3733 if (p_xyl .ge. p1st) nactive = nactive + 1
3734 if (p_csl .ge. p1st) nactive = nactive + 1
3735 if (p_to2 .ge. p1st) nactive = nactive + 1
3736 if (p_cro .ge. p1st) nactive = nactive + 1
3737 if (p_open .ge. p1st) nactive = nactive + 1
3738 if (p_onit .ge. p1st) nactive = nactive + 1
3739 if (p_op3 .ge. p1st) nactive = nactive + 1
3740 if (p_ro2 .ge. p1st) nactive = nactive + 1
3741 if (p_ano2 .ge. p1st) nactive = nactive + 1
3742 if (p_nap .ge. p1st) nactive = nactive + 1
3743 if (p_xo2 .ge. p1st) nactive = nactive + 1
3744 if (p_xpar .ge. p1st) nactive = nactive + 1
3745 if (nactive .le. 0) then
3746 igaschem_allowed_m2 = 0
3747 else if (nactive .eq. ndum) then
3748 igaschem_allowed_m2 = 2
3749 else
3750 msg = '*** subr set_gaschem_allowed_regimes'
3751 call peg_message( lunerr, msg )
3752 write(msg,90200) 2, nactive, ndum
3753 call peg_message( lunerr, msg )
3754 call peg_error_fatal( lunerr, msg )
3755 end if
3756
3757 ! determine if regime 3 is allowed
3758 nactive = 0
3759 ndum = 5
3760 if (p_iso .ge. p1st) nactive = nactive + 1
3761 if (p_isoprd .ge. p1st) nactive = nactive + 1
3762 if (p_isopp .ge. p1st) nactive = nactive + 1
3763 if (p_isopn .ge. p1st) nactive = nactive + 1
3764 if (p_isopo2 .ge. p1st) nactive = nactive + 1
3765 if (nactive .le. 0) then
3766 igaschem_allowed_m3 = 0
3767 else if (nactive .eq. ndum) then
3768 igaschem_allowed_m3 = 3
3769 else
3770 msg = '*** subr set_gaschem_allowed_regimes'
3771 call peg_message( lunerr, msg )
3772 write(msg,90200) 3, nactive, ndum
3773 call peg_message( lunerr, msg )
3774 call peg_error_fatal( lunerr, msg )
3775 end if
3776
3777 ! determine if regime 4 is allowed
3778 nactive = 0
3779 ndum = 11
3780 if (p_dms .ge. p1st) nactive = nactive + 1
3781 if (p_msa .ge. p1st) nactive = nactive + 1
3782 if (p_dmso .ge. p1st) nactive = nactive + 1
3783 if (p_dmso2 .ge. p1st) nactive = nactive + 1
3784 if (p_ch3so2h .ge. p1st) nactive = nactive + 1
3785 if (p_ch3sch2oo .ge. p1st) nactive = nactive + 1
3786 if (p_ch3so2 .ge. p1st) nactive = nactive + 1
3787 if (p_ch3so3 .ge. p1st) nactive = nactive + 1
3788 if (p_ch3so2oo .ge. p1st) nactive = nactive + 1
3789 if (p_ch3so2ch2oo .ge. p1st) nactive = nactive + 1
3790 if (p_mtf .ge. p1st) nactive = nactive + 1
3791 if (nactive .le. 0) then
3792 igaschem_allowed_m4 = 0
3793 else if (nactive .eq. ndum) then
3794 igaschem_allowed_m4 = 4
3795 else
3796 msg = '*** subr set_gaschem_allowed_regimes'
3797 call peg_message( lunerr, msg )
3798 write(msg,90200) 4, nactive, ndum
3799 call peg_message( lunerr, msg )
3800 call peg_error_fatal( lunerr, msg )
3801 end if
3802
3803 ! regime 1 must always be allowed
3804 if (igaschem_allowed_m1 .le. 0) then
3805 msg = '*** subr set_gaschem_allowed_regimes'
3806 call peg_message( lunerr, msg )
3807 write(msg,90300) 1
3808 call peg_message( lunerr, msg )
3809 call peg_error_fatal( lunerr, msg )
3810 end if
3811 90300 format( ' regime ', i1, ' must always be allowed' )
3812
3813 ! if regime 2 is allowed, then regime 1 must be allowed
3814 if (igaschem_allowed_m2 .gt. 0) then
3815 if (igaschem_allowed_m1 .le. 0) then
3816 msg = '*** subr set_gaschem_allowed_regimes'
3817 call peg_message( lunerr, msg )
3818 write(msg,90400) 2, 1
3819 call peg_message( lunerr, msg )
3820 call peg_error_fatal( lunerr, msg )
3821 end if
3822 end if
3823 90400 format( ' regime ', i1, ' allowed BUT regime ', i1, ' unallowed' )
3824
3825 ! if regime 3 is allowed, then regimes 1&2 must be allowed
3826 if (igaschem_allowed_m3 .gt. 0) then
3827 if (igaschem_allowed_m1 .le. 0) then
3828 msg = '*** subr set_gaschem_allowed_regimes'
3829 call peg_message( lunerr, msg )
3830 write(msg,90400) 3, 1
3831 call peg_message( lunerr, msg )
3832 call peg_error_fatal( lunerr, msg )
3833 else if (igaschem_allowed_m2 .le. 0) then
3834 msg = '*** subr set_gaschem_allowed_regimes'
3835 call peg_message( lunerr, msg )
3836 write(msg,90400) 3, 2
3837 call peg_message( lunerr, msg )
3838 call peg_error_fatal( lunerr, msg )
3839 end if
3840 end if
3841
3842 ! if regime 4 is allowed, then regime 1 must be allowed
3843 if (igaschem_allowed_m4 .gt. 0) then
3844 if (igaschem_allowed_m1 .le. 0) then
3845 msg = '*** subr set_gaschem_allowed_regimes'
3846 call peg_message( lunerr, msg )
3847 write(msg,90400) 4, 1
3848 call peg_message( lunerr, msg )
3849 call peg_error_fatal( lunerr, msg )
3850 end if
3851 end if
3852
3853 return
3854 end subroutine set_gaschem_allowed_regimes
3855
3856
3857
3858 !***********************************************************************
3859 ! <xx.> subr gasphotoconstants
3860 !
3861 ! purpose: copy photolytic rate constants from host arrays to local array
3862 !
3863 !-----------------------------------------------------------------------
3864 subroutine gasphotoconstants( rk_photo, &
3865 i_boxtest_units_convert, &
3866 it,jt,kt, ims,ime, jms,jme, kms,kme, &
3867 ph_o31d, ph_o33p, ph_no2, ph_no3o2, ph_no3o, ph_hno2, &
3868 ph_hno3, ph_hno4, ph_h2o2, ph_ch2or, ph_ch2om, &
3869 ph_ch3o2h, ph_n2o5 )
3870 !
3871 ! copies photolytic rate constants from host arrays to local arrays
3872 ! note1: currently 8 rate constants are scaled to other rate constants
3873 ! as is done in zz06gasphotolysis.f
3874 ! note2: currently the n2o5 rate is set to zero
3875 !
3876 use module_data_cbmz
3877 implicit none
3878
3879 ! subr arguments
3880 integer it,jt,kt, ims,ime, jms,jme, kms,kme
3881 integer i_boxtest_units_convert
3882 real rk_photo(nphoto)
3883 real, dimension( ims:ime, kms:kme, jms:jme ) :: &
3884 ph_o31d, ph_o33p, ph_no2, ph_no3o2, ph_no3o, ph_hno2, &
3885 ph_hno3, ph_hno4, ph_h2o2, ph_ch2or, ph_ch2om, &
3886 ph_ch3o2h, ph_n2o5
3887
3888 ! local variables
3889 real ft
3890
3891
3892 rk_photo(:) = 0.0
3893
3894 ! these from wrf/madronnich rate constants
3895 rk_photo(jphoto_no2) = ph_no2(it,kt,jt)
3896 rk_photo(jphoto_no3) = ph_no3o(it,kt,jt) &
3897 + ph_no3o2(it,kt,jt)
3898 rk_photo(jphoto_o3a) = ph_o33p(it,kt,jt)
3899 rk_photo(jphoto_o3b) = ph_o31d(it,kt,jt)
3900 rk_photo(jphoto_hono) = ph_hno2(it,kt,jt)
3901 rk_photo(jphoto_hno3) = ph_hno3(it,kt,jt)
3902 rk_photo(jphoto_hno4) = ph_hno4(it,kt,jt)
3903 rk_photo(jphoto_h2o2) = ph_h2o2(it,kt,jt)
3904 rk_photo(jphoto_ch3ooh) = ph_ch3o2h(it,kt,jt)
3905 rk_photo(jphoto_hchoa) = ph_ch2or(it,kt,jt)
3906 rk_photo(jphoto_hchob) = ph_ch2om(it,kt,jt)
3907 rk_photo(jphoto_n2o5) = ph_n2o5(it,kt,jt)
3908
3909 ! these scaled to other rate constants
3910 rk_photo(jphoto_ethooh) = 0.7 *rk_photo(jphoto_h2o2)
3911 rk_photo(jphoto_ald2) = 4.6e-4*rk_photo(jphoto_no2)
3912 rk_photo(jphoto_aone) = 7.8e-5*rk_photo(jphoto_no2)
3913 rk_photo(jphoto_mgly) = 9.64 *rk_photo(jphoto_hchoa)
3914 rk_photo(jphoto_open) = 9.04 *rk_photo(jphoto_hchoa)
3915 rk_photo(jphoto_rooh) = 0.7 *rk_photo(jphoto_h2o2)
3916 rk_photo(jphoto_onit) = 1.0e-4*rk_photo(jphoto_no2)
3917 rk_photo(jphoto_isoprd) = .025 *rk_photo(jphoto_hchob)
3918
3919 ! convert from (1/min) to (1/s)
3920 ! (except when i_boxtest_units_convert = 10 or 20)
3921 ft = 60.0
3922 if (i_boxtest_units_convert .eq. 10) ft = 1.0
3923 if (i_boxtest_units_convert .eq. 20) ft = 1.0
3924 if (ft .ne. 1.0) then
3925 rk_photo(:) = rk_photo(:)/ft
3926 end if
3927
3928
3929 return
3930 end subroutine gasphotoconstants
3931
3932
3933
3934 !-----------------------------------------------------------------------
3935 end module module_cbmz