source: trunk/LMDZ.MARS/libf/aeronomars/photochemistry.F @ 1119

Last change on this file since 1119 was 1036, checked in by emillour, 11 years ago

Mars GCM: (a first step towards using parallel dynamics)

  • IMPORTANT CHANGE: Implemented dynamic tracers. It is no longer necessary to compile the model with the '-t #' option, number of tracers is simply read from tracer.def file (as before). Adapted makegcm_* scripts (and co.) accordingly. Technical aspects of the switch to dynamic tracers are:
    • advtrac.h (in dyn3d) removed and replaced by module infotrac.F
    • tracer.h (in phymars) removed and replaced by module tracer_mod.F90 (which contains nqmx, the number of tracers, etc. and can be used anywhere in the physics).
  • Included some side cleanups: removed unused files (in dyn3d) anldoppler2.F, anl_mcdstats.F and anl_stats-diag.F, and all the unecessary dimensions.* files in grid/dimension.
  • Checked that changes are clean and that GCM yields identical results (in debug mode) to previous svn version.

EM

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