1 | !***************************************************************** |
---|
2 | ! |
---|
3 | ! Photochemical routine |
---|
4 | ! |
---|
5 | ! Author: Franck Lefevre |
---|
6 | ! ------ |
---|
7 | ! |
---|
8 | ! Version: 10/11/2015 |
---|
9 | ! |
---|
10 | !***************************************************************** |
---|
11 | |
---|
12 | subroutine photochemistry_asis(nlayer, nq, ngrid, & |
---|
13 | ig, lswitch, zycol, sza, fractcol, ptimestep, press, & |
---|
14 | temp, dens, zmmean, dist_sol, surfdust1d, & |
---|
15 | surfice1d, jo3, tau, iter) |
---|
16 | |
---|
17 | use callkeys_mod |
---|
18 | implicit none |
---|
19 | |
---|
20 | #include "chimiedata.h" |
---|
21 | |
---|
22 | !=================================================================== |
---|
23 | ! inputs: |
---|
24 | !=================================================================== |
---|
25 | |
---|
26 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
27 | integer, intent(in) :: nq ! number of tracers |
---|
28 | integer,intent(in) :: ngrid ! number of atmospheric columns |
---|
29 | integer :: ig ! grid point index |
---|
30 | |
---|
31 | real :: sza ! solar zenith angle (deg) |
---|
32 | real :: fractcol ! day fraction |
---|
33 | real :: ptimestep ! physics timestep (s) |
---|
34 | real :: press(nlayer) ! pressure (hpa) |
---|
35 | real :: temp(nlayer) ! temperature (k) |
---|
36 | real :: dens(nlayer) ! density (cm-3) |
---|
37 | real :: zmmean(nlayer) ! mean molar mass (g/mole) |
---|
38 | real :: dist_sol ! sun distance (au) |
---|
39 | real :: surfdust1d(nlayer) ! dust surface area (cm2/cm3) |
---|
40 | real :: surfice1d(nlayer) ! ice surface area (cm2/cm3) |
---|
41 | real :: tau ! optical depth at 7 hpa |
---|
42 | |
---|
43 | !=================================================================== |
---|
44 | ! input/output: |
---|
45 | !=================================================================== |
---|
46 | |
---|
47 | real :: zycol(nlayer,nq) ! chemical species volume mixing ratio |
---|
48 | |
---|
49 | !=================================================================== |
---|
50 | ! output: |
---|
51 | !=================================================================== |
---|
52 | |
---|
53 | integer :: iter(nlayer) ! iteration counter |
---|
54 | real :: jo3(nlayer) ! photodissociation rate o3 -> o1d |
---|
55 | |
---|
56 | !=================================================================== |
---|
57 | ! local: |
---|
58 | !=================================================================== |
---|
59 | |
---|
60 | integer :: phychemrat ! (physical timestep)/(nominal chemical timestep) |
---|
61 | integer :: j_o3_o1d, ilev |
---|
62 | integer :: iesp, nesp |
---|
63 | integer :: lswitch |
---|
64 | |
---|
65 | logical, save :: firstcall = .true. |
---|
66 | |
---|
67 | parameter (nesp = 17) ! number of species in the chemistry |
---|
68 | |
---|
69 | ! tracer indexes in the chemistry: |
---|
70 | |
---|
71 | integer,parameter :: i_co2 = 1 |
---|
72 | integer,parameter :: i_co = 2 |
---|
73 | integer,parameter :: i_o = 3 |
---|
74 | integer,parameter :: i_o1d = 4 |
---|
75 | integer,parameter :: i_o2 = 5 |
---|
76 | integer,parameter :: i_o3 = 6 |
---|
77 | integer,parameter :: i_h = 7 |
---|
78 | integer,parameter :: i_h2 = 8 |
---|
79 | integer,parameter :: i_oh = 9 |
---|
80 | integer,parameter :: i_ho2 = 10 |
---|
81 | integer,parameter :: i_h2o2 = 11 |
---|
82 | integer,parameter :: i_h2o = 12 |
---|
83 | integer,parameter :: i_n = 13 |
---|
84 | integer,parameter :: i_n2d = 14 |
---|
85 | integer,parameter :: i_no = 15 |
---|
86 | integer,parameter :: i_no2 = 16 |
---|
87 | integer,parameter :: i_n2 = 17 |
---|
88 | |
---|
89 | real :: stimestep ! standard timestep for the chemistry (s) |
---|
90 | real :: ctimestep ! real timestep for the chemistry (s) |
---|
91 | real :: dt_guess ! first-guess timestep (s) |
---|
92 | real :: dt_corrected ! corrected timestep (s) |
---|
93 | real :: dt_min = 1. ! minimum allowed timestep (s) |
---|
94 | real :: dtg ! correction factor for the timestep (s) |
---|
95 | real :: j(nlayer,nd) ! interpolated photolysis rates (s-1) |
---|
96 | real :: time ! internal time (between 0 and ptimestep, in s) |
---|
97 | |
---|
98 | real, dimension(nlayer,nesp) :: rm ! mixing ratios |
---|
99 | real (kind = 8), dimension(nesp) :: cold ! number densities at previous timestep (molecule.cm-3) |
---|
100 | real (kind = 8), dimension(nlayer,nesp) :: c ! number densities at current timestep (molecule.cm-3) |
---|
101 | real (kind = 8), dimension(nesp) :: cnew ! number densities at next timestep (molecule.cm-3) |
---|
102 | |
---|
103 | ! reaction rates |
---|
104 | |
---|
105 | real (kind = 8), dimension(nlayer, nb_phot_max) :: v_phot |
---|
106 | real (kind = 8), dimension(nlayer,nb_reaction_3_max) :: v_3 |
---|
107 | real (kind = 8), dimension(nlayer,nb_reaction_4_max) :: v_4 |
---|
108 | logical :: hetero_dust, hetero_ice |
---|
109 | |
---|
110 | ! matrix |
---|
111 | |
---|
112 | real (kind = 8), dimension(nesp,nesp) :: mat, mat1 |
---|
113 | integer, dimension(nesp) :: indx |
---|
114 | integer :: code |
---|
115 | |
---|
116 | ! production and loss terms (for first-guess solution only) |
---|
117 | |
---|
118 | real (kind = 8), dimension(nesp) :: prod, loss |
---|
119 | |
---|
120 | ! curvatures |
---|
121 | |
---|
122 | real :: ratio, curv, e, e1, e2, e3 |
---|
123 | |
---|
124 | !=================================================================== |
---|
125 | ! initialisation of the reaction indexes |
---|
126 | !=================================================================== |
---|
127 | |
---|
128 | if (firstcall) then |
---|
129 | print*,'photochemistry: initialize indexes' |
---|
130 | call indice(i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
131 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
132 | i_n, i_n2d, i_no, i_no2, i_n2) |
---|
133 | firstcall = .false. |
---|
134 | end if |
---|
135 | |
---|
136 | !=================================================================== |
---|
137 | ! initialisation of mixing ratios and densities |
---|
138 | !=================================================================== |
---|
139 | |
---|
140 | call gcmtochim(nlayer, nq, zycol, lswitch, nesp, & |
---|
141 | i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
142 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
143 | i_n, i_n2d, i_no, i_no2, i_n2, dens, rm, c) |
---|
144 | |
---|
145 | !=================================================================== |
---|
146 | ! interpolation of photolysis rates in the lookup table |
---|
147 | !=================================================================== |
---|
148 | |
---|
149 | call photolysis_asis(nlayer, ngrid, lswitch, press, temp, sza, fractcol, tau, zmmean, dist_sol, & |
---|
150 | rm(:,i_co2), rm(:,i_o3), v_phot) |
---|
151 | |
---|
152 | ! save o3 photolysis for output |
---|
153 | |
---|
154 | j_o3_o1d = 5 |
---|
155 | jo3(:) = v_phot(:,j_o3_o1d) |
---|
156 | |
---|
157 | !=================================================================== |
---|
158 | ! reaction rates |
---|
159 | !=================================================================== |
---|
160 | ! switches for heterogeneous chemistry |
---|
161 | ! hetero_ice : reactions on ice clouds |
---|
162 | ! hetero_dust : reactions on dust |
---|
163 | !=================================================================== |
---|
164 | |
---|
165 | hetero_dust = .false. |
---|
166 | hetero_ice = .false. |
---|
167 | |
---|
168 | call reactionrates(nlayer, lswitch, dens, c(:,i_co2), c(:,i_o2), & |
---|
169 | press, temp, hetero_dust, hetero_ice, & |
---|
170 | surfdust1d, surfice1d, v_phot, v_3, v_4) |
---|
171 | |
---|
172 | !=================================================================== |
---|
173 | ! stimestep : standard chemical timestep (s) |
---|
174 | ! ctimestep : real chemical timestep (s), |
---|
175 | ! taking into account the physical timestep |
---|
176 | !=================================================================== |
---|
177 | |
---|
178 | stimestep = 600. ! standard value : 10 mn |
---|
179 | |
---|
180 | phychemrat = nint(ptimestep/stimestep) |
---|
181 | phychemrat = 1 |
---|
182 | |
---|
183 | ctimestep = ptimestep/real(phychemrat) |
---|
184 | |
---|
185 | !print*, "stimestep = ", stimestep |
---|
186 | !print*, "ptimestep = ", ptimestep |
---|
187 | !print*, "phychemrat = ", phychemrat |
---|
188 | !print*, "ctimestep = ", ctimestep |
---|
189 | !stop |
---|
190 | |
---|
191 | !=================================================================== |
---|
192 | ! loop over levels |
---|
193 | !=================================================================== |
---|
194 | |
---|
195 | do ilev = 1,lswitch - 1 |
---|
196 | |
---|
197 | ! initializations |
---|
198 | |
---|
199 | time = 0. |
---|
200 | iter(ilev) = 0 |
---|
201 | dt_guess = ctimestep |
---|
202 | cold(:) = c(ilev,:) |
---|
203 | |
---|
204 | ! internal loop for the chemistry |
---|
205 | |
---|
206 | do while (time < ptimestep) |
---|
207 | |
---|
208 | iter(ilev) = iter(ilev) + 1 ! iteration counter |
---|
209 | |
---|
210 | ! first-guess: fill matrix |
---|
211 | |
---|
212 | call fill_matrix(ilev, mat1, prod, loss, c, nesp, nlayer, v_phot, v_3, v_4) |
---|
213 | |
---|
214 | ! adaptative evaluation of the sub time step |
---|
215 | |
---|
216 | call define_dt(nesp, dt_corrected, dt_guess, ctimestep, cold(:), c(ilev,:), & |
---|
217 | mat1, prod, loss, dens(ilev)) |
---|
218 | |
---|
219 | if (time + dt_corrected > ptimestep) then |
---|
220 | dt_corrected = ptimestep - time |
---|
221 | end if |
---|
222 | |
---|
223 | ! if (dt_corrected /= dt_guess) then ! the timestep has been modified |
---|
224 | |
---|
225 | ! form the matrix identity + mat*dt_corrected |
---|
226 | |
---|
227 | mat(:,:) = mat1(:,:)*dt_corrected |
---|
228 | do iesp = 1,nesp |
---|
229 | mat(iesp,iesp) = 1. + mat(iesp,iesp) |
---|
230 | end do |
---|
231 | |
---|
232 | ! solve the linear system M*Cn+1 = Cn (RHS in cnew, then replaced by solution) |
---|
233 | |
---|
234 | cnew(:) = c(ilev,:) |
---|
235 | |
---|
236 | #ifdef LAPACK |
---|
237 | call dgesv(nesp,1,mat,nesp,indx,cnew,nesp,code) |
---|
238 | #else |
---|
239 | write(*,*) "photochemistry_asis error, missing LAPACK routine dgesv" |
---|
240 | stop |
---|
241 | #endif |
---|
242 | |
---|
243 | ! end if |
---|
244 | |
---|
245 | ! eliminate small values |
---|
246 | |
---|
247 | where (cnew(:)/dens(ilev) < 1.e-30) |
---|
248 | cnew(:) = 0. |
---|
249 | end where |
---|
250 | |
---|
251 | ! update concentrations |
---|
252 | |
---|
253 | cold(:) = c(ilev,:) |
---|
254 | c(ilev,:) = cnew(:) |
---|
255 | cnew(:) = 0. |
---|
256 | |
---|
257 | ! increment internal time |
---|
258 | |
---|
259 | time = time + dt_corrected |
---|
260 | dt_guess = dt_corrected ! first-guess timestep for next iteration |
---|
261 | |
---|
262 | end do ! while (time < ptimestep) |
---|
263 | |
---|
264 | end do ! ilev |
---|
265 | |
---|
266 | !=================================================================== |
---|
267 | ! save chemical species for the gcm |
---|
268 | !=================================================================== |
---|
269 | |
---|
270 | call chimtogcm(nlayer, nq, zycol, lswitch, nesp, & |
---|
271 | i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
272 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
273 | i_n, i_n2d, i_no, i_no2, i_n2, dens, c) |
---|
274 | |
---|
275 | contains |
---|
276 | |
---|
277 | !================================================================ |
---|
278 | |
---|
279 | subroutine define_dt(nesp, dtnew, dtold, ctimestep, cold, ccur, mat1, & |
---|
280 | prod, loss, dens) |
---|
281 | |
---|
282 | !================================================================ |
---|
283 | ! iterative evaluation of the appropriate time step dtnew |
---|
284 | ! according to curvature criterion based on |
---|
285 | ! e = 2 Rtol [r Cn+1 -(1-r) Cn + Cn-1 ]/[(1+r) Cn] |
---|
286 | ! with r = (tn - tn-1)/(tn+1 - tn) |
---|
287 | !================================================================ |
---|
288 | |
---|
289 | implicit none |
---|
290 | |
---|
291 | ! input |
---|
292 | |
---|
293 | integer :: nesp ! number of species in the chemistry |
---|
294 | |
---|
295 | real :: dtold, ctimestep |
---|
296 | real (kind = 8), dimension(nesp) :: cold, ccur |
---|
297 | real (kind = 8), dimension(nesp,nesp) :: mat1 |
---|
298 | real (kind = 8), dimension(nesp) :: prod, loss |
---|
299 | real :: dens |
---|
300 | |
---|
301 | ! output |
---|
302 | |
---|
303 | real :: dtnew |
---|
304 | |
---|
305 | ! local |
---|
306 | |
---|
307 | real (kind = 8), dimension(nesp) :: cnew |
---|
308 | real (kind = 8), dimension(nesp,nesp) :: mat |
---|
309 | real (kind = 8) :: atol, ratio, e, es, coef |
---|
310 | |
---|
311 | integer :: code, iesp, iter |
---|
312 | integer, dimension(nesp) :: indx |
---|
313 | |
---|
314 | real :: dttest |
---|
315 | |
---|
316 | ! parameters |
---|
317 | |
---|
318 | real (kind = 8), parameter :: dtmin = 10. ! minimum time step (s) |
---|
319 | real (kind = 8), parameter :: vmrtol = 1.e-11 ! absolute tolerance on vmr |
---|
320 | real (kind = 8), parameter :: rtol = 1./0.05 ! 1/rtol recommended value : 0.1-0.02 |
---|
321 | integer, parameter :: niter = 3 ! number of iterations |
---|
322 | real (kind = 8), parameter :: coefmax = 2. |
---|
323 | real (kind = 8), parameter :: coefmin = 0.1 |
---|
324 | logical :: fast_guess = .true. |
---|
325 | |
---|
326 | dttest = dtold ! dttest = dtold = dt_guess |
---|
327 | |
---|
328 | atol = vmrtol*dens ! absolute tolerance in molecule.cm-3 |
---|
329 | |
---|
330 | do iter = 1,niter |
---|
331 | |
---|
332 | if (fast_guess) then |
---|
333 | |
---|
334 | ! first guess : fast semi-implicit method |
---|
335 | |
---|
336 | do iesp = 1, nesp |
---|
337 | cnew(iesp) = (ccur(iesp) + prod(iesp)*dttest)/(1. + loss(iesp)*dttest) |
---|
338 | end do |
---|
339 | |
---|
340 | else |
---|
341 | |
---|
342 | ! first guess : form the matrix identity + mat*dt_guess |
---|
343 | |
---|
344 | mat(:,:) = mat1(:,:)*dttest |
---|
345 | do iesp = 1,nesp |
---|
346 | mat(iesp,iesp) = 1. + mat(iesp,iesp) |
---|
347 | end do |
---|
348 | |
---|
349 | ! form right-hand side (RHS) of the system |
---|
350 | |
---|
351 | cnew(:) = ccur(:) |
---|
352 | |
---|
353 | ! solve the linear system M*Cn+1 = Cn (RHS in cnew, then replaced by solution) |
---|
354 | |
---|
355 | #ifdef LAPACK |
---|
356 | call dgesv(nesp,1,mat,nesp,indx,cnew,nesp,code) |
---|
357 | #else |
---|
358 | write(*,*) "photochemistry_asis error, missing LAPACK routine dgesv" |
---|
359 | stop |
---|
360 | #endif |
---|
361 | |
---|
362 | end if |
---|
363 | |
---|
364 | ! ratio old/new subtimestep |
---|
365 | |
---|
366 | ratio = dtold/dttest |
---|
367 | |
---|
368 | ! e : local error indicatocitr |
---|
369 | |
---|
370 | e = 0. |
---|
371 | |
---|
372 | do iesp = 1,nesp |
---|
373 | es = 2.*abs((ratio*cnew(iesp) - (1. + ratio)*ccur(iesp) + cold(iesp)) & |
---|
374 | /(1. + ratio)/max(ccur(iesp),atol)) |
---|
375 | |
---|
376 | if (es > e) then |
---|
377 | e = es |
---|
378 | end if |
---|
379 | end do |
---|
380 | e = rtol*e |
---|
381 | |
---|
382 | ! timestep correction |
---|
383 | |
---|
384 | coef = max(coefmin, min(coefmax,0.8/sqrt(e))) |
---|
385 | |
---|
386 | dttest = max(dtmin,dttest*coef) |
---|
387 | dttest = min(ctimestep,dttest) |
---|
388 | |
---|
389 | end do ! iter |
---|
390 | |
---|
391 | ! new timestep |
---|
392 | |
---|
393 | dtnew = dttest |
---|
394 | |
---|
395 | end subroutine define_dt |
---|
396 | |
---|
397 | |
---|
398 | |
---|
399 | |
---|
400 | !====================================================================== |
---|
401 | |
---|
402 | subroutine reactionrates(nlayer, & |
---|
403 | lswitch, dens, co2, o2, press, t, & |
---|
404 | hetero_dust, hetero_ice, & |
---|
405 | surfdust1d, surfice1d, & |
---|
406 | v_phot, v_3, v_4) |
---|
407 | |
---|
408 | !================================================================ |
---|
409 | ! compute reaction rates ! |
---|
410 | !---------------------------------------------------------------- |
---|
411 | ! reaction type array ! |
---|
412 | !---------------------------------------------------------------- |
---|
413 | ! A + B --> C + D bimolecular v_4 ! |
---|
414 | ! A + A --> B + C quadratic v_3 ! |
---|
415 | ! A + C --> B + C quenching v_phot ! |
---|
416 | ! A + ice --> B + C heterogeneous v_phot ! |
---|
417 | !================================================================ |
---|
418 | |
---|
419 | use comcstfi_mod |
---|
420 | |
---|
421 | implicit none |
---|
422 | |
---|
423 | #include "chimiedata.h" |
---|
424 | |
---|
425 | !---------------------------------------------------------------------- |
---|
426 | ! input |
---|
427 | !---------------------------------------------------------------------- |
---|
428 | |
---|
429 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
430 | integer :: lswitch ! interface level between lower |
---|
431 | ! atmosphere and thermosphere chemistries |
---|
432 | real, dimension(nlayer) :: dens ! total number density (molecule.cm-3) |
---|
433 | real, dimension(nlayer) :: press ! pressure (hPa) |
---|
434 | real, dimension(nlayer) :: t ! temperature (K) |
---|
435 | real, dimension(nlayer) :: surfdust1d ! dust surface area (cm2.cm-3) |
---|
436 | real, dimension(nlayer) :: surfice1d ! ice surface area (cm2.cm-3) |
---|
437 | real (kind = 8), dimension(nlayer) :: co2 ! co2 number density (molecule.cm-3) |
---|
438 | real (kind = 8), dimension(nlayer) :: o2 ! o2 number density (molecule.cm-3) |
---|
439 | logical :: hetero_dust, hetero_ice ! switches for heterogeneous chemistry |
---|
440 | |
---|
441 | !---------------------------------------------------------------------- |
---|
442 | ! output |
---|
443 | !---------------------------------------------------------------------- |
---|
444 | |
---|
445 | real (kind = 8), dimension(nlayer, nb_phot_max) :: v_phot |
---|
446 | real (kind = 8), dimension(nlayer,nb_reaction_3_max) :: v_3 |
---|
447 | real (kind = 8), dimension(nlayer,nb_reaction_4_max) :: v_4 |
---|
448 | |
---|
449 | !---------------------------------------------------------------------- |
---|
450 | ! local |
---|
451 | !---------------------------------------------------------------------- |
---|
452 | |
---|
453 | integer :: ilev |
---|
454 | integer :: nb_phot, nb_reaction_3, nb_reaction_4 |
---|
455 | real :: ak0, ak1, xpo, rate |
---|
456 | real :: k1a0, k1b0, k1ainf, k1a, k1b, fc, fx, x, y, gam |
---|
457 | real, dimension(nlayer) :: deq |
---|
458 | real, dimension(nlayer) :: a001, a002, a003, & |
---|
459 | b001, b002, b003, b004, b005, b006, b007, & |
---|
460 | b008, b009, & |
---|
461 | c001, c002, c003, c004, c005, c006, c007, & |
---|
462 | c008, c009, c010, c011, c012, c013, c014, & |
---|
463 | c015, c016, c017, c018, & |
---|
464 | d001, d002, d003, d004, d005, d006, d007, & |
---|
465 | d008, d009, & |
---|
466 | e001, e002, & |
---|
467 | h001, h002, h003, h004, h005 |
---|
468 | |
---|
469 | !---------------------------------------------------------------------- |
---|
470 | ! initialisation |
---|
471 | !---------------------------------------------------------------------- |
---|
472 | |
---|
473 | nb_phot = 11 ! jmars.20140930 reduit de 13 a 11 |
---|
474 | nb_reaction_3 = 0 |
---|
475 | nb_reaction_4 = 0 |
---|
476 | |
---|
477 | !---------------------------------------------------------------------- |
---|
478 | ! oxygen reactions |
---|
479 | !---------------------------------------------------------------------- |
---|
480 | |
---|
481 | !--- a001: o + o2 + co2 -> o3 + co2 |
---|
482 | |
---|
483 | ! jpl 2003 |
---|
484 | ! |
---|
485 | ! co2 efficiency as a third body (2.075) |
---|
486 | ! from sehested et al., j. geophys. res., 100, 1995. |
---|
487 | |
---|
488 | a001(:) = 2.075*6.0e-34*(t(:)/300.)**(-2.4)*dens(:) |
---|
489 | |
---|
490 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
491 | v_4(:,nb_reaction_4) = a001(:) |
---|
492 | |
---|
493 | !--- a002: o + o + co2 -> o2 + co2 |
---|
494 | |
---|
495 | ! Tsang and Hampson, J. Chem. Phys. Ref. Data, 15, 1087, 1986 |
---|
496 | |
---|
497 | ! a002(:) = 2.5*5.2e-35*exp(900./t(:))*dens(:) |
---|
498 | |
---|
499 | ! Campbell and Gray, Chem. Phys. Lett., 18, 607, 1973 |
---|
500 | |
---|
501 | ! a002(:) = 1.2e-32*(300./t(:))**(2.0)*dens(:) ! yung expression |
---|
502 | a002(:) = 2.5*9.46e-34*exp(485./t(:))*dens(:) ! nist expression |
---|
503 | |
---|
504 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
505 | v_3(:,nb_reaction_3) = a002(:) |
---|
506 | |
---|
507 | !--- a003: o + o3 -> o2 + o2 |
---|
508 | |
---|
509 | ! jpl 2003 |
---|
510 | |
---|
511 | a003(:) = 8.0e-12*exp(-2060./t(:)) |
---|
512 | |
---|
513 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
514 | v_4(:,nb_reaction_4) = a003(:) |
---|
515 | |
---|
516 | !---------------------------------------------------------------------- |
---|
517 | ! o(1d) reactions |
---|
518 | !---------------------------------------------------------------------- |
---|
519 | |
---|
520 | !--- b001: o(1d) + co2 -> o + co2 |
---|
521 | |
---|
522 | ! jpl 2006 |
---|
523 | |
---|
524 | b001(:) = 7.5e-11*exp(115./t(:)) |
---|
525 | |
---|
526 | nb_phot = nb_phot + 1 |
---|
527 | v_phot(:,nb_phot) = b001(:)*co2(:) |
---|
528 | |
---|
529 | !--- b002: o(1d) + h2o -> oh + oh |
---|
530 | |
---|
531 | ! jpl 2006 |
---|
532 | |
---|
533 | b002(:) = 1.63e-10*exp(60./t(:)) |
---|
534 | |
---|
535 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
536 | v_4(:,nb_reaction_4) = b002(:) |
---|
537 | |
---|
538 | !--- b003: o(1d) + h2 -> oh + h |
---|
539 | |
---|
540 | ! jpl 2011 |
---|
541 | |
---|
542 | b003(:) = 1.2e-10 |
---|
543 | |
---|
544 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
545 | v_4(:,nb_reaction_4) = b003(:) |
---|
546 | |
---|
547 | !--- b004: o(1d) + o2 -> o + o2 |
---|
548 | |
---|
549 | ! jpl 2006 |
---|
550 | |
---|
551 | b004(:) = 3.3e-11*exp(55./t(:)) |
---|
552 | |
---|
553 | nb_phot = nb_phot + 1 |
---|
554 | v_phot(:,nb_phot) = b004(:)*o2(:) |
---|
555 | |
---|
556 | !--- b005: o(1d) + o3 -> o2 + o2 |
---|
557 | |
---|
558 | ! jpl 2003 |
---|
559 | |
---|
560 | b005(:) = 1.2e-10 |
---|
561 | |
---|
562 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
563 | v_4(:,nb_reaction_4) = b005(:) |
---|
564 | |
---|
565 | !--- b006: o(1d) + o3 -> o2 + o + o |
---|
566 | |
---|
567 | ! jpl 2003 |
---|
568 | |
---|
569 | b006(:) = 1.2e-10 |
---|
570 | |
---|
571 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
572 | v_4(:,nb_reaction_4) = b006(:) |
---|
573 | |
---|
574 | !--- b007: o(1d) + ch4 -> ch3 + oh |
---|
575 | |
---|
576 | ! jpl 2003 |
---|
577 | |
---|
578 | b007(:) = 1.5e-10*0.75 |
---|
579 | |
---|
580 | !--- b008: o(1d) + ch4 -> ch3o + h |
---|
581 | |
---|
582 | ! jpl 2003 |
---|
583 | |
---|
584 | b008(:) = 1.5e-10*0.20 |
---|
585 | ! |
---|
586 | !--- b009: o(1d) + ch4 -> ch2o + h2 |
---|
587 | |
---|
588 | ! jpl 2003 |
---|
589 | |
---|
590 | b009(:) = 1.5e-10*0.05 |
---|
591 | |
---|
592 | !---------------------------------------------------------------------- |
---|
593 | ! hydrogen reactions |
---|
594 | !---------------------------------------------------------------------- |
---|
595 | |
---|
596 | !--- c001: o + ho2 -> oh + o2 |
---|
597 | |
---|
598 | ! jpl 2003 |
---|
599 | |
---|
600 | c001(:) = 3.0e-11*exp(200./t(:)) |
---|
601 | |
---|
602 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
603 | v_4(:,nb_reaction_4) = c001(:) |
---|
604 | |
---|
605 | !--- c002: o + oh -> o2 + h |
---|
606 | |
---|
607 | ! jpl 2011 |
---|
608 | |
---|
609 | c002(:) = 1.8e-11*exp(180./t(:)) |
---|
610 | |
---|
611 | ! robertson and smith, j. chem. phys. a 110, 6673, 2006 |
---|
612 | |
---|
613 | ! c002(:) = 11.2e-11*t(:)**(-0.32)*exp(177./t(:)) |
---|
614 | |
---|
615 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
616 | v_4(:,nb_reaction_4) = c002(:) |
---|
617 | |
---|
618 | !--- c003: h + o3 -> oh + o2 |
---|
619 | |
---|
620 | ! jpl 2003 |
---|
621 | |
---|
622 | c003(:) = 1.4e-10*exp(-470./t(:)) |
---|
623 | |
---|
624 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
625 | v_4(:,nb_reaction_4) = c003(:) |
---|
626 | |
---|
627 | !--- c004: h + ho2 -> oh + oh |
---|
628 | |
---|
629 | ! jpl 2006 |
---|
630 | |
---|
631 | c004(:) = 7.2e-11 |
---|
632 | |
---|
633 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
634 | v_4(:,nb_reaction_4) = c004(:) |
---|
635 | |
---|
636 | !--- c005: h + ho2 -> h2 + o2 |
---|
637 | |
---|
638 | ! jpl 2006 |
---|
639 | |
---|
640 | c005(:) = 6.9e-12 |
---|
641 | |
---|
642 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
643 | v_4(:,nb_reaction_4) = c005(:) |
---|
644 | |
---|
645 | !--- c006: h + ho2 -> h2o + o |
---|
646 | |
---|
647 | ! jpl 2006 |
---|
648 | |
---|
649 | c006(:) = 1.6e-12 |
---|
650 | |
---|
651 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
652 | v_4(:,nb_reaction_4) = c006(:) |
---|
653 | |
---|
654 | !--- c007: oh + ho2 -> h2o + o2 |
---|
655 | |
---|
656 | ! jpl 2003 |
---|
657 | |
---|
658 | ! canty et al., grl, 2006 suggest to increase this rate |
---|
659 | ! by 20%. not done here. |
---|
660 | |
---|
661 | c007(:) = 4.8e-11*exp(250./t(:)) |
---|
662 | |
---|
663 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
664 | v_4(:,nb_reaction_4) = c007(:) |
---|
665 | |
---|
666 | !--- c008: ho2 + ho2 -> h2o2 + o2 |
---|
667 | |
---|
668 | ! jpl 2006 |
---|
669 | |
---|
670 | ! c008(:) = 3.5e-13*exp(430./t(:)) |
---|
671 | |
---|
672 | ! christensen et al., grl, 13, 2002 |
---|
673 | |
---|
674 | c008(:) = 1.5e-12*exp(19./t(:)) |
---|
675 | |
---|
676 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
677 | v_3(:,nb_reaction_3) = c008(:) |
---|
678 | |
---|
679 | !--- c009: oh + h2o2 -> h2o + ho2 |
---|
680 | |
---|
681 | ! jpl 2006 |
---|
682 | |
---|
683 | c009(:) = 1.8e-12 |
---|
684 | |
---|
685 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
686 | v_4(:,nb_reaction_4) = c009(:) |
---|
687 | |
---|
688 | !--- c010: oh + h2 -> h2o + h |
---|
689 | |
---|
690 | ! jpl 2006 |
---|
691 | |
---|
692 | c010(:) = 2.8e-12*exp(-1800./t(:)) |
---|
693 | |
---|
694 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
695 | v_4(:,nb_reaction_4) = c010(:) |
---|
696 | |
---|
697 | !--- c011: h + o2 + co2 -> ho2 + co2 |
---|
698 | |
---|
699 | ! jpl 2011 |
---|
700 | |
---|
701 | do ilev = 1,lswitch-1 |
---|
702 | ak0 = 2.5*4.4e-32*(t(ilev)/300.)**(-1.3) |
---|
703 | ak1 = 7.5e-11*(t(ilev)/300.)**(0.2) |
---|
704 | |
---|
705 | rate = (ak0*dens(ilev))/(1. + ak0*dens(ilev)/ak1) |
---|
706 | xpo = 1./(1. + alog10((ak0*dens(ilev))/ak1)**2) |
---|
707 | c011(ilev) = rate*0.6**xpo |
---|
708 | end do |
---|
709 | |
---|
710 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
711 | v_4(:,nb_reaction_4) = c011(:) |
---|
712 | |
---|
713 | !--- c012: o + h2o2 -> oh + ho2 |
---|
714 | |
---|
715 | ! jpl 2003 |
---|
716 | |
---|
717 | c012(:) = 1.4e-12*exp(-2000./t(:)) |
---|
718 | |
---|
719 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
720 | v_4(:,nb_reaction_4) = c012(:) |
---|
721 | |
---|
722 | !--- c013: oh + oh -> h2o + o |
---|
723 | |
---|
724 | ! jpl 2006 |
---|
725 | |
---|
726 | c013(:) = 1.8e-12 |
---|
727 | |
---|
728 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
729 | v_3(:,nb_reaction_3) = c013(:) |
---|
730 | |
---|
731 | !--- c014: oh + o3 -> ho2 + o2 |
---|
732 | |
---|
733 | ! jpl 2003 |
---|
734 | |
---|
735 | c014(:) = 1.7e-12*exp(-940./t(:)) |
---|
736 | |
---|
737 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
738 | v_4(:,nb_reaction_4) = c014(:) |
---|
739 | |
---|
740 | !--- c015: ho2 + o3 -> oh + o2 + o2 |
---|
741 | |
---|
742 | ! jpl 2003 |
---|
743 | |
---|
744 | c015(:) = 1.0e-14*exp(-490./t(:)) |
---|
745 | |
---|
746 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
747 | v_4(:,nb_reaction_4) = c015(:) |
---|
748 | |
---|
749 | !--- c016: ho2 + ho2 + co2 -> h2o2 + o2 + co2 |
---|
750 | |
---|
751 | ! jpl 2011 |
---|
752 | |
---|
753 | c016(:) = 2.5*2.1e-33*exp(920./t(:))*dens(:) |
---|
754 | |
---|
755 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
756 | v_3(:,nb_reaction_3) = c016(:) |
---|
757 | |
---|
758 | !--- c017: oh + oh + co2 -> h2o2 + co2 |
---|
759 | |
---|
760 | ! jpl 2003 |
---|
761 | |
---|
762 | do ilev = 1,lswitch-1 |
---|
763 | ak0 = 2.5*6.9e-31*(t(ilev)/300.)**(-1.0) |
---|
764 | ak1 = 2.6e-11*(t(ilev)/300.)**(0.0) |
---|
765 | |
---|
766 | rate = (ak0*dens(ilev))/(1. + ak0*dens(ilev)/ak1) |
---|
767 | xpo = 1./(1. + alog10((ak0*dens(ilev))/ak1)**2) |
---|
768 | c017(ilev) = rate*0.6**xpo |
---|
769 | end do |
---|
770 | |
---|
771 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
772 | v_3(:,nb_reaction_3) = c017(:) |
---|
773 | |
---|
774 | !--- c018: h + h + co2 -> h2 + co2 |
---|
775 | |
---|
776 | ! baulch et al., 2005 |
---|
777 | |
---|
778 | c018(:) = 2.5*1.8e-30*(t(:)**(-1.0))*dens(:) |
---|
779 | |
---|
780 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
781 | v_3(:,nb_reaction_3) = c018(:) |
---|
782 | |
---|
783 | !---------------------------------------------------------------------- |
---|
784 | ! nitrogen reactions |
---|
785 | !---------------------------------------------------------------------- |
---|
786 | |
---|
787 | !--- d001: no2 + o -> no + o2 |
---|
788 | |
---|
789 | ! jpl 2006 |
---|
790 | |
---|
791 | d001(:) = 5.1e-12*exp(210./t(:)) |
---|
792 | |
---|
793 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
794 | v_4(:,nb_reaction_4) = d001(:) |
---|
795 | |
---|
796 | !--- d002: no + o3 -> no2 + o2 |
---|
797 | |
---|
798 | ! jpl 2006 |
---|
799 | |
---|
800 | d002(:) = 3.0e-12*exp(-1500./t(:)) |
---|
801 | |
---|
802 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
803 | v_4(:,nb_reaction_4) = d002(:) |
---|
804 | |
---|
805 | !--- d003: no + ho2 -> no2 + oh |
---|
806 | |
---|
807 | ! jpl 2011 |
---|
808 | |
---|
809 | d003(:) = 3.3e-12*exp(270./t(:)) |
---|
810 | |
---|
811 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
812 | v_4(:,nb_reaction_4) = d003(:) |
---|
813 | |
---|
814 | !--- d004: n + no -> n2 + o |
---|
815 | |
---|
816 | ! jpl 2011 |
---|
817 | |
---|
818 | d004(:) = 2.1e-11*exp(100./t(:)) |
---|
819 | |
---|
820 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
821 | v_4(:,nb_reaction_4) = d004(:) |
---|
822 | |
---|
823 | !--- d005: n + o2 -> no + o |
---|
824 | |
---|
825 | ! jpl 2011 |
---|
826 | |
---|
827 | d005(:) = 1.5e-11*exp(-3600./t(:)) |
---|
828 | |
---|
829 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
830 | v_4(:,nb_reaction_4) = d005(:) |
---|
831 | |
---|
832 | !--- d006: no2 + h -> no + oh |
---|
833 | |
---|
834 | ! jpl 2011 |
---|
835 | |
---|
836 | d006(:) = 4.0e-10*exp(-340./t(:)) |
---|
837 | |
---|
838 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
839 | v_4(:,nb_reaction_4) = d006(:) |
---|
840 | |
---|
841 | !--- d007: n + o -> no |
---|
842 | |
---|
843 | d007(:) = 2.8e-17*(300./t(:))**0.5 |
---|
844 | |
---|
845 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
846 | v_4(:,nb_reaction_4) = d007(:) |
---|
847 | |
---|
848 | !--- d008: n + ho2 -> no + oh |
---|
849 | |
---|
850 | ! brune et al., j. chem. phys., 87, 1983 |
---|
851 | |
---|
852 | d008(:) = 2.19e-11 |
---|
853 | |
---|
854 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
855 | v_4(:,nb_reaction_4) = d008(:) |
---|
856 | |
---|
857 | !--- d009: n + oh -> no + h |
---|
858 | |
---|
859 | ! atkinson et al., j. phys. chem. ref. data, 18, 881, 1989 |
---|
860 | |
---|
861 | d009(:) = 3.8e-11*exp(85./t(:)) |
---|
862 | |
---|
863 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
864 | v_4(:,nb_reaction_4) = d009(:) |
---|
865 | |
---|
866 | !---------------------------------------------------------------------- |
---|
867 | ! carbon reactions |
---|
868 | !---------------------------------------------------------------------- |
---|
869 | |
---|
870 | !--- e001: oh + co -> co2 + h |
---|
871 | |
---|
872 | ! jpl 2003 |
---|
873 | |
---|
874 | ! e001(:) = 1.5e-13*(1 + 0.6*press(:)/1013.) |
---|
875 | |
---|
876 | ! mccabe et al., grl, 28, 3135, 2001 |
---|
877 | |
---|
878 | ! e001(:) = 1.57e-13 + 3.54e-33*dens(:) |
---|
879 | |
---|
880 | ! jpl 2006 |
---|
881 | |
---|
882 | ! ak0 = 1.5e-13*(t(:)/300.)**(0.6) |
---|
883 | ! ak1 = 2.1e-9*(t(:)/300.)**(6.1) |
---|
884 | ! rate1 = ak0/(1. + ak0/(ak1/dens(:))) |
---|
885 | ! xpo1 = 1./(1. + alog10(ak0/(ak1/dens(:)))**2) |
---|
886 | |
---|
887 | ! ak0 = 5.9e-33*(t(:)/300.)**(-1.4) |
---|
888 | ! ak1 = 1.1e-12*(t(:)/300.)**(1.3) |
---|
889 | ! rate2 = (ak0*dens(:))/(1. + ak0*dens(:)/ak1) |
---|
890 | ! xpo2 = 1./(1. + alog10((ak0*dens(:))/ak1)**2) |
---|
891 | |
---|
892 | ! e001(:) = rate1*0.6**xpo1 + rate2*0.6**xpo2 |
---|
893 | |
---|
894 | ! joshi et al., 2006 |
---|
895 | |
---|
896 | do ilev = 1,lswitch-1 |
---|
897 | k1a0 = 1.34*2.5*dens(ilev) & |
---|
898 | *1/(1/(3.62e-26*t(ilev)**(-2.739)*exp(-20./t(ilev))) & |
---|
899 | + 1/(6.48e-33*t(ilev)**(0.14)*exp(-57./t(ilev)))) ! typo in paper corrected |
---|
900 | k1b0 = 1.17e-19*t(ilev)**(2.053)*exp(139./t(ilev)) & |
---|
901 | + 9.56e-12*t(ilev)**(-0.664)*exp(-167./t(ilev)) |
---|
902 | k1ainf = 1.52e-17*t(ilev)**(1.858)*exp(28.8/t(ilev)) & |
---|
903 | + 4.78e-8*t(ilev)**(-1.851)*exp(-318./t(ilev)) |
---|
904 | x = k1a0/(k1ainf - k1b0) |
---|
905 | y = k1b0/(k1ainf - k1b0) |
---|
906 | fc = 0.628*exp(-1223./t(ilev)) + (1. - 0.628)*exp(-39./t(ilev)) & |
---|
907 | + exp(-t(ilev)/255.) |
---|
908 | fx = fc**(1./(1. + (alog(x))**2)) ! typo in paper corrected |
---|
909 | k1a = k1a0*((1. + y)/(1. + x))*fx |
---|
910 | k1b = k1b0*(1./(1.+x))*fx |
---|
911 | |
---|
912 | e001(ilev) = k1a + k1b |
---|
913 | end do |
---|
914 | |
---|
915 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
916 | v_4(:,nb_reaction_4) = e001(:) |
---|
917 | |
---|
918 | !--- e002: o + co + m -> co2 + m |
---|
919 | |
---|
920 | ! tsang and hampson, 1986. |
---|
921 | |
---|
922 | e002(:) = 2.5*6.5e-33*exp(-2184./t(:))*dens(:) |
---|
923 | |
---|
924 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
925 | v_4(:,nb_reaction_4) = e002(:) |
---|
926 | |
---|
927 | !---------------------------------------------------------------------- |
---|
928 | ! heterogeneous chemistry |
---|
929 | !---------------------------------------------------------------------- |
---|
930 | |
---|
931 | if (hetero_ice) then |
---|
932 | |
---|
933 | ! k = (surface*v*gamma)/4 (s-1) |
---|
934 | ! v = 100*sqrt(8rt/(pi*m)) (cm s-1) |
---|
935 | |
---|
936 | !--- h001: ho2 + ice -> products |
---|
937 | |
---|
938 | ! cooper and abbatt, 1996: gamma = 0.025 |
---|
939 | |
---|
940 | gam = 0.025 |
---|
941 | h001(:) = surfice1d(:) & |
---|
942 | *100.*sqrt(8.*8.31*t(:)/(33.e-3*pi))*gam/4. |
---|
943 | |
---|
944 | ! h002: oh + ice -> products |
---|
945 | |
---|
946 | ! cooper and abbatt, 1996: gamma = 0.03 |
---|
947 | |
---|
948 | gam = 0.03 |
---|
949 | h002(:) = surfice1d(:) & |
---|
950 | *100.*sqrt(8.*8.31*t(:)/(17.e-3*pi))*gam/4. |
---|
951 | |
---|
952 | !--- h003: h2o2 + ice -> products |
---|
953 | |
---|
954 | ! gamma = 0. test value |
---|
955 | |
---|
956 | gam = 0. |
---|
957 | h003(:) = surfice1d(:) & |
---|
958 | *100.*sqrt(8.*8.31*t(:)/(34.e-3*pi))*gam/4. |
---|
959 | else |
---|
960 | h001(:) = 0. |
---|
961 | h002(:) = 0. |
---|
962 | h003(:) = 0. |
---|
963 | end if |
---|
964 | |
---|
965 | nb_phot = nb_phot + 1 |
---|
966 | v_phot(:,nb_phot) = h001(:) |
---|
967 | |
---|
968 | nb_phot = nb_phot + 1 |
---|
969 | v_phot(:,nb_phot) = h002(:) |
---|
970 | |
---|
971 | nb_phot = nb_phot + 1 |
---|
972 | v_phot(:,nb_phot) = h003(:) |
---|
973 | |
---|
974 | if (hetero_dust) then |
---|
975 | |
---|
976 | !--- h004: ho2 + dust -> products |
---|
977 | |
---|
978 | ! jacob, 2000: gamma = 0.2 |
---|
979 | ! see dereus et al., atm. chem. phys., 2005 |
---|
980 | |
---|
981 | gam = 0.2 |
---|
982 | h004(:) = surfdust1d(:) & |
---|
983 | *100.*sqrt(8.*8.31*t(:)/(33.e-3*pi))*gam/4. |
---|
984 | |
---|
985 | !--- h005: h2o2 + dust -> products |
---|
986 | |
---|
987 | ! gamma = 5.e-4 |
---|
988 | ! see dereus et al., atm. chem. phys., 2005 |
---|
989 | |
---|
990 | gam = 5.e-4 |
---|
991 | h005(:) = surfdust1d(:) & |
---|
992 | *100.*sqrt(8.*8.31*t(:)/(34.e-3*pi))*gam/4. |
---|
993 | else |
---|
994 | h004(:) = 0. |
---|
995 | h005(:) = 0. |
---|
996 | end if |
---|
997 | |
---|
998 | nb_phot = nb_phot + 1 |
---|
999 | v_phot(:,nb_phot) = h004(:) |
---|
1000 | |
---|
1001 | nb_phot = nb_phot + 1 |
---|
1002 | v_phot(:,nb_phot) = h005(:) |
---|
1003 | |
---|
1004 | end subroutine reactionrates |
---|
1005 | |
---|
1006 | !====================================================================== |
---|
1007 | |
---|
1008 | subroutine fill_matrix(ilev, mat, prod, loss, c, nesp, nlayer, v_phot, v_3, v_4) |
---|
1009 | |
---|
1010 | !====================================================================== |
---|
1011 | ! filling of the jacobian matrix |
---|
1012 | !====================================================================== |
---|
1013 | |
---|
1014 | use types_asis |
---|
1015 | |
---|
1016 | implicit none |
---|
1017 | |
---|
1018 | #include "chimiedata.h" |
---|
1019 | |
---|
1020 | ! input |
---|
1021 | |
---|
1022 | integer :: ilev ! level index |
---|
1023 | integer :: nesp ! number of species in the chemistry |
---|
1024 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
1025 | |
---|
1026 | real (kind = 8), dimension(nlayer,nesp) :: c ! number densities |
---|
1027 | real (kind = 8), dimension(nlayer, nb_phot_max) :: v_phot |
---|
1028 | real (kind = 8), dimension(nlayer,nb_reaction_3_max) :: v_3 |
---|
1029 | real (kind = 8), dimension(nlayer,nb_reaction_4_max) :: v_4 |
---|
1030 | |
---|
1031 | ! output |
---|
1032 | |
---|
1033 | real (kind = 8), dimension(nesp,nesp), intent(out) :: mat ! matrix |
---|
1034 | real (kind = 8), dimension(nesp), intent(out) :: prod, loss |
---|
1035 | |
---|
1036 | ! local |
---|
1037 | |
---|
1038 | integer :: iesp |
---|
1039 | integer :: ind_phot_2,ind_phot_4,ind_phot_6 |
---|
1040 | integer :: ind_3_2,ind_3_4,ind_3_6 |
---|
1041 | integer :: ind_4_2,ind_4_4,ind_4_6,ind_4_8 |
---|
1042 | integer :: iphot,i3,i4 |
---|
1043 | |
---|
1044 | real(kind = jprb) :: eps, eps_4 ! implicit/explicit coefficient |
---|
1045 | |
---|
1046 | ! initialisations |
---|
1047 | |
---|
1048 | mat(:,:) = 0. |
---|
1049 | prod(:) = 0. |
---|
1050 | loss(:) = 0. |
---|
1051 | |
---|
1052 | ! photodissociations |
---|
1053 | ! or reactions a + c -> b + c |
---|
1054 | ! or reactions a + ice -> b + c |
---|
1055 | |
---|
1056 | do iphot = 1,nb_phot_max |
---|
1057 | |
---|
1058 | ind_phot_2 = indice_phot(iphot)%z2 |
---|
1059 | ind_phot_4 = indice_phot(iphot)%z4 |
---|
1060 | ind_phot_6 = indice_phot(iphot)%z6 |
---|
1061 | |
---|
1062 | mat(ind_phot_2,ind_phot_2) = mat(ind_phot_2,ind_phot_2) + indice_phot(iphot)%z1*v_phot(ilev,iphot) |
---|
1063 | mat(ind_phot_4,ind_phot_2) = mat(ind_phot_4,ind_phot_2) - indice_phot(iphot)%z3*v_phot(ilev,iphot) |
---|
1064 | mat(ind_phot_6,ind_phot_2) = mat(ind_phot_6,ind_phot_2) - indice_phot(iphot)%z5*v_phot(ilev,iphot) |
---|
1065 | |
---|
1066 | loss(ind_phot_2) = loss(ind_phot_2) + indice_phot(iphot)%z1*v_phot(ilev,iphot) |
---|
1067 | prod(ind_phot_4) = prod(ind_phot_4) + indice_phot(iphot)%z3*v_phot(ilev,iphot)*c(ilev,ind_phot_2) |
---|
1068 | prod(ind_phot_6) = prod(ind_phot_6) + indice_phot(iphot)%z5*v_phot(ilev,iphot)*c(ilev,ind_phot_2) |
---|
1069 | |
---|
1070 | end do |
---|
1071 | |
---|
1072 | ! reactions a + a -> b + c |
---|
1073 | |
---|
1074 | do i3 = 1,nb_reaction_3_max |
---|
1075 | |
---|
1076 | ind_3_2 = indice_3(i3)%z2 |
---|
1077 | ind_3_4 = indice_3(i3)%z4 |
---|
1078 | ind_3_6 = indice_3(i3)%z6 |
---|
1079 | |
---|
1080 | mat(ind_3_2,ind_3_2) = mat(ind_3_2,ind_3_2) + indice_3(i3)%z1*v_3(ilev,i3)*c(ilev,ind_3_2) |
---|
1081 | mat(ind_3_4,ind_3_2) = mat(ind_3_4,ind_3_2) - indice_3(i3)%z3*v_3(ilev,i3)*c(ilev,ind_3_2) |
---|
1082 | mat(ind_3_6,ind_3_2) = mat(ind_3_6,ind_3_2) - indice_3(i3)%z5*v_3(ilev,i3)*c(ilev,ind_3_2) |
---|
1083 | |
---|
1084 | loss(ind_3_2) = loss(ind_3_2) + indice_3(i3)%z1*v_3(ilev,i3)*c(ilev,ind_3_2) |
---|
1085 | prod(ind_3_4) = prod(ind_3_4) + indice_3(i3)%z3*v_3(ilev,i3)*c(ilev,ind_3_2)*c(ilev,ind_3_2) |
---|
1086 | prod(ind_3_6) = prod(ind_3_6) + indice_3(i3)%z5*v_3(ilev,i3)*c(ilev,ind_3_2)*c(ilev,ind_3_2) |
---|
1087 | |
---|
1088 | end do |
---|
1089 | |
---|
1090 | ! reactions a + b -> c + d |
---|
1091 | |
---|
1092 | eps = 1.d-10 |
---|
1093 | |
---|
1094 | do i4 = 1,nb_reaction_4_max |
---|
1095 | |
---|
1096 | ind_4_2 = indice_4(i4)%z2 |
---|
1097 | ind_4_4 = indice_4(i4)%z4 |
---|
1098 | ind_4_6 = indice_4(i4)%z6 |
---|
1099 | ind_4_8 = indice_4(i4)%z8 |
---|
1100 | |
---|
1101 | eps_4 = abs(c(ilev,ind_4_2))/(abs(c(ilev,ind_4_2)) + abs(c(ilev,ind_4_4)) + eps) |
---|
1102 | eps_4 = min(eps_4,1.0_jprb) |
---|
1103 | |
---|
1104 | mat(ind_4_2,ind_4_2) = mat(ind_4_2,ind_4_2) + indice_4(i4)%z1*v_4(ilev,i4)*(1. - eps_4)*c(ilev,ind_4_4) |
---|
1105 | mat(ind_4_2,ind_4_4) = mat(ind_4_2,ind_4_4) + indice_4(i4)%z1*v_4(ilev,i4)*eps_4*c(ilev,ind_4_2) |
---|
1106 | mat(ind_4_4,ind_4_2) = mat(ind_4_4,ind_4_2) + indice_4(i4)%z3*v_4(ilev,i4)*(1. - eps_4)*c(ilev,ind_4_4) |
---|
1107 | mat(ind_4_4,ind_4_4) = mat(ind_4_4,ind_4_4) + indice_4(i4)%z3*v_4(ilev,i4)*eps_4*c(ilev,ind_4_2) |
---|
1108 | mat(ind_4_6,ind_4_2) = mat(ind_4_6,ind_4_2) - indice_4(i4)%z5*v_4(ilev,i4)*(1. - eps_4)*c(ilev,ind_4_4) |
---|
1109 | mat(ind_4_6,ind_4_4) = mat(ind_4_6,ind_4_4) - indice_4(i4)%z5*v_4(ilev,i4)*eps_4*c(ilev,ind_4_2) |
---|
1110 | mat(ind_4_8,ind_4_2) = mat(ind_4_8,ind_4_2) - indice_4(i4)%z7*v_4(ilev,i4)*(1. - eps_4)*c(ilev,ind_4_4) |
---|
1111 | mat(ind_4_8,ind_4_4) = mat(ind_4_8,ind_4_4) - indice_4(i4)%z7*v_4(ilev,i4)*eps_4*c(ilev,ind_4_2) |
---|
1112 | |
---|
1113 | loss(ind_4_2) = loss(ind_4_2) + indice_4(i4)%z1*v_4(ilev,i4)*c(ilev,ind_4_4) |
---|
1114 | loss(ind_4_4) = loss(ind_4_4) + indice_4(i4)%z3*v_4(ilev,i4)*c(ilev,ind_4_2) |
---|
1115 | prod(ind_4_6) = prod(ind_4_6) + indice_4(i4)%z5*v_4(ilev,i4)*c(ilev,ind_4_2)*c(ilev,ind_4_4) |
---|
1116 | prod(ind_4_8) = prod(ind_4_8) + indice_4(i4)%z7*v_4(ilev,i4)*c(ilev,ind_4_2)*c(ilev,ind_4_4) |
---|
1117 | |
---|
1118 | end do |
---|
1119 | |
---|
1120 | end subroutine fill_matrix |
---|
1121 | |
---|
1122 | !================================================================ |
---|
1123 | |
---|
1124 | subroutine indice(i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1125 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1126 | i_n, i_n2d, i_no, i_no2, i_n2) |
---|
1127 | |
---|
1128 | !================================================================ |
---|
1129 | ! set the "indice" arrays used to fill the jacobian matrix ! |
---|
1130 | !---------------------------------------------------------------- |
---|
1131 | ! reaction type ! |
---|
1132 | !---------------------------------------------------------------- |
---|
1133 | ! A + hv --> B + C photolysis indice_phot ! |
---|
1134 | ! A + B --> C + D bimolecular indice_4 ! |
---|
1135 | ! A + A --> B + C quadratic indice_3 ! |
---|
1136 | ! A + C --> B + C quenching indice_phot ! |
---|
1137 | ! A + ice --> B + C heterogeneous indice_phot ! |
---|
1138 | !================================================================ |
---|
1139 | |
---|
1140 | use types_asis |
---|
1141 | |
---|
1142 | implicit none |
---|
1143 | |
---|
1144 | #include "chimiedata.h" |
---|
1145 | |
---|
1146 | ! input |
---|
1147 | |
---|
1148 | integer :: i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1149 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1150 | i_n, i_n2d, i_no, i_no2, i_n2 |
---|
1151 | |
---|
1152 | ! local |
---|
1153 | |
---|
1154 | integer :: nb_phot, nb_reaction_3, nb_reaction_4 |
---|
1155 | integer :: i_dummy |
---|
1156 | |
---|
1157 | allocate (indice_phot(nb_phot_max)) |
---|
1158 | allocate (indice_3(nb_reaction_3_max)) |
---|
1159 | allocate (indice_4(nb_reaction_4_max)) |
---|
1160 | |
---|
1161 | i_dummy = 1 |
---|
1162 | |
---|
1163 | nb_phot = 0 |
---|
1164 | nb_reaction_3 = 0 |
---|
1165 | nb_reaction_4 = 0 |
---|
1166 | |
---|
1167 | !=========================================================== |
---|
1168 | ! O2 + hv -> O + O |
---|
1169 | !=========================================================== |
---|
1170 | |
---|
1171 | nb_phot = nb_phot + 1 |
---|
1172 | |
---|
1173 | indice_phot(nb_phot) = z3spec(1.0, i_o2, 2.0, i_o, 0.0, i_dummy) |
---|
1174 | |
---|
1175 | !=========================================================== |
---|
1176 | ! O2 + hv -> O + O(1D) |
---|
1177 | !=========================================================== |
---|
1178 | |
---|
1179 | nb_phot = nb_phot + 1 |
---|
1180 | |
---|
1181 | indice_phot(nb_phot) = z3spec(1.0, i_o2, 1.0, i_o, 1.0, i_o1d) |
---|
1182 | |
---|
1183 | !=========================================================== |
---|
1184 | ! CO2 + hv -> CO + O |
---|
1185 | !=========================================================== |
---|
1186 | |
---|
1187 | nb_phot = nb_phot + 1 |
---|
1188 | |
---|
1189 | indice_phot(nb_phot) = z3spec(1.0, i_co2, 1.0, i_co, 1.0, i_o) |
---|
1190 | |
---|
1191 | !=========================================================== |
---|
1192 | ! CO2 + hv -> CO + O(1D) |
---|
1193 | !=========================================================== |
---|
1194 | |
---|
1195 | nb_phot = nb_phot + 1 |
---|
1196 | |
---|
1197 | indice_phot(nb_phot) = z3spec(1.0, i_co2, 1.0, i_co, 1.0, i_o1d) |
---|
1198 | |
---|
1199 | !=========================================================== |
---|
1200 | ! O3 + hv -> O2 + O(1D) |
---|
1201 | !=========================================================== |
---|
1202 | |
---|
1203 | nb_phot = nb_phot + 1 |
---|
1204 | |
---|
1205 | indice_phot(nb_phot) = z3spec(1.0, i_o3, 1.0, i_o2, 1.0, i_o1d) |
---|
1206 | |
---|
1207 | !=========================================================== |
---|
1208 | ! O3 + hv -> O2 + O |
---|
1209 | !=========================================================== |
---|
1210 | |
---|
1211 | nb_phot = nb_phot + 1 |
---|
1212 | |
---|
1213 | indice_phot(nb_phot) = z3spec(1.0, i_o3, 1.0, i_o2, 1.0, i_o) |
---|
1214 | |
---|
1215 | !=========================================================== |
---|
1216 | ! H2O + hv -> H + OH |
---|
1217 | !=========================================================== |
---|
1218 | |
---|
1219 | nb_phot = nb_phot + 1 |
---|
1220 | |
---|
1221 | indice_phot(nb_phot) = z3spec(1.0, i_h2o, 1.0, i_h, 1.0, i_oh) |
---|
1222 | |
---|
1223 | !=========================================================== |
---|
1224 | ! H2O2 + hv -> OH + OH |
---|
1225 | !=========================================================== |
---|
1226 | |
---|
1227 | nb_phot = nb_phot + 1 |
---|
1228 | |
---|
1229 | indice_phot(nb_phot) = z3spec(1.0, i_h2o2, 2.0, i_oh, 0.0, i_dummy) |
---|
1230 | |
---|
1231 | !=========================================================== |
---|
1232 | ! HO2 + hv -> OH + O |
---|
1233 | !=========================================================== |
---|
1234 | |
---|
1235 | nb_phot = nb_phot + 1 |
---|
1236 | |
---|
1237 | indice_phot(nb_phot) = z3spec(1.0, i_ho2, 1.0, i_oh, 1.0, i_o) |
---|
1238 | |
---|
1239 | !=========================================================== |
---|
1240 | ! NO + hv -> N + O |
---|
1241 | !=========================================================== |
---|
1242 | |
---|
1243 | nb_phot = nb_phot + 1 |
---|
1244 | |
---|
1245 | indice_phot(nb_phot) = z3spec(1.0, i_no, 1.0, i_n, 1.0, i_o) |
---|
1246 | |
---|
1247 | !=========================================================== |
---|
1248 | ! NO2 + hv -> NO + O |
---|
1249 | !=========================================================== |
---|
1250 | |
---|
1251 | nb_phot = nb_phot + 1 |
---|
1252 | |
---|
1253 | indice_phot(nb_phot) = z3spec(1.0, i_no2, 1.0, i_no, 1.0, i_o) |
---|
1254 | |
---|
1255 | !=========================================================== |
---|
1256 | ! a001 : O + O2 + CO2 -> O3 + CO2 |
---|
1257 | !=========================================================== |
---|
1258 | |
---|
1259 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1260 | |
---|
1261 | indice_4(nb_reaction_4) = z4spec(1.0, i_o, 1.0, i_o2, 1.0, i_o3, 0.0, i_dummy) |
---|
1262 | |
---|
1263 | !=========================================================== |
---|
1264 | ! a002 : O + O + CO2 -> O2 + CO2 |
---|
1265 | !=========================================================== |
---|
1266 | |
---|
1267 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1268 | |
---|
1269 | indice_3(nb_reaction_3) = z3spec(2.0, i_o, 1.0, i_o2, 0.0, i_dummy) |
---|
1270 | |
---|
1271 | !=========================================================== |
---|
1272 | ! a003 : O + O3 -> O2 + O2 |
---|
1273 | !=========================================================== |
---|
1274 | |
---|
1275 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1276 | |
---|
1277 | indice_4(nb_reaction_4) = z4spec(1.0, i_o, 1.0, i_o3, 2.0, i_o2, 0.0, i_dummy) |
---|
1278 | |
---|
1279 | !=========================================================== |
---|
1280 | ! b001 : O(1D) + CO2 -> O + CO2 |
---|
1281 | !=========================================================== |
---|
1282 | |
---|
1283 | nb_phot = nb_phot + 1 |
---|
1284 | |
---|
1285 | indice_phot(nb_phot) = z3spec(1.0, i_o1d, 1.0, i_o, 0.0, i_dummy) |
---|
1286 | |
---|
1287 | !=========================================================== |
---|
1288 | ! b002 : O(1D) + H2O -> OH + OH |
---|
1289 | !=========================================================== |
---|
1290 | |
---|
1291 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1292 | |
---|
1293 | indice_4(nb_reaction_4) = z4spec(1.0, i_o1d, 1.0, i_h2o, 2.0, i_oh, 0.0, i_dummy) |
---|
1294 | |
---|
1295 | !=========================================================== |
---|
1296 | ! b003 : O(1D) + H2 -> OH + H |
---|
1297 | !=========================================================== |
---|
1298 | |
---|
1299 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1300 | |
---|
1301 | indice_4(nb_reaction_4) = z4spec(1.0, i_o1d, 1.0, i_h2, 1.0, i_oh, 1.0, i_h) |
---|
1302 | |
---|
1303 | !=========================================================== |
---|
1304 | ! b004 : O(1D) + O2 -> O + O2 |
---|
1305 | !=========================================================== |
---|
1306 | |
---|
1307 | nb_phot = nb_phot + 1 |
---|
1308 | |
---|
1309 | indice_phot(nb_phot) = z3spec(1.0, i_o1d, 1.0, i_o, 0.0, i_dummy) |
---|
1310 | |
---|
1311 | !=========================================================== |
---|
1312 | ! b005 : O(1D) + O3 -> O2 + O2 |
---|
1313 | !=========================================================== |
---|
1314 | |
---|
1315 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1316 | |
---|
1317 | indice_4(nb_reaction_4) = z4spec(1.0, i_o1d, 1.0, i_o3, 2.0, i_o2, 0.0, i_dummy) |
---|
1318 | |
---|
1319 | !=========================================================== |
---|
1320 | ! b006 : O(1D) + O3 -> O2 + O + O |
---|
1321 | !=========================================================== |
---|
1322 | |
---|
1323 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1324 | |
---|
1325 | indice_4(nb_reaction_4) = z4spec(1.0, i_o1d, 1.0, i_o3, 1.0, i_o2, 2.0, i_o) |
---|
1326 | |
---|
1327 | !=========================================================== |
---|
1328 | ! c001 : O + HO2 -> OH + O2 |
---|
1329 | !=========================================================== |
---|
1330 | |
---|
1331 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1332 | |
---|
1333 | indice_4(nb_reaction_4) = z4spec(1.0, i_o, 1.0, i_ho2, 1.0, i_oh, 1.0, i_o2) |
---|
1334 | |
---|
1335 | !=========================================================== |
---|
1336 | ! c002 : O + OH -> O2 + H |
---|
1337 | !=========================================================== |
---|
1338 | |
---|
1339 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1340 | |
---|
1341 | indice_4(nb_reaction_4) = z4spec(1.0, i_o, 1.0, i_oh, 1.0, i_o2, 1.0, i_h) |
---|
1342 | |
---|
1343 | !=========================================================== |
---|
1344 | ! c003 : H + O3 -> OH + O2 |
---|
1345 | !=========================================================== |
---|
1346 | |
---|
1347 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1348 | |
---|
1349 | indice_4(nb_reaction_4) = z4spec(1.0, i_h, 1.0, i_o3, 1.0, i_oh, 1.0, i_o2) |
---|
1350 | |
---|
1351 | !=========================================================== |
---|
1352 | ! c004 : H + HO2 -> OH + OH |
---|
1353 | !=========================================================== |
---|
1354 | |
---|
1355 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1356 | |
---|
1357 | indice_4(nb_reaction_4) = z4spec(1.0, i_h, 1.0, i_ho2, 2.0, i_oh, 0.0, i_dummy) |
---|
1358 | |
---|
1359 | !=========================================================== |
---|
1360 | ! c005 : H + HO2 -> H2 + O2 |
---|
1361 | !=========================================================== |
---|
1362 | |
---|
1363 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1364 | |
---|
1365 | indice_4(nb_reaction_4) = z4spec(1.0, i_h, 1.0, i_ho2, 1.0, i_h2, 1.0, i_o2) |
---|
1366 | |
---|
1367 | !=========================================================== |
---|
1368 | ! c006 : H + HO2 -> H2O + O |
---|
1369 | !=========================================================== |
---|
1370 | |
---|
1371 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1372 | |
---|
1373 | indice_4(nb_reaction_4) = z4spec(1.0, i_h, 1.0, i_ho2, 1.0, i_h2o, 1.0, i_o) |
---|
1374 | |
---|
1375 | !=========================================================== |
---|
1376 | ! c007 : OH + HO2 -> H2O + O2 |
---|
1377 | !=========================================================== |
---|
1378 | |
---|
1379 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1380 | |
---|
1381 | indice_4(nb_reaction_4) = z4spec(1.0, i_oh, 1.0, i_ho2, 1.0, i_h2o, 1.0, i_o2) |
---|
1382 | |
---|
1383 | !=========================================================== |
---|
1384 | ! c008 : HO2 + HO2 -> H2O2 + O2 |
---|
1385 | !=========================================================== |
---|
1386 | |
---|
1387 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1388 | |
---|
1389 | indice_3(nb_reaction_3) = z3spec(2.0, i_ho2, 1.0, i_h2o2, 1.0, i_o2) |
---|
1390 | |
---|
1391 | !=========================================================== |
---|
1392 | ! c009 : OH + H2O2 -> H2O + HO2 |
---|
1393 | !=========================================================== |
---|
1394 | |
---|
1395 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1396 | |
---|
1397 | indice_4(nb_reaction_4) = z4spec(1.0, i_oh, 1.0, i_h2o2, 1.0, i_h2o, 1.0, i_ho2) |
---|
1398 | |
---|
1399 | !=========================================================== |
---|
1400 | ! c010 : OH + H2 -> H2O + H |
---|
1401 | !=========================================================== |
---|
1402 | |
---|
1403 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1404 | |
---|
1405 | indice_4(nb_reaction_4) = z4spec(1.0, i_oh, 1.0, i_h2, 1.0, i_h2o, 1.0, i_h) |
---|
1406 | |
---|
1407 | !=========================================================== |
---|
1408 | ! c011 : H + O2 + CO2 -> HO2 + CO2 |
---|
1409 | !=========================================================== |
---|
1410 | |
---|
1411 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1412 | |
---|
1413 | indice_4(nb_reaction_4) = z4spec(1.0, i_h, 1.0, i_o2, 1.0, i_ho2, 0.0, i_dummy) |
---|
1414 | |
---|
1415 | !=========================================================== |
---|
1416 | ! c012 : O + H2O2 -> OH + HO2 |
---|
1417 | !=========================================================== |
---|
1418 | |
---|
1419 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1420 | |
---|
1421 | indice_4(nb_reaction_4) = z4spec(1.0, i_o, 1.0, i_h2o2, 1.0, i_oh, 1.0, i_ho2) |
---|
1422 | |
---|
1423 | !=========================================================== |
---|
1424 | ! c013 : OH + OH -> H2O + O |
---|
1425 | !=========================================================== |
---|
1426 | |
---|
1427 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1428 | |
---|
1429 | indice_3(nb_reaction_3) = z3spec(2.0, i_oh, 1.0, i_h2o, 1.0, i_o) |
---|
1430 | |
---|
1431 | !=========================================================== |
---|
1432 | ! c014 : OH + O3 -> HO2 + O2 |
---|
1433 | !=========================================================== |
---|
1434 | |
---|
1435 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1436 | |
---|
1437 | indice_4(nb_reaction_4) = z4spec(1.0, i_oh, 1.0, i_o3, 1.0, i_ho2, 1.0, i_o2) |
---|
1438 | |
---|
1439 | !=========================================================== |
---|
1440 | ! c015 : HO2 + O3 -> OH + O2 + O2 |
---|
1441 | !=========================================================== |
---|
1442 | |
---|
1443 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1444 | |
---|
1445 | indice_4(nb_reaction_4) = z4spec(1.0, i_ho2, 1.0, i_o3, 1.0, i_oh, 2.0, i_o2) |
---|
1446 | |
---|
1447 | !=========================================================== |
---|
1448 | ! c016 : HO2 + HO2 + CO2 -> H2O2 + O2 + CO2 |
---|
1449 | !=========================================================== |
---|
1450 | |
---|
1451 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1452 | |
---|
1453 | indice_3(nb_reaction_3) = z3spec(2.0, i_ho2, 1.0, i_h2o2, 1.0, i_o2) |
---|
1454 | |
---|
1455 | !=========================================================== |
---|
1456 | ! c017 : OH + OH + CO2 -> H2O2 + CO2 |
---|
1457 | !=========================================================== |
---|
1458 | |
---|
1459 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1460 | |
---|
1461 | indice_3(nb_reaction_3) = z3spec(2.0, i_oh, 1.0, i_h2o2, 0.0, i_dummy) |
---|
1462 | |
---|
1463 | !=========================================================== |
---|
1464 | ! c018 : H + H + CO2 -> H2 + CO2 |
---|
1465 | !=========================================================== |
---|
1466 | |
---|
1467 | nb_reaction_3 = nb_reaction_3 + 1 |
---|
1468 | |
---|
1469 | indice_3(nb_reaction_3) = z3spec(2.0, i_h, 1.0, i_h2, 0.0, i_dummy) |
---|
1470 | |
---|
1471 | !=========================================================== |
---|
1472 | ! d001 : NO2 + O -> NO + O2 |
---|
1473 | !=========================================================== |
---|
1474 | |
---|
1475 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1476 | |
---|
1477 | indice_4(nb_reaction_4) = z4spec(1.0, i_no2, 1.0, i_o, 1.0, i_no, 1.0, i_o2) |
---|
1478 | |
---|
1479 | !=========================================================== |
---|
1480 | ! d002 : NO + O3 -> NO2 + O2 |
---|
1481 | !=========================================================== |
---|
1482 | |
---|
1483 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1484 | |
---|
1485 | indice_4(nb_reaction_4) = z4spec(1.0, i_no, 1.0, i_o3, 1.0, i_no2, 1.0, i_o2) |
---|
1486 | |
---|
1487 | !=========================================================== |
---|
1488 | ! d003 : NO + HO2 -> NO2 + OH |
---|
1489 | !=========================================================== |
---|
1490 | |
---|
1491 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1492 | |
---|
1493 | indice_4(nb_reaction_4) = z4spec(1.0, i_no, 1.0, i_ho2, 1.0, i_no2, 1.0, i_oh) |
---|
1494 | |
---|
1495 | !=========================================================== |
---|
1496 | ! d004 : N + NO -> N2 + O |
---|
1497 | !=========================================================== |
---|
1498 | |
---|
1499 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1500 | |
---|
1501 | indice_4(nb_reaction_4) = z4spec(1.0, i_n, 1.0, i_no, 1.0, i_n2, 1.0, i_o) |
---|
1502 | |
---|
1503 | !=========================================================== |
---|
1504 | ! d005 : N + O2 -> NO + O |
---|
1505 | !=========================================================== |
---|
1506 | |
---|
1507 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1508 | |
---|
1509 | indice_4(nb_reaction_4) = z4spec(1.0, i_n, 1.0, i_o2, 1.0, i_no, 1.0, i_o) |
---|
1510 | |
---|
1511 | !=========================================================== |
---|
1512 | ! d006 : NO2 + H -> NO + OH |
---|
1513 | !=========================================================== |
---|
1514 | |
---|
1515 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1516 | |
---|
1517 | indice_4(nb_reaction_4) = z4spec(1.0, i_no2, 1.0, i_h, 1.0, i_no, 1.0, i_oh) |
---|
1518 | |
---|
1519 | !=========================================================== |
---|
1520 | ! d007 : N + O -> NO |
---|
1521 | !=========================================================== |
---|
1522 | |
---|
1523 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1524 | |
---|
1525 | indice_4(nb_reaction_4) = z4spec(1.0, i_n, 1.0, i_o, 1.0, i_no, 0.0, i_dummy) |
---|
1526 | |
---|
1527 | !=========================================================== |
---|
1528 | ! d008 : N + HO2 -> NO + OH |
---|
1529 | !=========================================================== |
---|
1530 | |
---|
1531 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1532 | |
---|
1533 | indice_4(nb_reaction_4) = z4spec(1.0, i_n, 1.0, i_ho2, 1.0, i_no, 1.0, i_oh) |
---|
1534 | |
---|
1535 | !=========================================================== |
---|
1536 | ! d009 : N + OH -> NO + H |
---|
1537 | !=========================================================== |
---|
1538 | |
---|
1539 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1540 | |
---|
1541 | indice_4(nb_reaction_4) = z4spec(1.0, i_n, 1.0, i_oh, 1.0, i_no, 1.0, i_h) |
---|
1542 | |
---|
1543 | !=========================================================== |
---|
1544 | ! e001 : CO + OH -> CO2 + H |
---|
1545 | !=========================================================== |
---|
1546 | |
---|
1547 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1548 | |
---|
1549 | indice_4(nb_reaction_4) = z4spec(1.0, i_co, 1.0, i_oh, 1.0, i_co2, 1.0, i_h) |
---|
1550 | |
---|
1551 | !=========================================================== |
---|
1552 | ! e002 : CO + O + M -> CO2 + M |
---|
1553 | !=========================================================== |
---|
1554 | |
---|
1555 | nb_reaction_4 = nb_reaction_4 + 1 |
---|
1556 | |
---|
1557 | indice_4(nb_reaction_4) = z4spec(1.0, i_co, 1.0, i_o, 1.0, i_co2, 0.0, i_dummy) |
---|
1558 | |
---|
1559 | !=========================================================== |
---|
1560 | ! h001: HO2 + ice -> products |
---|
1561 | ! treated as |
---|
1562 | ! HO2 -> 0.5 H2O + 0.75 O2 |
---|
1563 | !=========================================================== |
---|
1564 | |
---|
1565 | nb_phot = nb_phot + 1 |
---|
1566 | |
---|
1567 | indice_phot(nb_phot) = z3spec(1.0, i_ho2, 0.5, i_h2o, 0.75, i_o2) |
---|
1568 | |
---|
1569 | !=========================================================== |
---|
1570 | ! h002: OH + ice -> products |
---|
1571 | ! treated as |
---|
1572 | ! OH -> 0.5 H2O + 0.25 O2 |
---|
1573 | !=========================================================== |
---|
1574 | |
---|
1575 | nb_phot = nb_phot + 1 |
---|
1576 | |
---|
1577 | indice_phot(nb_phot) = z3spec(1.0, i_oh, 0.5, i_h2o, 0.25, i_o2) |
---|
1578 | |
---|
1579 | !=========================================================== |
---|
1580 | ! h003: H2O2 + ice -> products |
---|
1581 | ! treated as |
---|
1582 | ! H2O2 -> H2O + 0.5 O2 |
---|
1583 | !=========================================================== |
---|
1584 | |
---|
1585 | nb_phot = nb_phot + 1 |
---|
1586 | |
---|
1587 | indice_phot(nb_phot) = z3spec(1.0, i_h2o2, 1.0, i_h2o, 0.5, i_o2) |
---|
1588 | |
---|
1589 | !=========================================================== |
---|
1590 | ! h004: HO2 + dust -> products |
---|
1591 | ! treated as |
---|
1592 | ! HO2 -> 0.5 H2O + 0.75 O2 |
---|
1593 | !=========================================================== |
---|
1594 | |
---|
1595 | nb_phot = nb_phot + 1 |
---|
1596 | |
---|
1597 | indice_phot(nb_phot) = z3spec(1.0, i_ho2, 0.5, i_h2o, 0.75, i_o2) |
---|
1598 | |
---|
1599 | !=========================================================== |
---|
1600 | ! h005: H2O2 + dust -> products |
---|
1601 | ! treated as |
---|
1602 | ! H2O2 -> H2O + 0.5 O2 |
---|
1603 | !=========================================================== |
---|
1604 | |
---|
1605 | nb_phot = nb_phot + 1 |
---|
1606 | |
---|
1607 | indice_phot(nb_phot) = z3spec(1.0, i_h2o2, 1.0, i_h2o, 0.5, i_o2) |
---|
1608 | |
---|
1609 | !=========================================================== |
---|
1610 | ! check dimensions |
---|
1611 | !=========================================================== |
---|
1612 | |
---|
1613 | print*, 'nb_phot = ', nb_phot |
---|
1614 | print*, 'nb_reaction_4 = ', nb_reaction_4 |
---|
1615 | print*, 'nb_reaction_3 = ', nb_reaction_3 |
---|
1616 | |
---|
1617 | if ((nb_phot /= nb_phot_max) .or. & |
---|
1618 | (nb_reaction_3 /= nb_reaction_3_max) .or. & |
---|
1619 | (nb_reaction_4 /= nb_reaction_4_max)) then |
---|
1620 | print*, 'wrong dimensions in indice' |
---|
1621 | stop |
---|
1622 | end if |
---|
1623 | |
---|
1624 | end subroutine indice |
---|
1625 | |
---|
1626 | !***************************************************************** |
---|
1627 | |
---|
1628 | subroutine gcmtochim(nlayer, nq, zycol, lswitch, nesp, & |
---|
1629 | i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1630 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1631 | i_n, i_n2d, i_no, i_no2, i_n2, & |
---|
1632 | dens, rm, c) |
---|
1633 | |
---|
1634 | !***************************************************************** |
---|
1635 | |
---|
1636 | use tracer_h, only: igcm_co2, igcm_co, igcm_o, igcm_o1d, & |
---|
1637 | & igcm_o2, igcm_o3, igcm_h, igcm_h2, igcm_oh, & |
---|
1638 | & igcm_ho2, igcm_h2o2, igcm_h2o_vap, & |
---|
1639 | & igcm_n, igcm_n2d, igcm_no, igcm_no2, igcm_n2 |
---|
1640 | |
---|
1641 | use callkeys_mod |
---|
1642 | |
---|
1643 | implicit none |
---|
1644 | |
---|
1645 | |
---|
1646 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1647 | ! input: |
---|
1648 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1649 | |
---|
1650 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
1651 | integer, intent(in) :: nq ! number of tracers in the gcm |
---|
1652 | integer :: nesp ! number of species in the chemistry |
---|
1653 | integer :: lswitch ! interface level between chemistries |
---|
1654 | |
---|
1655 | integer :: i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1656 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1657 | i_n, i_n2d, i_no, i_no2, i_n2 |
---|
1658 | |
---|
1659 | real :: zycol(nlayer,nq) ! volume mixing ratios in the gcm |
---|
1660 | real :: dens(nlayer) ! total number density (molecule.cm-3) |
---|
1661 | |
---|
1662 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1663 | ! output: |
---|
1664 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1665 | |
---|
1666 | real, dimension(nlayer,nesp) :: rm ! volume mixing ratios |
---|
1667 | real (kind = 8), dimension(nlayer,nesp) :: c ! number densities (molecule.cm-3) |
---|
1668 | |
---|
1669 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1670 | ! local: |
---|
1671 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1672 | |
---|
1673 | integer :: l, iesp |
---|
1674 | logical,save :: firstcall = .true. |
---|
1675 | |
---|
1676 | |
---|
1677 | ! first call initializations |
---|
1678 | |
---|
1679 | if (firstcall) then |
---|
1680 | |
---|
1681 | ! identify the indexes of the tracers we need |
---|
1682 | |
---|
1683 | if (igcm_co2 == 0) then |
---|
1684 | write(*,*) "gcmtochim: Error; no CO2 tracer !!!" |
---|
1685 | stop |
---|
1686 | endif |
---|
1687 | if (igcm_co == 0) then |
---|
1688 | write(*,*) "gcmtochim: Error; no CO tracer !!!" |
---|
1689 | stop |
---|
1690 | end if |
---|
1691 | if (igcm_o == 0) then |
---|
1692 | write(*,*) "gcmtochim: Error; no O tracer !!!" |
---|
1693 | stop |
---|
1694 | end if |
---|
1695 | if (igcm_o1d == 0) then |
---|
1696 | write(*,*) "gcmtochim: Error; no O1D tracer !!!" |
---|
1697 | stop |
---|
1698 | end if |
---|
1699 | if (igcm_o2 == 0) then |
---|
1700 | write(*,*) "gcmtochim: Error; no O2 tracer !!!" |
---|
1701 | stop |
---|
1702 | end if |
---|
1703 | if (igcm_o3 == 0) then |
---|
1704 | write(*,*) "gcmtochim: Error; no O3 tracer !!!" |
---|
1705 | stop |
---|
1706 | end if |
---|
1707 | if (igcm_h == 0) then |
---|
1708 | write(*,*) "gcmtochim: Error; no H tracer !!!" |
---|
1709 | stop |
---|
1710 | end if |
---|
1711 | if (igcm_h2 == 0) then |
---|
1712 | write(*,*) "gcmtochim: Error; no H2 tracer !!!" |
---|
1713 | stop |
---|
1714 | end if |
---|
1715 | if (igcm_oh == 0) then |
---|
1716 | write(*,*) "gcmtochim: Error; no OH tracer !!!" |
---|
1717 | stop |
---|
1718 | end if |
---|
1719 | if (igcm_ho2 == 0) then |
---|
1720 | write(*,*) "gcmtochim: Error; no HO2 tracer !!!" |
---|
1721 | stop |
---|
1722 | end if |
---|
1723 | if (igcm_h2o2 == 0) then |
---|
1724 | write(*,*) "gcmtochim: Error; no H2O2 tracer !!!" |
---|
1725 | stop |
---|
1726 | end if |
---|
1727 | if (igcm_n == 0) then |
---|
1728 | write(*,*) "gcmtochim: Error; no N tracer !!!" |
---|
1729 | ! stop |
---|
1730 | end if |
---|
1731 | if (igcm_n2d == 0) then |
---|
1732 | write(*,*) "gcmtochim: Error; no N2D tracer !!!" |
---|
1733 | ! stop |
---|
1734 | end if |
---|
1735 | if (igcm_no == 0) then |
---|
1736 | write(*,*) "gcmtochim: Error; no NO tracer !!!" |
---|
1737 | ! stop |
---|
1738 | end if |
---|
1739 | if (igcm_no2 == 0) then |
---|
1740 | write(*,*) "gcmtochim: Error; no NO2 tracer !!!" |
---|
1741 | ! stop |
---|
1742 | end if |
---|
1743 | if (igcm_n2 == 0) then |
---|
1744 | write(*,*) "gcmtochim: Error; no N2 tracer !!!" |
---|
1745 | stop |
---|
1746 | end if |
---|
1747 | if (igcm_h2o_vap == 0) then |
---|
1748 | write(*,*) "gcmtochim: Error; no water vapor tracer !!!" |
---|
1749 | stop |
---|
1750 | end if |
---|
1751 | firstcall = .false. |
---|
1752 | end if ! of if (firstcall) |
---|
1753 | |
---|
1754 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1755 | ! initialise mixing ratios |
---|
1756 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1757 | |
---|
1758 | do l = 1,lswitch-1 |
---|
1759 | rm(l,i_co2) = zycol(l, igcm_co2) |
---|
1760 | rm(l,i_co) = zycol(l, igcm_co) |
---|
1761 | rm(l,i_o) = zycol(l, igcm_o) |
---|
1762 | rm(l,i_o1d) = zycol(l, igcm_o1d) |
---|
1763 | rm(l,i_o2) = zycol(l, igcm_o2) |
---|
1764 | rm(l,i_o3) = zycol(l, igcm_o3) |
---|
1765 | rm(l,i_h) = zycol(l, igcm_h) |
---|
1766 | rm(l,i_h2) = zycol(l, igcm_h2) |
---|
1767 | rm(l,i_oh) = zycol(l, igcm_oh) |
---|
1768 | rm(l,i_ho2) = zycol(l, igcm_ho2) |
---|
1769 | rm(l,i_h2o2) = zycol(l, igcm_h2o2) |
---|
1770 | rm(l,i_h2o) = zycol(l, igcm_h2o_vap) |
---|
1771 | rm(l,i_n) = zycol(l, igcm_n) |
---|
1772 | rm(l,i_n2d) = zycol(l, igcm_n2d) |
---|
1773 | rm(l,i_no) = zycol(l, igcm_no) |
---|
1774 | rm(l,i_no2) = zycol(l, igcm_no2) |
---|
1775 | rm(l,i_n2) = zycol(l, igcm_n2) |
---|
1776 | end do |
---|
1777 | |
---|
1778 | where (rm(:,:) < 1.e-30) |
---|
1779 | rm(:,:) = 0. |
---|
1780 | end where |
---|
1781 | |
---|
1782 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1783 | ! initialise number densities |
---|
1784 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1785 | |
---|
1786 | do iesp = 1,nesp |
---|
1787 | do l = 1,lswitch-1 |
---|
1788 | c(l,iesp) = rm(l,iesp)*dens(l) |
---|
1789 | end do |
---|
1790 | end do |
---|
1791 | |
---|
1792 | end subroutine gcmtochim |
---|
1793 | |
---|
1794 | !***************************************************************** |
---|
1795 | |
---|
1796 | subroutine chimtogcm(nlayer, nq, zycol, lswitch, nesp, & |
---|
1797 | i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1798 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1799 | i_n, i_n2d, i_no, i_no2, i_n2, dens, c) |
---|
1800 | |
---|
1801 | !***************************************************************** |
---|
1802 | |
---|
1803 | use tracer_h, only: igcm_co2, igcm_co, igcm_o, igcm_o1d, & |
---|
1804 | igcm_o2, igcm_o3, igcm_h, igcm_h2, igcm_oh, & |
---|
1805 | igcm_ho2, igcm_h2o2, igcm_h2o_vap, & |
---|
1806 | igcm_n, igcm_n2d, igcm_no, igcm_no2, igcm_n2 |
---|
1807 | |
---|
1808 | use callkeys_mod |
---|
1809 | |
---|
1810 | implicit none |
---|
1811 | |
---|
1812 | |
---|
1813 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1814 | ! inputs: |
---|
1815 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1816 | |
---|
1817 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
1818 | integer, intent(in) :: nq ! number of tracers in the gcm |
---|
1819 | integer :: nesp ! number of species in the chemistry |
---|
1820 | integer :: lswitch ! interface level between chemistries |
---|
1821 | integer :: i_co2, i_co, i_o, i_o1d, i_o2, i_o3, i_h, & |
---|
1822 | i_h2, i_oh, i_ho2, i_h2o2, i_h2o, & |
---|
1823 | i_n, i_n2d, i_no, i_no2, i_n2 |
---|
1824 | |
---|
1825 | real :: dens(nlayer) ! total number density (molecule.cm-3) |
---|
1826 | real (kind = 8), dimension(nlayer,nesp) :: c ! number densities (molecule.cm-3) |
---|
1827 | |
---|
1828 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1829 | ! output: |
---|
1830 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1831 | |
---|
1832 | real zycol(nlayer,nq) ! volume mixing ratios in the gcm |
---|
1833 | |
---|
1834 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1835 | ! local: |
---|
1836 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1837 | |
---|
1838 | integer l |
---|
1839 | |
---|
1840 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1841 | ! save mixing ratios for the gcm |
---|
1842 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1843 | |
---|
1844 | do l = 1,lswitch-1 |
---|
1845 | zycol(l, igcm_co2) = c(l,i_co2)/dens(l) |
---|
1846 | zycol(l, igcm_co) = c(l,i_co)/dens(l) |
---|
1847 | zycol(l, igcm_o) = c(l,i_o)/dens(l) |
---|
1848 | zycol(l, igcm_o1d) = c(l,i_o1d)/dens(l) |
---|
1849 | zycol(l, igcm_o2) = c(l,i_o2)/dens(l) |
---|
1850 | zycol(l, igcm_o3) = c(l,i_o3)/dens(l) |
---|
1851 | zycol(l, igcm_h) = c(l,i_h)/dens(l) |
---|
1852 | zycol(l, igcm_h2) = c(l,i_h2)/dens(l) |
---|
1853 | zycol(l, igcm_oh) = c(l,i_oh)/dens(l) |
---|
1854 | zycol(l, igcm_ho2) = c(l,i_ho2)/dens(l) |
---|
1855 | zycol(l, igcm_h2o2) = c(l,i_h2o2)/dens(l) |
---|
1856 | zycol(l, igcm_h2o_vap) = c(l,i_h2o)/dens(l) |
---|
1857 | zycol(l, igcm_n) = c(l,i_n)/dens(l) |
---|
1858 | zycol(l, igcm_n2d) = c(l,i_n2d)/dens(l) |
---|
1859 | zycol(l, igcm_no) = c(l,i_no)/dens(l) |
---|
1860 | zycol(l, igcm_no2) = c(l,i_no2)/dens(l) |
---|
1861 | zycol(l, igcm_n2) = c(l,i_n2)/dens(l) |
---|
1862 | end do |
---|
1863 | |
---|
1864 | end subroutine chimtogcm |
---|
1865 | |
---|
1866 | end subroutine photochemistry_asis |
---|
1867 | |
---|