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