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