source: trunk/LMDZ.MARS/libf/aeronomars/paramfoto.F @ 461

Last change on this file since 461 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

File size: 17.1 KB
Line 
1c**********************************************************************
2     
3      subroutine paramfoto
4     $(lswitch,zdens,tx,timestep,zenit,nz,co2x,o2x,o3px,cox,hx,ohx,
5     & ho2x,h2x,h2ox,h2o2x,o1dx)
6
7c     oct 2002      fgg
8c**********************************************************************
9
10      implicit none
11
12      include 'param.h'
13      include 'param_v3.h'
14
15c     arguments
16
17      integer nz,lswitch
18      real zdens(nz)
19      real tx(nz)
20      real co2x(nz),o2x(nz),o3px(nz),cox(nz),hx(nz),ohx(nz)
21      real    h2x(nz),h2ox(nz),h2o2x(nz),o1dx(nz),ho2x(nz)
22      real    zenit
23      real    timestep
24
25c     local variables
26
27      real*8    dt,deltat
28      integer nt, ntmax
29      real*8    tauco2(nzmax,nreact)
30      real*8    tauo2(nzmax,nreact)
31      real*8    tauo3p(nzmax,nreact)
32      real*8    tauco(nzmax,nreact)
33      real*8    tauh(nzmax,nreact)
34      real*8    tauoh(nzmax,nreact)
35      real*8    tauho2(nzmax,nreact)
36      real*8    tauh2(nzmax,nreact)
37      real*8    tauh2o(nzmax,nreact)
38      real*8    tauo1d(nzmax,nreact)
39      real*8    tauh2o2(nzmax,nreact)
40      real*8    co2x8,o2x8,o3px8,cox8,hx8,ohx8,ho2x8,h2x8,h2ox8
41      real*8    h2o2x8,o1dx8
42      real*8    o1dx8aux,ohx8aux,ho2x8aux,dif_o1d,dif_oh,dif_ho2
43      real*8    timestep8
44      real*8    jdistot8(nabs,nzmax)
45      real*8    jdistot8_b(2,nzmax)
46      real*8    tx8
47
48      real*8    tmin(nzmax)
49      integer i,j,k,nada,n
50
51      real*8    auxp,auxl
52     
53      character*1 o1d_eq(nzmax)
54      character*1 oh_eq(nzmax)
55      character*1 ho2_eq(nzmax)
56      character*1 dn
57
58      logical firstcall
59      save firstcall
60      data firstcall /.true./
61
62      if (firstcall) then
63        firstcall= .false.
64        print*,'FOTO'
65      endif
66
67c 123  format(1x,f10.7,11(tr3,e12.7))
68
69c**********************************************************************
70
71C     Start: altitude loop
72      timestep8=dble(timestep)
73
74      do i=nz,lswitch,-1
75         
76
77c     Concentrations to real*8
78
79         co2x8=dble(co2x(i))
80         o2x8=dble(o2x(i))
81         o3px8=dble(o3px(i))
82         cox8=dble(cox(i))
83         hx8=dble(hx(i))
84         ohx8=dble(ohx(i))
85         ho2x8=dble(ho2x(i))
86         h2x8=dble(h2x(i))
87         h2ox8=dble(h2ox(i))
88         h2o2x8=dble(h2o2x(i))
89         o1dx8=dble(o1dx(i))
90         tx8=dble(tx(i))
91
92
93c     reaction rates         
94
95         call getch(tx8)
96
97c         print*,i,co2x(i),o2x(i),o3px(i),cox(i),hx(i),ohx(i),ho2x(i),
98c     &            h2x(i),h2ox(i),h2o2x(i),o1dx(i),tx(i),zdens(i)
99
100
101c     photodissotiation rates
102
103         jdistot(1,i) = 0.
104         jdistot(2,i) = 0.
105         jdistot(3,i) = 0.
106         jdistot(4,i) = 0.
107         jdistot(5,i) = 0.
108         jdistot(6,i) = 0.
109         jdistot_b(1,i) = 0.
110         jdistot_b(2,i) = 0.
111
112         call phdisrate(zenit,i)
113
114         do j=1,6
115            jdistot8(j,i) = dble(jdistot(j,i))
116         end do
117         jdistot8_b(1,i) = dble(jdistot_b(1,i))
118         jdistot8_b(2,i) = dble(jdistot_b(2,i))
119         
120         
121c     Lifetimes calculation
122
123         do j=1,nreact
124            tauco2(i,j) = 1.d30
125            tauo2(i,j) = 1.d30
126            tauo3p(i,j) = 1.d30
127            tauco(i,j) = 1.d30
128            tauh(i,j) = 1.d30
129            tauoh(i,j) = 1.d30
130            tauho2(i,j) = 1.d30
131            tauh2(i,j) = 1.d30
132            tauh2o(i,j) = 1.d30
133            tauo1d(i,j) = 1.d30
134            tauh2o2(i,j) = 1.d30
135         end do
136
137         if(jdistot8(1,i).gt.1d-30) tauco2(i,1) = 1.d0 / jdistot8(1,i)
138
139         if(ch2*o2x8*co2x8.gt.1d-30)
140     $      tauh(i,2) = 1.d0 / (ch2 * o2x8 * co2x8) 
141         if(ch2*hx8*co2x8.gt.1d-30)
142     $      tauo2(i,2) = 1.d0 / (ch2 * hx8 * co2x8)
143 
144         if(ch3*o3px8.gt.1d-30) tauho2(i,3) = 1.d0 / (ch3 * o3px8)
145         if(ch3*ho2x8.gt.1d-30) tauo3p(i,3) = 1.d0 / (ch3 * ho2x8)
146 
147         if(ch4*cox8.gt.1d-30) tauoh(i,4) = 1.d0 / (ch4 * cox8)
148         if(ch4*ohx8.gt.1d-30) tauco(i,4) = 1.d0 / (ch4 * ohx8)
149 
150         if(ch5*ho2x8.gt.1d-30)tauho2(i,5) = 1.d0 / (2.d0 * ch5 * ho2x8)
151
152         if(jdistot8(6,i).gt.1d-30) tauh2o2(i,6) = 1.d0 / jdistot8(6,i)
153
154         if(ch7*ohx8.gt.1d-30) tauho2(i,7) = 1.d0 / (ch7 * ohx8)
155         if(ch7*ho2x8.gt.1d-30) tauoh(i,7) = 1.d0 / (ch7 * ho2x8)
156 
157         if(jdistot8(4,i).gt.1d-30) tauh2o(i,8) = 1.d0 / jdistot8(4,i)
158
159         if(ch9*o1dx8.gt.1d-30) tauh2o(i,9) = 1.d0 / (ch9 * o1dx8)
160         if(ch9*h2ox8.gt.1d-30) tauo1d(i,9) = 1.d0 / (ch9 * h2ox8)
161
162         if(ch10*o3px8*co2x8.gt.1d-30)
163     $     tauo3p(i,10) = 1.d0 / (2.d0 * ch10 * o3px8 * co2x8)
164
165         if(ch11*o3px8.gt.1d-30) tauoh(i,11) = 1.d0 / (ch11 * o3px8)
166         if(ch11*ohx8.gt.1d-30) tauo3p(i,11) = 1.d0 / (ch11 * ohx8)
167
168         if(jdistot8(2,i).gt.1d-30) tauo2(i,12) = 1.d0 / jdistot8(2,i)
169
170         if(ch13*hx8.gt.1d-30) tauho2(i,13) = 1.d0 / (ch13 * hx8)
171         if(ch13*ho2x8.gt.1d-30) tauh(i,13) = 1.d0 / (ch13 * ho2x8)
172
173         if(ch14*o1dx8.gt.1d-30) tauh2(i,14) = 1.d0 / (ch14 * o1dx8)
174         if(ch14*h2x8.gt.1d-30) tauo1d(i,14) = 1.d0 / (ch14 * h2x8)
175
176         if(ch15*ohx8.gt.1d-30) tauh2(i,15) = 1.d0 / (ch15 * ohx8)
177         if(ch15*h2x8.gt.1d-30) tauoh(i,15) = 1.d0 / (ch15 * h2x8)
178
179         if(jdistot8_b(1,i).gt.1d-30) tauco2(i,16) =1.d0/jdistot8_b(1,i)
180
181         if(jdistot8_b(2,i).gt.1d-30) tauo2(i,17) =1.d0/jdistot8_b(2,i)
182
183         if(ch18*ohx8.gt.1d-30) tauh2o2(i,18) = 1.d0 / (ch18 * ohx8)
184         if(ch18*h2o2x8.gt.1d-30) tauoh(i,18) = 1.d0 / (ch18 * h2o2x8)
185
186         if(ch19*co2x8.gt.1d-30) tauo1d(i,19) = 1.d0 / (ch19 * co2x8)
187
188         if(ch20*o2x8.gt.1d-30) tauo1d(i,20) = 1.d0 / (ch20 * o2x8)
189
190         if(ch21*o2x8*co2x8.gt.1d-30) tauo3p(i,21) = 1.d0 /
191     $        (ch21 * o2x8 * co2x8)
192         if(ch21*o3px8*co2x8.gt.1d-30) tauo2(i,21) = 1.d0 /
193     $        (ch21 * o3px8 * co2x8)
194
195         if(jdistot8(5,i).gt.1d-30) tauh2(i,22) = 1.d0 / jdistot8(5,i)
196
197         !minimum lifetime
198         tminco2(i) = 1.d30
199         tmino2(i) = 1.d30
200         tmino3p(i) = 1.d30
201         tminco(i) = 1.d30
202         tminh(i) = 1.d30
203         tminoh(i) = 1.d30
204         tminho2(i) = 1.d30
205         tminh2(i) = 1.d30
206         tminh2o(i) = 1.d30
207         tmino1d(i) = 1.d30
208         tminh2o2(i) = 1.d30
209         
210         do j=1,nreact
211            tminco2(i) = min(tminco2(i),tauco2(i,j))
212            tmino2(i) = min(tmino2(i),tauo2(i,j))
213            tmino3p(i) = min(tmino3p(i),tauo3p(i,j))
214            tminco(i) = min(tminco(i),tauco(i,j))
215            tminh(i) = min(tminh(i),tauh(i,j))
216            tminoh(i) = min(tminoh(i),tauoh(i,j))
217            tminho2(i) = min(tminho2(i),tauho2(i,j))
218            tminh2(i) = min(tminh2(i),tauh2(i,j))
219            tminh2o(i) = min(tminh2o(i),tauh2o(i,j))
220            tmino1d(i) = min(tmino1d(i),tauo1d(i,j))
221            tminh2o2(i) = min(tminh2o2(i),tauh2o2(i,j))
222         end do
223
224         tmin(i) = 1.d30
225         tmin(i) = min(tmin(i),tminco2(i))
226         tmin(i) = min(tmin(i),tmino2(i))
227         tmin(i) = min(tmin(i),tmino3p(i))
228         tmin(i) = min(tmin(i),tminco(i))
229         tmin(i) = min(tmin(i),tminh(i))
230         tmin(i) = min(tmin(i),tminh2(i))
231         tmin(i) = min(tmin(i),tminh2o(i))
232         tmin(i) = min(tmin(i),tminh2o2(i))
233
234
235         !Internal timestep
236         if(tmin(i)/10.d0.gt.timestep8*3600.d0) then
237            deltat = timestep8 * 3600.d0
238         else
239            deltat = tmin(i)/10.d0
240         end if
241         if(deltat.ne.timestep8*3600.d0) then
242            n=int(timestep8*3600.d0/deltat)+1
243            deltat=timestep8*3600.d0/n
244         endif
245
246
247         !Photochemical equilibrium for O1D,OH and HO2?
248         if(tmino1d(i).lt.deltat/5.d0) then
249            o1d_eq(i)='Y'
250         else
251            o1d_eq(i)='N'
252c            o1d_eq(i)='Y'
253         endif
254         if(tminoh(i).lt.deltat/5.d0) then
255           oh_eq(i)='Y'
256         else
257           oh_eq(i)='N'
258c           oh_eq(i)='Y'
259         end if
260         if(tminho2(i).lt.deltat/5.d0) then
261            ho2_eq(i)='Y'
262         else
263            ho2_eq(i)='N'
264c            ho2_eq(i)='Y'
265         endif
266
267
268C     Start: temporal loop
269c        do dt=deltat,timestep8*3600.d0,deltat
270         do nt=1,int(timestep8*3600.d0/deltat)
271            dt=nt*deltat
272
273            !Productions and losses
274            do j=1,nreact
275               Lco2(i,j) = 0.d0
276               Pco2(i,j) = 0.d0
277               Lo2(i,j) = 0.d0
278               Po2(i,j) = 0.d0
279               Lo3p(i,j) = 0.d0
280               Po3p(i,j) = 0.d0
281               Lco(i,j) = 0.d0
282               Pco(i,j) = 0.d0
283               Ph(i,j) = 0.d0
284               Lh(i,j) = 0.d0
285               Poh(i,j) = 0.d0
286               Loh(i,j) = 0.d0
287               Pho2(i,j) = 0.d0
288               Lho2(i,j) = 0.d0
289               Ph2(i,j) = 0.d0
290               Lh2(i,j) = 0.d0
291               Ph2o(i,j) = 0.d0
292               Lh2o(i,j) = 0.d0
293               Po1d(i,j) = 0.d0
294               Lo1d(i,j) = 0.d0
295               Ph2o2(i,j) = 0.d0
296               Lh2o2(i,j) = 0.d0
297            end do
298            Pco2tot(i) = 0.d0
299            Lco2tot(i) = 0.d0
300            Po2tot(i) = 0.d0
301            Lo2tot(i) = 0.d0
302            Po3ptot(i) = 0.d0
303            Lo3ptot(i) = 0.d0
304            Pcotot(i) = 0.d0
305            Lcotot(i) = 0.d0
306            Phtot(i) = 0.d0
307            Lhtot(i) = 0.d0
308            Pohtot(i) = 0.d0
309            Lohtot(i) = 0.d0
310            Pho2tot(i) = 0.d0
311            Lho2tot(i) = 0.d0
312            Ph2tot(i) = 0.d0
313            Lh2tot(i) = 0.d0
314            Ph2otot(i) = 0.d0
315            Lh2otot(i) = 0.d0
316            Po1dtot(i) = 0.d0
317            Lo1dtot(i) = 0.d0
318            Ph2o2tot(i) = 0.d0
319            Lh2o2tot(i) = 0.d0
320
321
322c           R1: CO2 + hv -> CO + O
323
324            Lco2(i,1) = jdistot8(1,i)
325            Pco(i,1) = co2x8 * jdistot8(1,i)
326            Po3p(i,1) = co2x8 * jdistot8(1,i)
327
328
329c           R16(1b): CO2 + hv -> CO + O1D
330
331            Lco2(i,16) = jdistot8_b(1,i)
332            Pco(i,16) = co2x8 * jdistot8_b(1,i)
333            Po1d(i,16) = co2x8 * jdistot8_b(1,i)
334
335
336c           R2: H + O2 + CO2 -> HO2 + CO2
337
338            Lh(i,2) = ch2 * o2x8 * co2x8
339            Lo2(i,2) = ch2 * hx8 * co2x8
340            Pho2(i,2) = ch2 * hx8 * o2x8 * co2x8
341
342
343c           R3: O + HO2 -> OH + O2
344
345            Lo3p(i,3) = ch3 * ho2x8
346            Lho2(i,3) = ch3 * o3px8
347            Poh(i,3) = ch3 * o3px8 * ho2x8
348            Po2(i,3) = ch3 * o3px8 * ho2x8
349
350
351c           R4: CO + OH -> CO2 + H
352
353            Lco(i,4) = ch4 * ohx8
354            Loh(i,4) = ch4 * cox8
355            Pco2(i,4) = ch4 * cox8 * ohx8
356            Ph(i,4) = ch4 * cox8 * ohx8
357
358
359c           R5: 2HO2 -> H2O2 + O2
360
361            Lho2(i,5) = 2.d0 * ch5 * ho2x8
362            Po2(i,5) = ch5 * ho2x8 * ho2x8
363            Ph2o2(i,5) = ch5 * ho2x8 * ho2x8
364
365
366c           R6: H2O2 + hv -> 2OH
367
368            Lh2o2(i,6) = jdistot8(6,i)
369            Poh(i,6) = jdistot8(6,i) * h2o2x8
370
371
372c           R7: OH + HO2 -> H2O + O2
373
374            Loh(i,7) = ch7 * ho2x8
375            Lho2(i,7) = ch7 * ohx8
376            Po2(i,7) = ch7 * ohx8 * ho2x8
377            Ph2o(i,7) = ch7 * ohx8 * ho2x8
378
379
380c           R8: H20 + hv -> H + OH
381
382            Lh2o(i,8) = jdistot8(4,i)
383            Ph(i,8) = jdistot8(4,i) * h2ox8
384            Poh(i,8) = jdistot8(4,i) * h2ox8
385
386
387c           R9: O(1D) + H2O -> 2OH
388
389            Lo1d(i,9) = ch9 * h2ox8
390            Lh2o(i,9) = ch9 * o1dx8
391            Poh(i,9) = 2.d0 * ch9 * o1dx8 * h2ox8
392         
393
394c           R10: 2O + CO2 -> O2 + CO2
395
396            Lo3p(i,10) = 2.d0 * ch10 * o3px8 * co2x8 
397            Po2(i,10) = ch10 * o3px8 * o3px8 * co2x8
398
399
400c           R11: O + OH -> O2 + H
401
402            Lo3p(i,11) = ch11 * ohx8
403            Loh(i,11) = ch11 * o3px8
404            Po2(i,11) = ch11 * o3px8 * ohx8
405            Ph(i,11) = ch11 * o3px8 * ohx8
406
407
408c           R12: O2 + hv -> 2O
409
410            Lo2(i,12) = jdistot8(2,i)
411            Po3p(i,12) = 2.d0 * jdistot8(2,i) * o2x8
412
413
414c           R17(12b): O2 + hv -> O + O1D
415
416            Lo2(i,17) = jdistot8_b(2,i)
417            Po3p(i,17) = jdistot8_b(2,i) * o2x8
418            Po1d(i,17) = jdistot8_b(2,i) * o2x8
419
420
421c           R13: H + HO2 -> H2 + O2
422
423            Lh(i,13) = ch13 * ho2x8
424            Lho2(i,13) = ch13 * hx8
425            Ph2(i,13) = ch13 * hx8 * ho2x8
426            Po2(i,13) = ch13 * hx8 * ho2x8
427
428
429c           R14: O(1D) + H2 -> H + OH
430
431            Lo1d(i,14) = ch14 * h2x8
432            Lh2(i,14) = ch14 * o1dx8
433            Ph(i,14) = ch14 * o1dx8 * h2x8
434            Poh(i,14) = ch14 * o1dx8 * h2x8
435
436
437c           R15: OH + H2 -> H + H20
438
439            Loh(i,15) = ch15 * h2x8
440            Lh2(i,15) = ch15 * ohx8
441            Ph(i,15) = ch15 * ohx8 * h2x8
442            Ph2o(i,15) = ch15 * ohx8 * h2x8
443
444
445c           R18: OH + H2O2 -> H2O + HO2
446
447            Loh(i,18) = ch18 * h2o2x8
448            Lh2o2(i,18) = ch18 * ohx8
449            Ph2o(i,18) = ch18 * ohx8 * h2o2x8
450            Pho2(i,18) = ch18 * ohx8 * h2o2x8
451
452
453c           R19: O(1D) + CO2 -> O + CO2
454
455            Lo1d(i,19) = ch19 * co2x8
456            Po3p(i,19) = ch19 * o1dx8 * co2x8
457
458
459c           R20: O(1D) + O2 -> O + O2
460
461            Lo1d(i,20) = ch20 * o2x8
462            Po3p(i,20) = ch20 * o1dx8 * o2x8
463
464
465c           R21: O + O2 + CO2 -> O3 + CO2
466
467            Lo3p(i,21) = ch21 * o2x8 * co2x8
468            Lo2(i,21) = ch21 * o3px8 * co2x8
469
470
471c           R22: H2 + hv -> 2H
472
473            Lh2(i,22) = jdistot8(5,i)
474            Ph(i,22) = 2.d0 * jdistot8(5,i) * h2x8
475
476
477            do j=1,nreact
478               Pco2tot(i) = Pco2tot(i) + Pco2(i,j)
479               Lco2tot(i) = Lco2tot(i) + Lco2(i,j)
480               Po2tot(i) = Po2tot(i) + Po2(i,j)
481               Lo2tot(i) = Lo2tot(i) + Lo2(i,j)
482               Po3ptot(i) = Po3ptot(i) + Po3p(i,j)
483               Lo3ptot(i) = Lo3ptot(i) + Lo3p(i,j)
484               Pcotot(i) = Pcotot(i) + Pco(i,j)
485               Lcotot(i) = Lcotot(i) + Lco(i,j)
486               Phtot(i) = Phtot(i) + Ph(i,j)
487               Lhtot(i) = Lhtot(i) + Lh(i,j)
488               Pohtot(i) = Pohtot(i) + Poh(i,j)
489               Lohtot(i) = Lohtot(i) + Loh(i,j)
490               Pho2tot(i) = Pho2tot(i) + Pho2(i,j)
491               Lho2tot(i) = Lho2tot(i) + Lho2(i,j)
492               Ph2tot(i) = Ph2tot(i) + Ph2(i,j)
493               Lh2tot(i) = Lh2tot(i) + Lh2(i,j)
494               Ph2otot(i) = Ph2otot(i) + Ph2o(i,j)
495               Lh2otot(i) = Lh2otot(i) + Lh2o(i,j)
496               Po1dtot(i) = Po1dtot(i) + Po1d(i,j)
497               Lo1dtot(i) = Lo1dtot(i) + Lo1d(i,j)
498               Ph2o2tot(i) = Ph2o2tot(i) + Ph2o2(i,j)
499               Lh2o2tot(i) = Lh2o2tot(i) + Lh2o2(i,j)
500            end do
501
502
503            !New concentrations, implicit scheme
504            co2x8 = (co2x8 + Pco2tot(i) * deltat) /
505     $                (1.d0 + Lco2tot(i) * deltat)
506            o2x8 = (o2x8 + Po2tot(i) * deltat) /
507     $               (1.d0 + Lo2tot(i) * deltat)
508            o3px8 = (o3px8 + Po3ptot(i) * deltat) /
509     $                (1.d0 + Lo3ptot(i) * deltat)
510            cox8 = (cox8 + Pcotot(i) * deltat) /
511     $               (1.d0 + Lcotot(i) * deltat)
512            hx8 = (hx8 + Phtot(i) * deltat) /
513     $              (1.d0 + Lhtot(i) * deltat)
514             h2x8 = (h2x8 + Ph2tot(i) * deltat) /
515     $               (1.d0 + Lh2tot(i) * deltat)
516            h2ox8 = (h2ox8 + Ph2otot(i) * deltat) /
517     $                (1.d0 + Lh2otot(i) * deltat)
518            h2o2x8 = (h2o2x8 + Ph2o2tot(i) * deltat) /
519     $                (1.d0 + Lh2o2tot(i) * deltat)
520
521            if(o1d_eq(i).eq.'Y') then
522               auxp=jdistot8_b(1,i)*co2x8+jdistot8_b(2,i)*o2x8
523               auxl=ch9*h2ox8+ch14*h2x8+ch19*co2x8+ch20*o2x8
524               if(auxl.ne.0.) then
525                  o1dx8=auxp/auxl
526               else
527                  o1dx8=o1dx8+Po1dtot(i)*deltat
528               endif
529            else
530               o1dx8 = (o1dx8 + Po1dtot(i) * deltat) /
531     $                   (1.d0 + Lo1dtot(i) * deltat)
532            end if
533
534            if(oh_eq(i).eq.'Y') then
535               auxp=ch3*o3px8*ho2x8+jdistot(4,i)*h2ox8+
536     $              2.d0*ch9*o1dx8*h2ox8+ch14*o1dx8*h2x8+
537     $              2.d0*jdistot(6,i)*h2o2x8
538               auxl=ch4*cox8+ch7*ho2x8+ch11*o3px8+ch15*h2x8+
539     $              ch18*h2o2x8
540               if(auxl.ne.0.) then
541                  ohx8=auxp/auxl
542               else
543                  ohx8=ohx8+Pohtot(i)*deltat
544               end if
545            else
546               ohx8 = (ohx8 + Pohtot(i) * deltat) /
547     $                  (1.d0 + Lohtot(i) * deltat)
548            end if
549
550            if(ho2_eq(i).eq.'Y') then
551               auxp=ch2*hx8*o2x8*co2x8+ch18*ohx8*h2o2x8
552               auxl=ch3*o3px8+2.d0*ch5*ho2x8+ch7*ohx8+ch13*hx8
553               if(auxl.ne.0.) then
554                  ho2x8=auxp/auxl
555               else
556                  ho2x8=ho2x8+Pho2tot(i)*deltat
557               end if
558            else
559               ho2x8 = (ho2x8 + Pho2tot(i) * deltat) /
560     $                  (1.d0 + Lho2tot(i) * deltat)
561            end if
562
563
564C        End: temporal loop
565         end do
566         
567         
568         co2x(i)  = real(co2x8)
569         o2x(i)   = real(o2x8)
570         o3px(i)  = real(o3px8)
571         cox(i)   = real(cox8)
572         hx(i)    = real(hx8)
573         ohx(i)   = real(ohx8)
574         ho2x(i)  = real(ho2x8)
575         h2x(i)   = real(h2x8)
576         h2ox(i)  = real(h2ox8)
577         h2o2x(i) = real(h2o2x8)
578         o1dx(i)  = real(o1dx8)
579
580c         co2x(i)=(co2x8)
581c         o2x(i)=(o2x8)
582c         o3px(i)=(o3px8)
583c         cox(i)=(cox8)
584c         hx(i)=(hx8)
585c         ohx(i)=(ohx8)
586c         ho2x(i)=(ho2x8)
587c         h2x(i)=(h2x8)
588c         h2ox(i)=(h2ox8)
589c         h2o2x(i)=(h2o2x8)
590c         o1dx(i)=(o1dx8)
591
592C     End: altitude loop
593      end do
594
595      return
596
597      end
598
599
Note: See TracBrowser for help on using the repository browser.