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