source: trunk/MESOSCALE/LMDZ.MARS/libf_gcm/aeronomars/photochemist_B.F @ 1242

Last change on this file since 1242 was 57, checked in by aslmd, 14 years ago

mineur LMD_MM_MARS: ajout du GCM ancienne physique, systeme maintenant complet sur SVN (ne manque que la base de donnees d'etats initiaux)

File size: 50.0 KB
Line 
1c*****************************************************************
2c
3c     Photochemical routines for Scheme B
4c
5c     Author: Franck Lefevre
6c     ------
7c
8c     update 12/06/2003: compatibility with dust and iceparty
9c        (S. Lebonnois)
10c
11c*****************************************************************
12c
13      subroutine photochemist_B(lswitch,zycol, sza, ptimestep, press,
14c     $                          temp, dens, dist_sol)
15     $                          temp, dens, dist_sol, jo3)
16c
17c
18      implicit none
19c
20#include "dimensions.h"
21#include "dimphys.h"
22#include "chimiedata.h"
23#include "callkeys.h"
24c
25cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
26c     input/output:
27cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
28c     
29      real zycol(nlayermx,nqmx)  ! chemical species volume mixing ratio
30c if iceparty, zycol(l,nqmx-1) is ice surface area (in microns^2/cm^3)
31c
32cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
33c     inputs:
34cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
35c     
36      real sza                   ! solar zenith angle (deg)
37      real ptimestep             ! physics timestep (s)
38      real press(nlayermx)       ! pressure (hpa)
39      real temp(nlayermx)        ! temperature (k)
40      real dens(nlayermx)        ! density (cm-3)
41      real dist_sol              ! sun distance (au)
42c
43cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
44c     output:
45cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
46c     
47      real jo3(nlayermx)  ! photodissociation rate O3 -> O1D
48c
49cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
50c     local:
51cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
52c
53      integer phychemrat         ! ratio physics timestep/chemistry timestep
54      integer istep
55      integer i_o3,lev
56      integer nesp
57      integer lswitch
58c
59      parameter (nesp = 15)      ! number of species in the chemistry
60c
61      real ctimestep             ! chemistry timestep (s)
62      real rm(nlayermx,nesp)     ! species volume mixing ratio
63      real j(nlayermx,nd)        ! interpolated photolysis rates (s-1)
64      real icesurf(nlayermx)     ! ice surface area (cm^2/cm^3)
65      real rmo3(nlayermx)        ! ozone mixing ratio
66c
67c     reaction rates
68c
69      real a001(nlayermx), a002(nlayermx), a003(nlayermx)
70      real b001(nlayermx), b002(nlayermx), b003(nlayermx),
71     $     b004(nlayermx), b005(nlayermx), b006(nlayermx)
72      real c001(nlayermx), c002(nlayermx), c003(nlayermx),
73     $     c004(nlayermx), c005(nlayermx), c006(nlayermx),
74     $     c007(nlayermx), c008(nlayermx), c009(nlayermx),
75     $     c010(nlayermx), c011(nlayermx), c012(nlayermx),
76     $     c013(nlayermx), c014(nlayermx), c015(nlayermx),
77     $     c016(nlayermx), c017(nlayermx), c018(nlayermx)
78      real d001(nlayermx), d002(nlayermx), d003(nlayermx)
79      real e001(nlayermx), e002(nlayermx)
80      real h001(nlayermx), h002(nlayermx)
81c
82cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
83c     ctimestep : chemistry timestep (s)                         c
84cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
85c
86      ctimestep = 10.*60.
87c
88      phychemrat = int(ptimestep/ctimestep)
89c     print*,"PHYCHEMRAT = ",phychemrat
90c
91cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
92c     initialisation of chemical species and families            c
93cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
94c
95      call gcmtochim(zycol, lswitch, nesp, rm)
96c
97cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
98c     ice surface area for heterogenous chemistry
99cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
100c
101      if (iceparty) then
102       do lev=1,lswitch-1
103        icesurf(lev)= zycol(lev,nqmx-1)*1.e-8 ! microns^2 -> cm^2
104       enddo
105      else
106       do lev=1,lswitch-1
107        icesurf(lev)=0.0e0
108       enddo
109      endif
110c
111cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
112c     compute reaction rates                                     c
113cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
114c
115      call chemrates(lswitch, dens, press, temp, icesurf,
116     $               a001, a002, a003,
117     $               b001, b002, b003, b004, b005, b006,
118     $               c001, c002, c003, c004, c005, c006,
119     $               c007, c008, c009, c010, c011, c012,
120     $               c013, c014, c015, c016, c017, c018,
121     $               d001, d002, d003,
122     $               e001, e002,
123     $               h001, h002)
124c
125cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
126c     interpolation of photolysis rates in the lookup table      c
127cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
128c
129      i_o3 = 5
130c
131      do lev=1,lswitch-1
132        rmo3(lev)=rm(lev,i_o3)
133      enddo
134
135      call phot(lswitch, press, temp, sza, dist_sol, rmo3, j)
136c
137      do lev=1,lswitch-1
138        jo3(lev) = j(lev,5)
139      enddo
140
141cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
142c     loop over time within the physical timestep                c
143cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
144c
145      do istep = 1,phychemrat
146         call chimie_B(lswitch,nesp, rm, j, dens, ctimestep,
147     $                 press, temp, sza, dist_sol,
148     $                 a001, a002, a003,
149     $                 b001, b002, b003, b004, b005, b006,
150     $                 c001, c002, c003, c004, c005, c006,
151     $                 c007, c008, c009, c010, c011, c012,
152     $                 c013, c014, c015, c016, c017, c018,
153     $                 d001, d002, d003,
154     $                 e001, e002,
155     $                 h001, h002)
156c
157c        print*,'appel chimie ok iteration ', istep
158c
159      end do
160c
161cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
162c     save chemical species for the gcm                          c
163cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
164c
165      call chimtogcm(zycol, lswitch, nesp, rm)
166c
167      return
168      end
169c
170c*****************************************************************
171c
172      subroutine chimie_B(lswitch, nesp, rm, j, dens, dt,
173     $                    press, t, sza, dist_sol,
174     $                    a001, a002, a003,
175     $                    b001, b002, b003, b004, b005, b006,
176     $                    c001, c002, c003, c004, c005, c006,
177     $                    c007, c008, c009, c010, c011, c012,
178     $                    c013, c014, c015, c016, c017, c018,
179     $                    d001, d002, d003,
180     $                    e001, e002,
181     $                    h001, h002)
182c
183c*****************************************************************
184c
185      implicit none
186c
187#include "dimensions.h"
188#include "dimphys.h"
189#include "chimiedata.h"
190#include "callkeys.h"
191c
192cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
193c     input/output:
194cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
195c     
196      integer lswitch         ! interface level between chemistries
197      integer nesp            ! number of species
198c
199      real rm(nlayermx,nesp)  ! volume mixing ratios
200c
201cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
202c     inputs:
203cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
204c
205      real dens(nlayermx)     ! density (cm-3)
206      real dt                 ! chemistry timestep (s)
207      real j(nlayermx,nd)     ! interpolated photolysis rates (s-1)
208      real press(nlayermx)    ! pressure (hpa)
209      real t(nlayermx)        ! temperature (k)
210      real sza                ! solar zenith angle (deg)
211      real dist_sol           ! sun distance (au)
212c
213c     reaction rates
214c
215      real a001(nlayermx), a002(nlayermx), a003(nlayermx)
216      real b001(nlayermx), b002(nlayermx), b003(nlayermx),
217     $     b004(nlayermx), b005(nlayermx), b006(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)
226      real h001(nlayermx), h002(nlayermx)
227c
228cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
229c     local:
230cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
231c
232      integer HETERO
233      parameter (HETERO=0)
234
235      integer iesp, iter, l, niter
236      integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox,
237     $        i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2
238      integer j_o2_o, j_o2_o1d, j_co2_o, j_co2_o1d,
239     $        j_o3_o1d, j_o3_o, j_h2o, j_hdo, j_h2o2,
240     $        j_ho2, j_no2
241c
242      parameter (niter = 5)          ! iterations in the chemical scheme
243c
244      real*8 cc0(nlayermx,nesp)        ! initial number density (cm-3)
245      real*8 cc(nlayermx,nesp)         ! final number density (cm-3)
246      real*8 co2(nlayermx)             ! co2 number density (cm-3)
247      real*8 nox(nlayermx)             ! nox number density (cm-3)
248      real*8 no(nlayermx)              ! no  number density (cm-3)
249      real*8 no2(nlayermx)             ! no2 number density (cm-3)
250      real*8 production(nlayermx,nesp) ! production rate
251      real*8 loss(nlayermx,nesp)       ! loss rate
252c
253      real c007h(nlayermx)           ! c007 rate modified
254c                                    ! by heterogenous reactions h001 and h002
255c
256      real ro_o3, rh_ho2, roh_ho2, rno2_no
257c
258cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
259c     tracer numbering in the chemistry
260cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
261c
262      i_co2    =  1
263      i_co     =  2
264      i_o      =  3
265      i_o1d    =  4
266      i_o2     =  5
267      i_o3     =  6
268      i_h      =  7
269      i_h2     =  8
270      i_oh     =  9
271      i_ho2    = 10
272      i_h2o2   = 11
273      i_h2o    = 12
274      i_n2     = 13
275      i_hox    = 14
276      i_ox     = 15
277c
278cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
279c     numbering of photolysis rates
280cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
281c
282      j_o2_o     =  1      ! o2 + hv   -> o + o
283      j_o2_o1d   =  2      ! o2 + hv   -> o + o(1d)
284      j_co2_o    =  3      ! co2 + hv  -> co + o
285      j_co2_o1d  =  4      ! co2 + hv  -> co + o(1d)
286      j_o3_o1d   =  5      ! o3 + hv   -> o2 + o(1d)
287      j_o3_o     =  6      ! o3 + hv   -> o2 + o
288      j_h2o      =  7      ! h2o + hv  -> h + oh
289      j_hdo      =  8      ! hdo + hv  -> d + oh
290      j_h2o2     =  9      ! h2o2 + hv -> oh + oh
291      j_ho2      =  10     ! ho2 + hv  -> oh + o
292      j_no2      =  11     ! no2 + hv  -> no + o
293c
294cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
295c     volume mixing ratio -> density conversion
296cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
297c
298      do iesp = 1,nesp
299         do l = 1,lswitch-1
300            cc0(l,iesp) = rm(l,iesp)*dens(l)
301            cc(l,iesp)  = cc0(l,iesp)
302         end do
303      end do
304c
305cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
306c     co2 and nox number densities (cm-3)   
307c     
308c     co2 mixing ratio: 0.953
309c     nox mixing ratio: 6.e-10 (yung and demore, 1999)
310cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
311c
312      do l = 1,lswitch-1
313         co2(l) = cc(l,i_co2)
314         nox(l) = 6.e-10*dens(l)
315      end do
316c
317cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
318c     loop over iterations in the chemical scheme
319cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
320c
321      do iter = 1,niter
322c
323cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
324c     modification of c007
325c        by heterogenous reactions h001 and h002
326cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
327c
328      if (HETERO.eq.1) then
329       do l = 1,lswitch-1
330         c007h(l) = c007(l)
331     $            + h001(l)/max(cc(l,i_oh), 1.e-30*dens(l))
332     $            + h002(l)/max(cc(l,i_ho2),1.e-30*dens(l))
333c         write(*,*) "C007 = ",c007(l)," / C007H = ",c007h(l)
334       end do
335      else
336       do l = 1,lswitch-1
337         c007h(l) = c007(l)
338       end do
339      endif
340c
341cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
342c        nox species
343cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
344c        no2/no
345cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
346c
347      do l = 1,lswitch-1
348c
349         rno2_no = (d002(l)*cc(l,i_o3) + d003(l)*cc(l,i_ho2))
350     $            /(j(l,j_no2) +
351     $              d001(l)*max(cc(l,i_o),1.e-30*dens(l)))
352c
353cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
354c        no
355cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
356c
357         no(l) = nox(l)/(1. + rno2_no)
358c
359cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
360c        no2
361cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
362c
363         no2(l) = no(l)*rno2_no
364c
365      end do
366c
367cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
368c        hox species
369cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
370c        photochemical equilibrium for oh and ho2
371cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
372c        h/ho2
373cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
374c
375      do l = 1,lswitch-1
376c
377         rh_ho2 = (c001(l)*cc(l,i_o)
378     $           + c004(l)*cc(l,i_h)
379     $           + c005(l)*cc(l,i_h)
380     $           + c006(l)*cc(l,i_h)
381     $           + c007h(l)*cc(l,i_oh)
382     $           + 2.*c008(l)*cc(l,i_ho2)
383     $           + c015(l)*cc(l,i_o3)
384     $           + 2.*c016(l)*cc(l,i_ho2)
385     $           + j(l,j_ho2))
386     $          /(c011(l)*cc(l,i_o2))
387c
388cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
389c        oh/ho2
390cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
391c
392         roh_ho2 = (c001(l)*cc(l,i_o)
393     $            + c003(l)*cc(l,i_o3)*rh_ho2
394     $            + 2.*c004(l)*cc(l,i_h)
395     $            + 2.*c008(l)*cc(l,i_ho2)
396     $            + c015(l)*cc(l,i_o3)
397     $            + d003(l)*no(l)
398     $            + j(l,j_ho2)
399     $            + j(l,j_h2o)*cc(l,i_h2o)
400     $             /max(cc(l,i_ho2),dens(l)*1.e-30))
401     $            /(c002(l)*cc(l,i_o)
402     $            + c007h(l)*cc(l,i_ho2)
403     $            + e001(l)*cc(l,i_co))
404c
405cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
406c        h
407cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
408c
409         cc(l,i_h) = cc(l,i_hox)
410     $                 /(1. + (1. + roh_ho2)/rh_ho2)
411c
412cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
413c        ho2
414cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
415c
416         cc(l,i_ho2) = cc(l,i_h)/rh_ho2
417c
418cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
419c        oh
420cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
421c
422         cc(l,i_oh) = cc(l,i_ho2)*roh_ho2
423c
424      end do
425c
426cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
427c     modification of c007
428c        by heterogenous reactions h001 and h002
429c           after update of cc(l,i_oh) and cc(l,i_ho2)
430cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
431c
432      if (HETERO.eq.1) then
433       do l = 1,lswitch-1
434         c007h(l) = c007(l)
435     $            + h001(l)/max(cc(l,i_oh), 1.e-30*dens(l))
436     $            + h002(l)/max(cc(l,i_ho2),1.e-30*dens(l))
437c         write(*,*) "C007 = ",c007(l)," / C007H = ",c007h(l)
438       end do
439      endif
440c
441cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
442c        ox species
443cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
444c        day:
445c        - o1d at photochemical equilibrium
446c        - o3 at photochemical equilibrium
447c        - continuity equation for ox
448c        night:
449c        - o1d = 0
450c        - continuity equation for o3
451c        - continuity equation for o
452cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
453c       
454      if (sza .le. 95.) then
455c
456         do l = 1,lswitch-1
457c
458cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
459c        o(1d)
460cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
461c
462            cc(l,i_o1d) = (j(l,j_co2_o1d)*co2(l)
463     $                   + j(l,j_o2_o1d)*cc(l,i_o2)
464     $                   + j(l,j_o3_o1d)*cc(l,i_o3))
465     $                   /( b001(l)*co2(l)
466     $                   + b002(l)*cc(l,i_h2o)
467     $                   + b003(l)*cc(l,i_h2)
468     $                   + b004(l)*cc(l,i_o2)
469     $                   + b005(l)*cc(l,i_o3)
470     $                   + b006(l)*cc(l,i_o3))
471c
472cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
473c        o/o3
474cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
475c
476            ro_o3 = (j(l,j_o3_o1d) + j(l,j_o3_o)
477     $              + a003(l)*cc(l,i_o)
478     $              + c003(l)*cc(l,i_h)
479     $              + c014(l)*cc(l,i_oh)
480     $              + c015(l)*cc(l,i_ho2))
481     $              /(a001(l)*cc(l,i_o2))
482c
483cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
484c        o3
485cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
486c
487            cc(l,i_o3) = cc(l,i_ox)/(1. + ro_o3)
488c
489cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
490c        o
491cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
492c
493            cc(l,i_o) = cc(l,i_o3)*ro_o3
494c
495cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
496c        ox = o + o3
497cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
498c
499            production(l,i_ox) =
500     $                   + j(l,j_co2_o)*co2(l)
501     $                   + j(l,j_co2_o1d)*co2(l)
502     $                   + j(l,j_ho2)*cc(l,i_ho2)
503     $                   + 2.*j(l,j_o2_o)*cc(l,i_o2)
504     $                   + 2.*j(l,j_o2_o1d)*cc(l,i_o2)
505     $                   + c006(l)*cc(l,i_h)*cc(l,i_ho2)
506     $                   + c013(l)*cc(l,i_oh)*cc(l,i_oh)
507     $                   + d003(l)*cc(l,i_ho2)*no(l)
508c
509            loss(l,i_ox) = 2.*a002(l)*cc(l,i_o)*cc(l,i_o)
510     $                   + 2.*a003(l)*cc(l,i_o)*cc(l,i_o3)
511     $                   + c001(l)*cc(l,i_ho2)*cc(l,i_o)
512     $                   + c002(l)*cc(l,i_oh)*cc(l,i_o)
513     $                   + c003(l)*cc(l,i_h)*cc(l,i_o3)
514     $                   + c012(l)*cc(l,i_o)*cc(l,i_h2o2)
515     $                   + c014(l)*cc(l,i_o3)*cc(l,i_oh)
516     $                   + c015(l)*cc(l,i_o3)*cc(l,i_ho2)
517     $                   + d001(l)*cc(l,i_o)*no2(l)
518     $                   + e002(l)*cc(l,i_o)*cc(l,i_co)
519c
520            loss(l,i_ox) = loss(l,i_ox)/cc(l,i_ox)
521c
522         end do
523c
524      else
525c
526         do l = 1,lswitch-1
527c
528cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
529c        o(1d)
530cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
531c
532            cc(l,i_o1d) = 0.
533c
534cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
535c        o3
536cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
537c
538            production(l,i_o3) = a001(l)*cc(l,i_o2)*cc(l,i_o)
539c
540            loss(l,i_o3) = a003(l)*cc(l,i_o)
541     $                   + c003(l)*cc(l,i_h)
542     $                   + c014(l)*cc(l,i_oh)
543     $                   + c015(l)*cc(l,i_ho2)
544c
545
546cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
547c        o
548cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
549c
550            production(l,i_o) = c006(l)*cc(l,i_h)*cc(l,i_ho2)
551     $                        + c013(l)*cc(l,i_oh)*cc(l,i_oh)
552c
553            loss(l,i_o)  =  a001(l)*cc(l,i_o2)
554     $                   + 2.*a002(l)*cc(l,i_o)
555     $                   + a003(l)*cc(l,i_o3)
556     $                   + c001(l)*cc(l,i_ho2)
557     $                   + c002(l)*cc(l,i_oh)
558     $                   + c012(l)*cc(l,i_h2o2)
559     $                   + e002(l)*cc(l,i_co)
560c
561         end do
562      end if
563c
564cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
565c     other species
566cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
567c
568      do l = 1,lswitch-1
569c
570cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
571c        co
572cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
573c
574         production(l,i_co) = j(l,j_co2_o)*co2(l)
575     $                      + j(l,j_co2_o1d)*co2(l)
576c
577         loss(l,i_co) = e001(l)*cc(l,i_oh)
578     $                + e002(l)*cc(l,i_o)
579c
580cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
581c        o2
582cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
583c
584         production(l,i_o2) =
585     $                  j(l,j_o3_o)*cc(l,i_o3)
586     $                + j(l,j_o3_o1d)*cc(l,i_o3)
587     $                + a002(l)*cc(l,i_o)*cc(l,i_o)
588     $                + 2.*a003(l)*cc(l,i_o)*cc(l,i_o3)
589     $                + 2.*b005(l)*cc(l,i_o1d)*cc(l,i_o3)
590     $                + b006(l)*cc(l,i_o1d)*cc(l,i_o3)
591     $                + c001(l)*cc(l,i_o)*cc(l,i_ho2)
592     $                + c002(l)*cc(l,i_o)*cc(l,i_oh)
593     $                + c003(l)*cc(l,i_h)*cc(l,i_o3)
594     $                + c005(l)*cc(l,i_h)*cc(l,i_ho2)
595     $                + c007h(l)*cc(l,i_oh)*cc(l,i_ho2)
596     $                + c008(l)*cc(l,i_ho2)*cc(l,i_ho2)
597     $                + c014(l)*cc(l,i_o3)*cc(l,i_oh)
598     $                + 2.*c015(l)*cc(l,i_o3)*cc(l,i_ho2)
599     $                + c016(l)*cc(l,i_ho2)*cc(l,i_ho2)
600     $                + d001(l)*cc(l,i_o)*no2(l)
601c
602         loss(l,i_o2) = j(l,j_o2_o)
603     $                + j(l,j_o2_o1d)
604     $                + a001(l)*cc(l,i_o)
605     $                + c011(l)*cc(l,i_h)
606c
607cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
608c        h2
609cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
610c
611         production(l,i_h2) = c005(l)*cc(l,i_h)*cc(l,i_ho2)
612     $                      + c018(l)*cc(l,i_h)*cc(l,i_h)
613c
614         loss(l,i_h2) = b003(l)*cc(l,i_o1d)
615     $                + c010(l)*cc(l,i_oh)
616c
617cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
618c        h2o
619cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
620c
621         production(l,i_h2o) =
622     $                   c006(l)*cc(l,i_h)*cc(l,i_ho2)
623     $                 + c007h(l)*cc(l,i_oh)*cc(l,i_ho2)
624     $                 + c009(l)*cc(l,i_oh)*cc(l,i_h2o2)
625     $                 + c010(l)*cc(l,i_oh)*cc(l,i_h2)
626     $                 + c013(l)*cc(l,i_oh)*cc(l,i_oh)
627c
628         loss(l,i_h2o) = j(l,j_h2o)
629     $                 + b002(l)*cc(l,i_o1d)
630c
631cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
632c        h2o2
633cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
634c
635         production(l,i_h2o2) =
636     $                    c008(l)*cc(l,i_ho2)*cc(l,i_ho2)
637     $                  + c016(l)*cc(l,i_ho2)*cc(l,i_ho2)
638     $                  + c017(l)*cc(l,i_oh)*cc(l,i_oh)
639c
640         loss(l,i_h2o2) = j(l,j_h2o2)
641     $                  + c009(l)*cc(l,i_oh)
642     $                  + c012(l)*cc(l,i_o)
643c
644cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
645c        hox = h + oh + ho2
646cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
647c
648         production(l,i_hox) =
649     $                   2.*j(l,j_h2o)*cc(l,i_h2o)
650     $                 + 2.*j(l,j_h2o2)*cc(l,i_h2o2)
651     $                 + 2.*b002(l)*cc(l,i_o1d)*cc(l,i_h2o)
652     $                 + 2.*b003(l)*cc(l,i_o1d)*cc(l,i_h2)
653     $                 + 2.*c012(l)*cc(l,i_o)*cc(l,i_h2o2)
654c
655         loss(l,i_hox) = 2.*c005(l)*cc(l,i_h)*cc(l,i_ho2)
656     $                 + 2.*c006(l)*cc(l,i_h)*cc(l,i_ho2)
657     $                 + 2.*c007h(l)*cc(l,i_oh)*cc(l,i_ho2)
658     $                 + 2.*c008(l)*cc(l,i_ho2)*cc(l,i_ho2)
659     $                 + 2.*c013(l)*cc(l,i_oh)*cc(l,i_oh)
660     $                 + 2.*c016(l)*cc(l,i_ho2)*cc(l,i_ho2)
661     $                 + 2.*c017(l)*cc(l,i_oh)*cc(l,i_oh)
662     $                 + 2.*c018(l)*cc(l,i_h)*cc(l,i_h)
663c
664         loss(l,i_hox) = loss(l,i_hox)/cc(l,i_hox)
665c
666      end do
667c
668cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
669c     update number densities
670cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
671c
672c     long-lived species
673c
674      do l = 1,lswitch-1
675         cc(l,i_co)  = (cc0(l,i_co) + production(l,i_co)*dt)
676     $              /(1. + loss(l,i_co)*dt)
677         cc(l,i_o2)  = (cc0(l,i_o2) + production(l,i_o2)*dt)
678     $              /(1. + loss(l,i_o2)*dt)
679         cc(l,i_h2)  = (cc0(l,i_h2) + production(l,i_h2)*dt)
680     $              /(1. + loss(l,i_h2)*dt)
681         cc(l,i_h2o2)= (cc0(l,i_h2o2) + production(l,i_h2o2)*dt)
682     $              /(1. + loss(l,i_h2o2)*dt)
683         cc(l,i_h2o) = (cc0(l,i_h2o) + production(l,i_h2o)*dt)
684     $              /(1. + loss(l,i_h2o)*dt)
685         cc(l,i_hox) = (cc0(l,i_hox) + production(l,i_hox)*dt)
686     $              /(1. + loss(l,i_hox)*dt)
687      end do
688c
689c     ox species
690c
691      if (sza .le. 95.) then
692         do l = 1,lswitch-1
693            cc(l,i_ox) = (cc0(l,i_ox) + production(l,i_ox)*dt)
694     $                   /(1. + loss(l,i_ox)*dt)
695         end do
696      else
697         do l = 1,lswitch-1
698            cc(l,i_o)  = (cc0(l,i_o) + production(l,i_o)*dt)
699     $                   /(1. + loss(l,i_o)*dt)
700            cc(l,i_o3) = (cc0(l,i_o3) + production(l,i_o3)*dt)
701     $                   /(1. + loss(l,i_o3)*dt)
702            cc(l,i_ox) = cc(l,i_o) + cc(l,i_o3)
703         end do
704      end if
705c
706cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
707c     end of loop over iterations
708cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
709c
710      end do
711c
712cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
713c     density -> volume mixing ratio conversion
714cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
715c
716      do iesp = 1,nesp
717         do l = 1,lswitch-1
718            rm(l,iesp) = max(cc(l,iesp)/dens(l), 1.e-30)
719         end do
720      end do
721c
722      return
723      end
724c
725c*****************************************************************
726c
727      subroutine phot(lswitch, press, temp, sza, dist_sol, rmo3, j)
728c
729c*****************************************************************
730c
731      implicit none
732c
733#include "dimensions.h"
734#include "dimphys.h"
735#include "chimiedata.h"
736c
737cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
738c     inputs:
739cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
740c       
741      integer lswitch         ! interface level between chemistries
742      real press(nlayermx)             ! pressure (hPa)
743      real temp(nlayermx)              ! temperature (K)
744      real sza                         ! solar zenith angle (deg)
745      real dist_sol                    ! sun distance (AU)
746      real rmo3(nlayermx)              ! ozone mixing ratio
747c
748cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
749c     output:
750cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
751c
752      real j(nlayermx,nd)              ! interpolated photolysis rates (s-1)
753c
754cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
755c     local:
756cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
757c
758      integer icol, ij, indsza, indcol, indozo, indtemp,
759     $        iozo, isza, it, iztop, l
760c
761      real col(nlayermx)               ! overhead air column   (molecule cm-2)
762      real colo3(nlayermx)             ! overhead ozone column (molecule cm-2)
763      real poids(2,2,2,2)              ! 4D interpolation weights
764      real tref                        ! temperature  at 1.9 hPa in the gcm (K)
765      real table_temp(ntemp)           ! temperatures at 1.9 hPa in jmars   (K)
766      real cinf, csup, cicol,
767     $     ciozo, cisza, citemp
768      real colo3min, dp, scaleh
769      real ratio_o3(nlayermx)
770c
771ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
772c     day/night criterion
773ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
774c
775      if (sza .le. 95.) then
776c
777ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
778c     temperatures at 1.9 hPa in jmars
779ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
780c     
781      table_temp(1) = 226.2
782      table_temp(2) = 206.2
783      table_temp(3) = 186.2
784      table_temp(4) = 169.8
785c
786ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
787c     interpolation in solar zenith angle
788ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
789c
790      indsza = nsza - 1
791      do isza = 1,nsza
792         if (szatab(isza) .ge. sza) then
793            indsza = min(indsza,isza - 1)
794            indsza = max(indsza, 1)
795         end if
796      end do
797      cisza = (sza - szatab(indsza))
798     $       /(szatab(indsza + 1) - szatab(indsza))
799c
800ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
801c     air and ozone columns
802ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
803c
804c     air column at model top: assign standard value
805c
806      scaleh = 10.
807      iztop = nint(scaleh*log(press(1)/press(lswitch-1)))
808      iztop = min(iztop,200)
809      col(lswitch-1) = colairtab(iztop)
810c
811c     ozone column at model top
812c
813      colo3(lswitch-1) = 0.
814c
815c     air and ozone columns for other levels
816c
817      do l = lswitch-2,1,-1
818         dp = (press(l) - press(l+1))*100.
819         col(l) = col(l+1) + factor*dp
820         col(l) = min(col(l), colairtab(0))
821         colo3(l) = colo3(l+1)
822     $            + factor*(rmo3(l+1) + rmo3(l))*0.5*dp
823      end do
824c
825c     ratio ozone column/minimal theoretical column (ro3 = 7.171e-10)
826c
827      do l = 1,lswitch-1
828         colo3min    = col(l)*7.171e-10
829         ratio_o3(l) = colo3(l)/colo3min
830         ratio_o3(l) = min(ratio_o3(l), table_ozo(nozo)*10.)
831         ratio_o3(l) = max(ratio_o3(l), table_ozo(1)*10.)
832      end do
833c
834ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
835c     temperature dependence
836ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
837c
838c     1) search for temperature at 1.9 hPa (tref): vertical interpolation
839c
840c      print*,'****** recherche de T a 1.9 hPa ******'
841      tref = temp(1)
842      do l = (lswitch-1)-1,1,-1
843         if (press(l) .gt. 1.9) then
844            cinf = (press(l) - 1.9)
845     $            /(press(l) - press(l+1))
846            csup = 1. - cinf
847            tref = cinf*temp(l+1) + csup*temp(l)
848c
849c            print*, press(l+1), temp(l+1)
850c            print*, press(l  ), temp(l  )
851c            print*, 'tref = ', tref
852            go to 10
853         end if
854      end do
855 10   continue
856c
857c     2) interpolation in temperature
858c
859      tref = min(tref,table_temp(1))
860      tref = max(tref,table_temp(ntemp))
861c
862      do it = 2, ntemp
863         if (table_temp(it) .le. tref) then
864            citemp = (log(tref) - log(table_temp(it)))
865     $              /(log(table_temp(it-1)) - log(table_temp(it)))
866            indtemp = it - 1
867c            print*,'interpolation entre ', table_temp(it), ' et ',
868c     $              table_temp(it-1)
869c            print*,'indtemp = ', indtemp, ' citemp = ', citemp
870            goto 20
871         end if
872      end do
873  20  continue
874c      print*,'****** fin interpolation en T ********'
875c
876ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
877c     loop over vertical levels
878ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
879c
880      do l = 1,lswitch-1
881c
882c     interpolation in air column
883c
884         do icol = 0,200
885            if (colairtab(icol) .lt. col(l)) then
886               cicol = (log(col(l)) - log(colairtab(icol)))
887     $                /(log(colairtab(icol-1)) - log(colairtab(icol)))
888               indcol = icol - 1
889               goto 30
890            end if
891         end do
892 30      continue
893c
894cc    interpolation in ozone column
895c
896         indozo = nozo - 1
897         do iozo = 1,nozo
898            if (table_ozo(iozo)*10. .ge. ratio_o3(l)) then
899               indozo = min(indozo, iozo - 1)
900               indozo = max(indozo, 1)
901            end if
902         end do
903         ciozo = (ratio_o3(l) - table_ozo(indozo)*10.)
904     $          /(table_ozo(indozo + 1)*10. - table_ozo(indozo)*10.)
905c
906c         write(6,'(i2,4e13.4,5f10.4)')
907c    $         l, press(l), col(l),
908c    $         colairtab(indcol), colairtab(indcol+1),
909c    $         colo3(l)/2.687e15,
910c    $         ratio_o3(l), table_ozo(indozo), table_ozo(indozo+1),
911c    $         ciozo
912c
913cc    4-dimensional interpolation weights
914c
915         poids(1,1,1,1) = citemp*(1.-cisza)*cicol*(1.-ciozo)
916         poids(1,1,1,2) = citemp*(1.-cisza)*cicol*ciozo
917         poids(1,1,2,1) = citemp*(1.-cisza)*(1.-cicol)*(1.-ciozo)
918         poids(1,1,2,2) = citemp*(1.-cisza)*(1.-cicol)*ciozo
919         poids(1,2,1,1) = citemp*cisza*cicol*(1.-ciozo)
920         poids(1,2,1,2) = citemp*cisza*cicol*ciozo
921         poids(1,2,2,1) = citemp*cisza*(1.-cicol)*(1.-ciozo)
922         poids(1,2,2,2) = citemp*cisza*(1.-cicol)*ciozo
923         poids(2,1,1,1) = (1.-citemp)*(1.-cisza)*cicol*(1.-ciozo)
924         poids(2,1,1,2) = (1.-citemp)*(1.-cisza)*cicol*ciozo
925         poids(2,1,2,1) = (1.-citemp)*(1.-cisza)*(1.-cicol)*(1.-ciozo)
926         poids(2,1,2,2) = (1.-citemp)*(1.-cisza)*(1.-cicol)*ciozo
927         poids(2,2,1,1) = (1.-citemp)*cisza*cicol*(1.-ciozo)
928         poids(2,2,1,2) = (1.-citemp)*cisza*cicol*ciozo
929         poids(2,2,2,1) = (1.-citemp)*cisza*(1.-cicol)*(1.-ciozo)
930         poids(2,2,2,2) = (1.-citemp)*cisza*(1.-cicol)*ciozo
931c
932cc    4-dimensional interpolation in the lookup table
933c
934         do ij = 1,nd
935           j(l,ij) =
936     $     poids(1,1,1,1)*jphot(indtemp,indsza,indcol,indozo,ij)
937     $   + poids(1,1,1,2)*jphot(indtemp,indsza,indcol,indozo+1,ij)
938     $   + poids(1,1,2,1)*jphot(indtemp,indsza,indcol+1,indozo,ij)
939     $   + poids(1,1,2,2)*jphot(indtemp,indsza,indcol+1,indozo+1,ij)
940     $   + poids(1,2,1,1)*jphot(indtemp,indsza+1,indcol,indozo,ij)
941     $   + poids(1,2,1,2)*jphot(indtemp,indsza+1,indcol,indozo+1,ij)
942     $   + poids(1,2,2,1)*jphot(indtemp,indsza+1,indcol+1,indozo,ij)
943     $   + poids(1,2,2,2)*jphot(indtemp,indsza+1,indcol+1,indozo+1,ij)
944     $   + poids(2,1,1,1)*jphot(indtemp+1,indsza,indcol,indozo,ij)
945     $   + poids(2,1,1,2)*jphot(indtemp+1,indsza,indcol,indozo+1,ij)
946     $   + poids(2,1,2,1)*jphot(indtemp+1,indsza,indcol+1,indozo,ij)
947     $   + poids(2,1,2,2)*jphot(indtemp+1,indsza,indcol+1,indozo+1,ij)
948     $   + poids(2,2,1,1)*jphot(indtemp+1,indsza+1,indcol,indozo,ij)
949     $   + poids(2,2,1,2)*jphot(indtemp+1,indsza+1,indcol,indozo+1,ij)
950     $   + poids(2,2,2,1)*jphot(indtemp+1,indsza+1,indcol+1,indozo,ij)
951     $   + poids(2,2,2,2)*jphot(indtemp+1,indsza+1,indcol+1,indozo+1,ij)
952         end do
953c
954cc    correction for sun distance
955c
956         do ij = 1,nd
957            j(l,ij) = j(l,ij)*(1.52/dist_sol)**2.
958         end do
959c
960ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
961c     end of loop over vertical levels
962ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
963c
964      end do
965c
966      else
967c
968ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
969c     night
970ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
971c
972      do ij = 1,nd
973         do l = 1,lswitch-1
974            j(l,ij) = 0.
975         end do
976      end do
977c
978      end if
979c
980      return
981      end
982c
983c*****************************************************************
984c
985      subroutine gcmtochim(zycol, lswitch, nesp, rm)
986c
987c*****************************************************************
988c
989      implicit none
990c
991#include "dimensions.h"
992#include "dimphys.h"
993#include "callkeys.h"
994c
995cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
996c     inputs:
997cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
998c     
999      real zycol(nlayermx,nqmx)! species volume mixing ratio in the gcm
1000c
1001      integer nesp            ! number of species in the chemistry
1002      integer lswitch         ! interface level between chemistries
1003c
1004cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1005c     outputs:
1006cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1007c
1008      real rm(nlayermx,nesp)   ! species volume mixing ratio
1009c     
1010cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1011c     local:
1012cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1013c
1014      integer l
1015      integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox,
1016     $        i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2
1017c     
1018cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1019c     tracers numbering in the gcm
1020cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1021c
1022c     co2      =  nqchem_min           
1023c     co       =  nqchem_min + 1         
1024c     o        =  nqchem_min + 2
1025c     o(1d)    =  nqchem_min + 3
1026c     o2       =  nqchem_min + 4
1027c     o3       =  nqchem_min + 5
1028c     h        =  nqchem_min + 6
1029c     h2       =  nqchem_min + 7
1030c     oh       =  nqchem_min + 8
1031c     ho2      =  nqchem_min + 9
1032c     h2o2     =  nqchem_min + 10
1033c     n2       =  nqchem_min + 11
1034c     ar       =  nqchem_min + 12
1035c     h2o      =  nqmx
1036c     
1037cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1038c     tracers numbering in the chemistry
1039cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1040c
1041      i_co2    =  1
1042      i_co     =  2
1043      i_o      =  3
1044      i_o1d    =  4
1045      i_o2     =  5
1046      i_o3     =  6
1047      i_h      =  7
1048      i_h2     =  8
1049      i_oh     =  9
1050      i_ho2    = 10
1051      i_h2o2   = 11
1052      i_h2o    = 12
1053      i_n2     = 13
1054      i_hox    = 14
1055      i_ox     = 15
1056c
1057cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1058c     initialise chemical species
1059cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1060c
1061      do l = 1,lswitch-1
1062       rm(l,i_co2)  = max(zycol(l, nqchem_min),      1.e-30)
1063       rm(l,i_co)   = max(zycol(l, nqchem_min + 1),  1.e-30)
1064       rm(l,i_o)    = max(zycol(l, nqchem_min + 2),  1.e-30)
1065       rm(l,i_o1d)  = max(zycol(l, nqchem_min + 3),  1.e-30)
1066       rm(l,i_o2)   = max(zycol(l, nqchem_min + 4),  1.e-30)
1067       rm(l,i_o3)   = max(zycol(l, nqchem_min + 5),  1.e-30)
1068       rm(l,i_h)    = max(zycol(l, nqchem_min + 6),  1.e-30)
1069       rm(l,i_h2)   = max(zycol(l, nqchem_min + 7),  1.e-30)
1070       rm(l,i_oh)   = max(zycol(l, nqchem_min + 8),  1.e-30)
1071       rm(l,i_ho2)  = max(zycol(l, nqchem_min + 9),  1.e-30)
1072       rm(l,i_h2o2) = max(zycol(l, nqchem_min + 10), 1.e-30)
1073       rm(l,i_n2)   = max(zycol(l, nqchem_min + 11), 1.e-30)
1074       rm(l,i_h2o)  = max(zycol(l, nqmx),            1.e-30)
1075c
1076      end do
1077c
1078cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1079c     initialise chemical families                     c
1080cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1081c
1082      do l = 1,lswitch-1
1083         rm(l,i_hox) = rm(l,i_h)
1084     $               + rm(l,i_oh)
1085     $               + rm(l,i_ho2)
1086         rm(l,i_ox)  = rm(l,i_o)
1087     $               + rm(l,i_o3)
1088      end do
1089c
1090      return
1091      end
1092c
1093c*****************************************************************
1094c
1095      subroutine chimtogcm(zycol, lswitch, nesp, rm)
1096c
1097c*****************************************************************
1098c
1099      implicit none
1100c
1101#include "dimensions.h"
1102#include "dimphys.h"
1103#include "callkeys.h"
1104c
1105cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1106c     inputs:
1107cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1108c
1109      integer nesp               ! number of species in the chemistry
1110      integer lswitch            ! interface level between chemistries
1111c     
1112      real rm(nlayermx,nesp)     ! species volume mixing ratio
1113c
1114cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1115c     output:
1116cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1117c     
1118      real zycol(nlayermx,nqmx)  ! species volume mixing ratio in the gcm
1119c
1120cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1121c     local:
1122cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1123c
1124      integer l
1125      integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox,
1126     $        i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2
1127c     
1128cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1129c     tracers numbering in the gcm
1130cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1131c
1132c     co2      =  nqchem_min
1133c     co       =  nqchem_min + 1
1134c     o        =  nqchem_min + 2
1135c     o(1d)    =  nqchem_min + 3
1136c     o2       =  nqchem_min + 4
1137c     o3       =  nqchem_min + 5
1138c     h        =  nqchem_min + 6
1139c     h2       =  nqchem_min + 7
1140c     oh       =  nqchem_min + 8
1141c     ho2      =  nqchem_min + 9
1142c     h2o2     =  nqchem_min + 10
1143c     n2       =  nqchem_min + 11
1144c     ar       =  nqchem_min + 12
1145c     h2o      =  nqmx
1146c     
1147cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1148c     tracers numbering in the chemistry
1149cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1150c
1151      i_co2    =  1
1152      i_co     =  2
1153      i_o      =  3
1154      i_o1d    =  4
1155      i_o2     =  5
1156      i_o3     =  6
1157      i_h      =  7
1158      i_h2     =  8
1159      i_oh     =  9
1160      i_ho2    = 10
1161      i_h2o2   = 11
1162      i_h2o    = 12
1163      i_n2     = 13
1164      i_hox    = 14
1165      i_ox     = 15
1166c
1167cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1168c     save mixing ratios for the gcm
1169cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1170c
1171      do l = 1,lswitch-1
1172         zycol(l, nqchem_min)      = rm(l,i_co2)
1173         zycol(l, nqchem_min + 1)  = rm(l,i_co)
1174         zycol(l, nqchem_min + 2)  = rm(l,i_o)
1175         zycol(l, nqchem_min + 3)  = rm(l,i_o1d)
1176         zycol(l, nqchem_min + 4)  = rm(l,i_o2)
1177         zycol(l, nqchem_min + 5)  = rm(l,i_o3)
1178         zycol(l, nqchem_min + 6)  = rm(l,i_h) 
1179         zycol(l, nqchem_min + 7)  = rm(l,i_h2)
1180         zycol(l, nqchem_min + 8)  = rm(l,i_oh)
1181         zycol(l, nqchem_min + 9)  = rm(l,i_ho2)
1182         zycol(l, nqchem_min + 10) = rm(l,i_h2o2)
1183         zycol(l, nqchem_min + 11) = rm(l,i_n2)
1184         zycol(l, nqmx)            = rm(l,i_h2o)
1185      end do
1186c
1187c  water ice: zycol(nqmx-1) has not been changed if iceparty on
1188
1189      return
1190      end
1191c
1192c*****************************************************************
1193c
1194      subroutine chemrates(lswitch, dens, press, t, icesurf,
1195     $                     a001, a002, a003,
1196     $                     b001, b002, b003, b004, b005, b006,
1197     $                     c001, c002, c003, c004, c005, c006,
1198     $                     c007, c008, c009, c010, c011, c012,
1199     $                     c013, c014, c015, c016, c017, c018,
1200     $                     d001, d002, d003,
1201     $                     e001, e002,
1202     $                     h001, h002)
1203c
1204c*****************************************************************
1205c
1206      implicit none
1207c
1208#include "dimensions.h"
1209#include "dimphys.h"
1210c
1211cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1212c     inputs:                                                    c
1213cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1214c
1215      integer lswitch                  ! interface level between chemistries
1216
1217      real dens(nlayermx)              ! density (cm-3)
1218      real press(nlayermx)             ! pressure (hpa)
1219      real t(nlayermx)                 ! temperature (k)
1220      real icesurf(nlayermx)           ! ice surface area (cm^2/cm^3)
1221c
1222cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1223c     outputs:                                                   c
1224cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1225c
1226      real a001(nlayermx), a002(nlayermx), a003(nlayermx)
1227      real b001(nlayermx), b002(nlayermx), b003(nlayermx),
1228     $     b004(nlayermx), b005(nlayermx), b006(nlayermx)
1229      real c001(nlayermx), c002(nlayermx), c003(nlayermx),
1230     $     c004(nlayermx), c005(nlayermx), c006(nlayermx),
1231     $     c007(nlayermx), c008(nlayermx), c009(nlayermx),
1232     $     c010(nlayermx), c011(nlayermx), c012(nlayermx),
1233     $     c013(nlayermx), c014(nlayermx), c015(nlayermx),
1234     $     c016(nlayermx), c017(nlayermx), c018(nlayermx)
1235      real d001(nlayermx), d002(nlayermx), d003(nlayermx)
1236      real e001(nlayermx), e002(nlayermx)
1237      real h001(nlayermx), h002(nlayermx)
1238c
1239cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1240c     local:                                                     c
1241cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1242c
1243      real ak0, ak1, rate, xpo
1244c
1245      integer l
1246c
1247cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1248c     compute reaction rates from JPL 1997, otherwise mentioned
1249cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1250c
1251      do l = 1,lswitch-1
1252c
1253cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1254c        oxygen compounds
1255cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1256c
1257ccc      a001: o + o2 + co2 -> o3 + co2
1258c
1259c        jpl 2003
1260c
1261         a001(l) = 2.5
1262     $             *6.0e-34*(t(l)/300.)**(-2.4)*dens(l)
1263c
1264c        nair et al., 1994
1265c
1266c        a001(l) = 1.3e-34*exp(724./t(l))*dens(l)
1267c
1268ccc      a002: o + o + co2 -> o2 + co2
1269c
1270c        Tsang and Hampson, J. Chem. Phys. Ref. Data, 15, 1087, 1986
1271c
1272         a002(l) = 2.5*5.2e-35*exp(900./t(l))*dens(l)
1273c
1274c        Campbell and Gray, Chem. Phys. Lett., 18, 607, 1973
1275c
1276c        a002(l) = 1.2e-32*(300./t(l))**(2.0)*dens(l)
1277c
1278ccc      a003: o + o3 -> o2 + o2
1279c
1280c        jpl 2003
1281c
1282         a003(l) = 8.0e-12*exp(-2060./t(l))
1283c
1284cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1285c        reactions with o(1d)
1286cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1287c
1288ccc      b001: o(1d) + co2  -> o + co2
1289c
1290c        jpl 2003
1291c
1292         b001(l) = 7.4e-11*exp(120./t(l))
1293c
1294ccc      b002: o(1d) + h2o  -> oh + oh
1295c
1296c        jpl 2003
1297c
1298         b002(l) = 2.2e-10
1299c   
1300ccc      b003: o(1d) + h2  -> oh + h
1301c
1302c        jpl 2003
1303c
1304         b003(l) = 1.1e-10
1305c
1306c        nair et al., 1994
1307c
1308c        b003(l) = 1.0e-10
1309c   
1310ccc      b004: o(1d) + o2  -> o + o2
1311c
1312c        jpl 2003
1313c
1314         b004(l) = 3.2e-11*exp(70./t(l))
1315c   
1316ccc      b005: o(1d) + o3  -> o2 + o2
1317c
1318c        jpl 2003
1319c
1320         b005(l) = 1.2e-10
1321c   
1322ccc      b006: o(1d) + o3  -> o2 + o + o
1323c
1324c        jpl 2003
1325c
1326         b006(l) = 1.2e-10
1327c
1328cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1329c        hydrogen compounds   
1330cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1331c
1332ccc      c001: o + ho2 -> oh + o2
1333c
1334c        jpl 2003
1335c
1336         c001(l) = 3.0e-11*exp(200./t(l))
1337c        c001(l) = 0.75*3.0e-11*exp(200./t(l))
1338c
1339ccc      c002: o + oh -> o2 + h
1340c
1341c        jpl 2003
1342c
1343         c002(l) = 2.2e-11*exp(120./t(l))
1344c
1345ccc      c003: h + o3 -> oh + o2
1346c
1347c        jpl 2003
1348c
1349         c003(l) = 1.4e-10*exp(-470./t(l))
1350c
1351ccc      c004: h + ho2 -> oh + oh
1352c
1353c        jpl 2003
1354c
1355         c004(l) = 8.1e-11*0.90
1356c
1357ccc      c005: h + ho2 -> h2 + o2
1358c
1359c        jpl 2003
1360c
1361         c005(l) = 8.1e-11*0.08
1362c
1363ccc      c006: h + ho2 -> h2o + o
1364c
1365c        jpl 2003
1366c
1367         c006(l) = 8.1e-11*0.02
1368c
1369ccc      c007: oh + ho2 -> h2o + o2
1370c
1371c        jpl 2003
1372c
1373         c007(l) = 4.8e-11*exp(250./t(l))
1374c        c007(l) = 0.75*4.8e-11*exp(250./t(l))
1375c
1376ccc      c008: ho2 + ho2 -> h2o2 + o2
1377c
1378c        jpl 2003
1379c
1380c        c008(l) = 2.3e-13*exp(600./t(l))
1381c
1382c        christensen et al., grl, 13, 2002
1383c
1384         c008(l) = 1.5e-12*exp(19./t(l))
1385c
1386ccc      c009: oh + h2o2 -> h2o + ho2
1387c
1388c        jpl 2003
1389c
1390         c009(l) = 2.9e-12*exp(-160./t(l))
1391c
1392ccc      c010: oh + h2 -> h2o + h
1393c
1394c        jpl 2003
1395c
1396         c010(l) = 5.5e-12*exp(-2000./t(l))
1397c
1398ccc      c011: h + o2 + co2 -> ho2 + co2
1399c
1400c        jpl 2003
1401c
1402         ak0 = 2.5*5.7e-32*(t(l)/300.)**(-1.6)
1403         ak1 = 7.5e-11*(t(l)/300.)**(0.0)
1404c
1405         rate = (ak0*dens(l))/(1. + ak0*dens(l)/ak1)
1406         xpo = 1./(1. + alog10((ak0*dens(l))/ak1)**2)
1407         c011(l) = rate*0.6**xpo
1408c
1409ccc      c012: o + h2o2 -> oh + ho2
1410c
1411c        jpl 2003
1412c
1413         c012(l) = 1.4e-12*exp(-2000./t(l))
1414c
1415ccc      c013: oh + oh -> h2o + o
1416c
1417c        jpl 2003
1418c
1419         c013(l) = 4.2e-12*exp(-240./t(l))
1420c
1421ccc      c014: oh + o3 -> ho2 + o2
1422c
1423c        jpl 2003
1424c
1425         c014(l) = 1.7e-12*exp(-940./t(l))
1426c
1427c        jpl 2000
1428c
1429c        c014(l) = 1.5e-12*exp(-880./t(l))
1430c
1431c        nair et al., 1994 (jpl 1997)
1432c
1433c        c014(l) = 1.6e-12*exp(-940./t(l))
1434c
1435ccc      c015: ho2 + o3 -> oh + o2 + o2
1436c
1437c        jpl 2003
1438c
1439         c015(l) = 1.0e-14*exp(-490./t(l))
1440c
1441c        jpl 2000
1442c
1443c        c015(l) = 2.0e-14*exp(-680./t(l))
1444c
1445c        nair et al., 1994 (jpl 1997)
1446c
1447c        c015(l) = 1.1e-14*exp(-500./t(l))
1448c
1449ccc      c016: ho2 + ho2 + co2 -> h2o2 + o2 + co2
1450c
1451c        jpl 2003
1452c
1453         c016(l) = 2.5*1.7e-33
1454     $              *exp(1000./t(l))*dens(l)
1455c
1456ccc      c017: oh + oh + co2 -> h2o2 + co2
1457c
1458c        jpl 2003
1459c
1460         ak0 = 2.5*6.9e-31*(t(l)/300.)**(-1.0)
1461         ak1 = 2.6e-11*(t(l)/300.)**(0.0)
1462c
1463c        jpl 1997
1464c
1465c        ak0 = 2.5*6.2e-31*(t(l)/300.)**(-1.0)
1466c        ak1 = 2.6e-11*(t(l)/300.)**(0.0)
1467c
1468c        nair et al., 1994
1469c
1470c        ak0 = 2.5*7.1e-31*(t(l)/300.)**(-0.8)
1471c        ak1 = 1.5e-11*(t(l)/300.)**(0.0)
1472c
1473         rate = (ak0*dens(l))/(1. + ak0*dens(l)/ak1)
1474         xpo = 1./(1. + alog10((ak0*dens(l))/ak1)**2)
1475         c017(l) = rate*0.6**xpo
1476c
1477ccc      c018: h + h + co2 -> h2 + co2
1478c
1479c        baulch et al., 1992
1480c
1481         c018(l) = 2.5*8.85e-33*(t(l)/298.)**(-0.6)*dens(l)
1482c
1483cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1484c        nitrogen compounds
1485cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1486c
1487ccc      d001: no2 + o -> no + o2
1488c
1489c        jpl 2003
1490c
1491         d001(l) = 5.6e-12*exp(180./t(l))
1492c
1493ccc      d002: no + o3 -> no2 + o2
1494c
1495c        jpl 2003
1496c
1497         d002(l) = 3.0e-12*exp(-1500./t(l))
1498c
1499ccc      d003: no + ho2 -> no2 + oh
1500c
1501c        jpl 2003
1502c
1503         d003(l) = 3.5e-12*exp(250./t(l))
1504c
1505cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1506c        carbon compounds
1507cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1508c
1509ccc      e001: oh + co -> co2 + h
1510c
1511c        jpl 2003
1512c
1513c        e001(l) = 1.5e-13*(1 + 0.6*press(l)/1013.)
1514c
1515c        mccabe et al., grl, 28, 3135, 2001
1516c
1517         e001(l) = 1.57e-13 + 3.54e-33*dens(l)
1518c
1519ccc      e002: o + co + m -> co2 + m
1520c
1521c        tsang and hampson, 1986.
1522c
1523         e002(l) = 2.5*6.5e-33*exp(-2184./t(l))*dens(l)
1524c
1525c        baulch et al., butterworths, 1976.
1526c
1527c        e002(l) = 1.6e-32*exp(-2184./t(l))*dens(l)
1528c
1529cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1530c        heterogenous chemistry
1531cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
1532c
1533ccc      h001: ho2 + ice surface -> (ho2)ads
1534ccc           (ho2)ads + oh -> h2o + o2
1535ccc      total : oh + ho2 -> h2o + o2, first order
1536ccc         K = h001*[ho2],
1537ccc      h001 = (S*v*gamma1)/4. (s-1), v=100*sqrt(3kT/m_ho2) (cm s-1)
1538c
1539c        cooper and abbatt, 1996
1540c
1541         h001(l) = icesurf(l) * 25.*sqrt(3*8.31*t(l)/33.) * 0.025
1542c         h001(l) = 0.0e0
1543c
1544ccc      h002: oh + ice surface -> (oh)ads
1545ccc           (oh)ads + ho2 -> h2o + o2
1546ccc      total : oh + ho2 -> h2o + o2, first order
1547ccc         K = h002*[oh],
1548ccc      h002 = (S*v*gamma2)/4. (s-1), v=100*sqrt(3kT/m_oh) (cm s-1)
1549c
1550c        cooper and abbatt, 1996
1551c
1552         h002(l) = icesurf(l) * 25.*sqrt(3*8.31*t(l)/17.) * 0.03
1553c         h002(l) = 0.0e0
1554c
1555c        write(*,*) "level ",l," / icesurf=",icesurf(l),
1556c     $       " / 0.25vg1=",(25.*sqrt(3*8.31*t(l)/33.) * 0.025),
1557c     $       " / 0.25vg2=",(25.*sqrt(3*8.31*t(l)/17.) * 0.03)
1558c
1559      end do
1560c
1561      return
1562      end
Note: See TracBrowser for help on using the repository browser.