source: trunk/LMDZ.MARS/libf/aeronomars/calchim_mod.F90 @ 2170

Last change on this file since 2170 was 2170, checked in by flefevre, 5 years ago
  • la chimie des ions est désormais optionnelle (FGG)
  • regroupement des options de la photochimie dans calchim_mod (FL)
  • mise à jour cinétique CO + OH, inactive pour le moment (FL)
File size: 28.6 KB
RevLine 
[2162]1MODULE calchim_mod
2
3  IMPLICIT NONE
4
[2164]5  INTEGER, SAVE :: ichemistry ! compute chemistry every ichemistry physics step
[2162]6  REAL,SAVE,ALLOCATABLE :: zdqchim(:,:,:) ! Tendancy on pq due to photochemistry
7  REAL,SAVE,ALLOCATABLE :: zdqschim(:,:) ! Tendancy on qsurf due to photochemistry
8
9  CONTAINS
10
[2158]11      subroutine calchim(ngrid,nlayer,nq,                           &
[1495]12                         ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0,         &
13                         zzlev,zzlay,zday,pq,pdq,dqchim,dqschim,dqcloud,    &
[2031]14                         dqscloud,tau,co2ice,                               &
[1495]15                         pu,pdu,pv,pdv,surfdust,surfice)
16
17      use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_o1d, igcm_o2, &
18                            igcm_o3, igcm_h, igcm_h2, igcm_oh, igcm_ho2,  &
19                            igcm_h2o2, igcm_ch4, igcm_n2, igcm_h2o_vap,   &
20                            igcm_no, igcm_n, igcm_no2, igcm_n2d,          &
21                            igcm_o2plus, igcm_co2plus, igcm_oplus,        &
22                            igcm_coplus, igcm_cplus, igcm_nplus,          &
23                            igcm_noplus, igcm_n2plus, igcm_hplus,         &
24                            igcm_hco2plus, igcm_elec, mmol
25
26      use conc_mod, only: mmean ! mean molecular mass of the atmosphere
[2162]27      use comcstfi_h, only: pi
[2170]28      use photolysis_mod, only: init_photolysis, nphot
[2158]29      use iono_h, only: temp_elect
[1495]30
31      implicit none
32
33!=======================================================================
34!
35!   subject:
36!   --------
37!
38!  Prepare the call for the photochemical module, and send back the
39!  tendencies from photochemistry in the chemical species mass mixing ratios
40!
41!   Arguments:
42!   ----------
43!
44!  Input:
45!
[2031]46!    ptimestep              timestep (s)
[1495]47!    pplay(ngrid,nlayer)    Pressure at the middle of the layers (Pa)
48!    pplev(ngrid,nlayer+1)  Intermediate pressure levels (Pa)
49!    pt(ngrid,nlayer)       Temperature (K)
50!    pdt(ngrid,nlayer)      Temperature tendency (K)
51!    pu(ngrid,nlayer)       u component of the wind (ms-1)
52!    pdu(ngrid,nlayer)      u component tendency (K)
53!    pv(ngrid,nlayer)       v component of the wind (ms-1)
54!    pdv(ngrid,nlayer)      v component tendency (K)
[2031]55!    dist_sol               distance of the sun (AU)
56!    mu0(ngrid)             cos of solar zenith angle (=1 when sun at zenith)
57!    pq(ngrid,nlayer,nq)    advected fields, ie chemical species here
58!    pdq(ngrid,nlayer,nq)   previous tendencies on pq
59!    tau(ngrid)             dust optical depth
60!    co2ice(ngrid)          co2 ice surface layer (kg.m-2)
[1495]61!    surfdust(ngrid,nlayer) dust surface area (m2/m3)
62!    surfice(ngrid,nlayer)  ice surface area (m2/m3)
63!
64!  Output:
65!
[2031]66!    dqchim(ngrid,nlayer,nq) tendencies on pq due to chemistry
67!    dqschim(ngrid,nq)       tendencies on qsurf
[1495]68!
69!=======================================================================
70
[2162]71include "callkeys.h"
[1495]72
73!     input:
74
75      integer,intent(in) :: ngrid    ! number of atmospheric columns
76      integer,intent(in) :: nlayer   ! number of atmospheric layers
77      integer,intent(in) :: nq       ! number of tracers
78      real :: ptimestep
79      real :: pplay(ngrid,nlayer)    ! pressure at the middle of the layers
80      real :: zzlay(ngrid,nlayer)    ! pressure at the middle of the layers
81      real :: pplev(ngrid,nlayer+1)  ! intermediate pressure levels
82      real :: zzlev(ngrid,nlayer+1)  ! altitude at layer boundaries
83      real :: pt(ngrid,nlayer)       ! temperature
84      real :: pdt(ngrid,nlayer)      ! temperature tendency
85      real :: pu(ngrid,nlayer)       ! u component of the wind (m.s-1)
86      real :: pdu(ngrid,nlayer)      ! u component tendency
87      real :: pv(ngrid,nlayer)       ! v component of the wind (m.s-1)
88      real :: pdv(ngrid,nlayer)      ! v component tendency
89      real :: dist_sol               ! distance of the sun (AU)
90      real :: mu0(ngrid)             ! cos of solar zenith angle (=1 when sun at zenith)
91      real :: pq(ngrid,nlayer,nq)    ! tracers mass mixing ratio
92      real :: pdq(ngrid,nlayer,nq)   ! previous tendencies
93      real :: zday                   ! date (time since Ls=0, in martian days)
[2031]94      real :: tau(ngrid)             ! dust optical depth
[1495]95      real :: co2ice(ngrid)          ! co2 ice surface layer (kg.m-2)
96      real :: surfdust(ngrid,nlayer) ! dust surface area (m2/m3)
[2031]97      real :: surfice(ngrid,nlayer)  ! ice surface area (m2/m3)
[1495]98
99!     output:
100
101      real :: dqchim(ngrid,nlayer,nq)   ! tendencies on pq due to chemistry
102      real :: dqschim(ngrid,nq)         ! tendencies on qsurf
103      real :: dqcloud(ngrid,nlayer,nq)  ! tendencies on pq due to condensation
104      real :: dqscloud(ngrid,nq)        ! tendencies on qsurf
105
106!     local variables:
107
108      integer,save :: nbq                   ! number of tracers used in the chemistry
109      integer,allocatable,save :: niq(:)    ! array storing the indexes of the tracers
110      integer :: iloc(1)                    ! index of major species
111      integer :: ig,l,i,iq,iqmax
112      integer :: foundswitch, lswitch
113      integer,save :: chemthermod
114
115      integer,save :: i_co2  = 0
116      integer,save :: i_co   = 0
117      integer,save :: i_o    = 0
118      integer,save :: i_o1d  = 0
119      integer,save :: i_o2   = 0
120      integer,save :: i_o3   = 0
121      integer,save :: i_h    = 0
122      integer,save :: i_h2   = 0
123      integer,save :: i_oh   = 0
124      integer,save :: i_ho2  = 0
125      integer,save :: i_h2o2 = 0
126      integer,save :: i_ch4  = 0
127      integer,save :: i_n2   = 0
128      integer,save :: i_h2o  = 0
129      integer,save :: i_n    = 0
130      integer,save :: i_no   = 0
131      integer,save :: i_no2  = 0
132      integer,save :: i_n2d  = 0
133      integer,save :: i_co2plus=0
134      integer,save :: i_oplus=0
135      integer,save :: i_o2plus=0
136      integer,save :: i_coplus=0
137      integer,save :: i_cplus=0
138      integer,save :: i_nplus=0
139      integer,save :: i_noplus=0
140      integer,save :: i_n2plus=0
141      integer,save :: i_hplus=0
142      integer,save :: i_hco2plus=0
143      integer,save :: i_elec=0
144
145      integer :: ig_vl1
146
[2170]147      integer :: nb_reaction_3_max     ! number of quadratic reactions
148      integer :: nb_reaction_4_max     ! number of bimolecular reactions
149      integer :: nquench               ! number of quenching + heterogeneous reactions
150      integer :: nphotion              ! number of photoionizations
151      integer :: nb_phot_max           ! total number of photolysis+photoionizations+quenching reactions
152
153
[1495]154      real    :: latvl1, lonvl1
155      real    :: zq(ngrid,nlayer,nq)   ! pq+pdq*ptimestep before chemistry
156                                       ! new mole fraction after
157      real    :: zt(ngrid,nlayer)      ! temperature
158      real    :: zu(ngrid,nlayer)      ! u component of the wind
159      real    :: zv(ngrid,nlayer)      ! v component of the wind
[2027]160      real    :: taucol                ! dust optical depth at the surface
[2030]161      real    :: kb                    ! boltzmann constant
[1495]162
[2170]163      logical, save :: firstcall = .true.
164      logical, save :: depos       ! switch for dry deposition
165      logical, save :: ionchem     ! switch for ionospheric chemistry
166      logical, save :: jonline     ! switch for online photodissociation rates or lookup table
167      logical, save :: unichim     ! only one unified chemical scheme at all
168                                   ! layers (default), or upper atmospheric
169                                   ! scheme in the thermosphere
[1495]170
171!     for each column of atmosphere:
172
173      real :: zpress(nlayer)       ! Pressure (mbar)
174      real :: zdens(nlayer)        ! Density  (cm-3)
175      real :: ztemp(nlayer)        ! Temperature (K)
[2158]176      real :: ztelec(nlayer)       ! Electronic temperature (K)
[1495]177      real :: zlocal(nlayer)       ! Altitude (km)
178      real :: zycol(nlayer,nq)     ! Composition (mole fractions)
179      real :: zmmean(nlayer)       ! Mean molecular mass (g.mole-1)
180      real :: szacol               ! Solar zenith angle
181      real :: surfice1d(nlayer)    ! Ice surface area (cm2/cm3)
182      real :: surfdust1d(nlayer)   ! Dust surface area (cm2/cm3)
183      real :: jo3(nlayer)          ! Photodissociation rate O3->O1D (s-1)
[2030]184      real :: jh2o(nlayer)         ! Photodissociation rate H2O->H+OH (s-1)
[2170]185      real :: em_no(nlayer)        ! NO nightglow emission rate
186      real :: em_o2(nlayer)        ! O2 nightglow emission rate     
[1495]187
[2170]188      integer :: iter(nlayer)      ! Number of chemical iterations
189                                   ! within one physical timestep
[1495]190
191!     for output:
192
193      logical :: output             ! to issue calls to writediagfi and stats
194      real :: jo3_3d(ngrid,nlayer)  ! Photodissociation rate O3->O1D (s-1)
[2030]195      real :: jh2o_3d(ngrid,nlayer)  ! Photodissociation rate H2O->H+OH (s-1)
[1888]196      real :: emission_no(ngrid,nlayer) !NO emission rate
197      real :: emission_o2(ngrid,nlayer) !O2 emission rate
[1495]198      real :: iter_3d(ngrid,nlayer) ! Number of chemical iterations
199                                    !  within one physical timestep
200
[2170]201!=======================================================================
202!     main dashboard for the chemistry
203!=======================================================================
[2158]204
[2170]205      unichim = .true.     ! true : unified chemistry ! false : separate models in lower and upper atmosphere
206      jonline = .true.     ! true : on-line calculation of photodissociation rates ! false : lookup table
207      ionchem = .false.    ! switch for ionospheric chemistry
208      depos   = .false.    ! switch for dry deposition
209      output  = .false.    ! true : triggers writing of specific outputs (photolysis rates, emission rates, etc)
210
[1495]211!=======================================================================
212!     initialization of the chemistry (first call only)
213!=======================================================================
214
215      if (firstcall) then
216
217         if (photochem) then
[2044]218            if (jonline) then
219               print*,'calchim: Read UV absorption cross-sections'
220               call init_photolysis     ! on-line photolysis
221            else
222               print*,'calchim: Read photolysis lookup table'
223               call read_phototable     ! off-line photolysis
224            end if
[1495]225         end if
[2158]226
227         if(.not.unichim) then
228            !Read reaction rates from external file if the upper atmospheric
229            !chemistry is called
230            call chemthermos_readini
231         endif
232
[1495]233         ! find index of chemical tracers to use
234         allocate(niq(nq))
235         ! Listed here are all tracers that can go into photochemistry
236         nbq = 0 ! to count number of tracers
[2158]237         ! Species ALWAYS present if photochem=.T.
[1495]238         i_co2 = igcm_co2
239         if (i_co2 == 0) then
240            write(*,*) "calchim: Error; no CO2 tracer !!!"
241            stop
242         else
243            nbq = nbq + 1
244            niq(nbq) = i_co2
245         end if
246         i_co = igcm_co
247         if (i_co == 0) then
248            write(*,*) "calchim: Error; no CO tracer !!!"
249            stop
250         else
251            nbq = nbq + 1
252            niq(nbq) = i_co
253         end if
254         i_o = igcm_o
255         if (i_o == 0) then
256            write(*,*) "calchim: Error; no O tracer !!!"
257            stop
258         else
259            nbq = nbq + 1
260            niq(nbq) = i_o
261         end if
262         i_o1d = igcm_o1d
263         if (i_o1d == 0) then
264            write(*,*) "calchim: Error; no O1D tracer !!!"
265            stop
266         else
267            nbq = nbq + 1
268            niq(nbq) = i_o1d
269         end if
270         i_o2 = igcm_o2
271         if (i_o2 == 0) then
272            write(*,*) "calchim: Error; no O2 tracer !!!"
273            stop
274         else
275            nbq = nbq + 1
276            niq(nbq) = i_o2
277         end if
278         i_o3 = igcm_o3
279         if (i_o3 == 0) then
280            write(*,*) "calchim: Error; no O3 tracer !!!"
281            stop
282         else
283            nbq = nbq + 1
284            niq(nbq) = i_o3
285         end if
286         i_h = igcm_h
287         if (i_h == 0) then
288            write(*,*) "calchim: Error; no H tracer !!!"
289            stop
290         else
291            nbq = nbq + 1
292            niq(nbq) = i_h
293         end if
294         i_h2 = igcm_h2
295         if (i_h2 == 0) then
296            write(*,*) "calchim: Error; no H2 tracer !!!"
297            stop
298         else
299            nbq = nbq + 1
300            niq(nbq) = i_h2
301         end if
302         i_oh = igcm_oh
303         if (i_oh == 0) then
304            write(*,*) "calchim: Error; no OH tracer !!!"
305            stop
306         else
307            nbq = nbq + 1
308            niq(nbq) = i_oh
309         end if
310         i_ho2 = igcm_ho2
311         if (i_ho2 == 0) then
312            write(*,*) "calchim: Error; no HO2 tracer !!!"
313            stop
314         else
315            nbq = nbq + 1
316            niq(nbq) = i_ho2
317         end if
318         i_h2o2 = igcm_h2o2
319         if (i_h2o2 == 0) then
320            write(*,*) "calchim: Error; no H2O2 tracer !!!"
321            stop
322         else
323            nbq = nbq + 1
324            niq(nbq) = i_h2o2
325         end if
326         i_ch4 = igcm_ch4
327         if (i_ch4 == 0) then
328            write(*,*) "calchim: Error; no CH4 tracer !!!"
329            write(*,*) "CH4 will be ignored in the chemistry"
330         else
331            nbq = nbq + 1
332            niq(nbq) = i_ch4
333         end if
334         i_n2 = igcm_n2
335         if (i_n2 == 0) then
336            write(*,*) "calchim: Error; no N2 tracer !!!"
337            stop
338         else
339            nbq = nbq + 1
340            niq(nbq) = i_n2
341         end if
342         i_n = igcm_n
343         if (i_n == 0) then
344            if (photochem) then
345               write(*,*) "calchim: Error; no N tracer !!!"
346               stop
347            end if
348         else
349            nbq = nbq + 1
350            niq(nbq) = i_n
351         end if
352         i_n2d = igcm_n2d
353         if (i_n2d == 0) then
354            if (photochem) then
355               write(*,*) "calchim: Error; no N2D tracer !!!"
356               stop
357            end if
358         else
359            nbq = nbq + 1
360            niq(nbq) = i_n2d
361         end if
362         i_no = igcm_no
363         if (i_no == 0) then
364            if (photochem) then
365               write(*,*) "calchim: Error; no NO tracer !!!"
366               stop
367            end if
368         else
369            nbq = nbq + 1
370            niq(nbq) = i_no
371         end if
372         i_no2 = igcm_no2
373         if (i_no2 == 0) then
374            if (photochem) then
375               write(*,*) "calchim: Error; no NO2 tracer !!!"
376               stop
377            end if
378         else
379            nbq = nbq + 1
380            niq(nbq) = i_no2
381         end if
382         i_h2o = igcm_h2o_vap
383         if (i_h2o == 0) then
384            write(*,*) "calchim: Error; no water vapor tracer !!!"
385            stop
386         else
387            nbq = nbq + 1
388            niq(nbq) = i_h2o
389         end if
[2158]390         i_o2plus = igcm_o2plus
391         if (i_o2plus == 0) then
[2170]392            write(*,*) "calchim: no O2+ tracer "
393            write(*,*) "Only neutral chemistry"
[2158]394         else
395            nbq = nbq + 1
396            niq(nbq) = i_o2plus
[2170]397            ionchem = .true.
398            write(*,*) "calchim: O2+ tracer found in traceur.def"
399            write(*,*) "Ion chemistry included"
[2158]400         end if
401         i_co2plus = igcm_co2plus
[2170]402         if(ionchem) then
403            if (i_co2plus == 0) then
404               write(*,*) "calchim: Error, no CO2+ tracer !!!"
405               write(*,*) "CO2+ is needed if O2+ is in traceur.def"
406               stop
407            else
408               nbq = nbq + 1
409               niq(nbq) = i_co2plus
410            end if
[2158]411         else
[2170]412            if (i_co2plus /= 0) then
413               write(*,*) "calchim: Error: CO2+ is present, but O2+ is not!!!"
414               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
415               stop
416            endif
417         endif
[2158]418         i_oplus=igcm_oplus
[2170]419         if(ionchem) then
420            if (i_oplus == 0) then
421               write(*,*) "calchim: Error, no O+ tracer !!!"
422               write(*,*) "O+ is needed if O2+ is in traceur.def"
423               stop
424            else
425               nbq = nbq + 1
426               niq(nbq) = i_oplus
427            end if
[2158]428         else
[2170]429            if (i_oplus /= 0) then
430               write(*,*) "calchim: Error: O+ is present, but O2+ is not!!!"
431               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
432               stop
433            endif
434         endif
[2158]435         i_noplus=igcm_noplus
[2170]436         if(ionchem) then
437            if (i_noplus == 0) then
438               write(*,*) "calchim: Error, no NO+ tracer !!!"
439               write(*,*) "NO+ is needed if O2+ is in traceur.def"
440               stop
441            else
442               nbq = nbq + 1
443               niq(nbq) = i_noplus
444            end if
[2158]445         else
[2170]446            if (i_noplus /= 0) then
447               write(*,*) "calchim: Error: NO+ is present, but O2+ is not!!!"
448               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
449            endif
450         endif
[2158]451         i_coplus=igcm_coplus
[2170]452         if(ionchem) then
453            if (i_coplus == 0) then
454               write(*,*) "calchim: Error, no CO+ tracer !!!"
455               write(*,*) "CO+ is needed if O2+ is in traceur.def"
456               stop               
457            else
458               nbq = nbq + 1
459               niq(nbq) = i_coplus
460            end if
[2158]461         else
[2170]462            if (i_coplus /= 0) then
463               write(*,*) "calchim: Error: CO+ is present, but O2+ is not!!!"
464               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
465            endif
466         endif
[2158]467         i_cplus=igcm_cplus
[2170]468         if(ionchem) then
469            if (i_cplus == 0) then
470               write(*,*) "calchim: Error, no C+ tracer !!!"
471               write(*,*) "C+ is needed if O2+ is in traceur.def"
472               stop
473            else
474               nbq = nbq + 1
475               niq(nbq) = i_cplus
476            end if
[2158]477         else
[2170]478            if (i_cplus /= 0) then
479               write(*,*) "calchim: Error: C+ is present, but O2+ is not!!!"
480               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
481            endif
482         endif
[2158]483         i_n2plus=igcm_n2plus
[2170]484         if(ionchem) then
485            if (i_n2plus == 0) then
486               write(*,*) "calchim: Error, no N2+ tracer !!!"
487               write(*,*) "N2+ is needed if O2+ is in traceur.def"
488               stop
489            else
490               nbq = nbq + 1
491               niq(nbq) = i_n2plus
492            end if
[2158]493         else
[2170]494            if (i_n2plus /= 0) then
495               write(*,*) "calchim: Error: N2+ is present, but O2+ is not!!!"
496               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
497            endif
498         endif
[2158]499         i_nplus=igcm_nplus
[2170]500         if(ionchem) then
501            if (i_nplus == 0) then
502               write(*,*) "calchim: Error, no N+ tracer !!!"
503               write(*,*) "N+ is needed if O2+ is in traceur.def"
504               stop
505            else
506               nbq = nbq + 1
507               niq(nbq) = i_nplus
508            end if
[2158]509         else
[2170]510            if (i_nplus /= 0) then
511               write(*,*) "calchim: Error: N+ is present, but O2+ is not!!!"
512               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
513            endif
514         endif
[2158]515         i_hplus=igcm_hplus
[2170]516         if(ionchem) then
517            if (i_hplus == 0) then
518               write(*,*) "calchim: Error, no H+ tracer !!!"
519               write(*,*) "H+ is needed if O2+ is in traceur.def"
520               stop
521            else
522               nbq = nbq + 1
523               niq(nbq) = i_hplus
524            end if
[2158]525         else
[2170]526            if (i_hplus /= 0) then
527               write(*,*) "calchim: Error: H+ is present, but O2+ is not!!!"
528               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
529            endif
530         endif
[2158]531         i_hco2plus=igcm_hco2plus
[2170]532         if(ionchem) then
533            if (i_hco2plus == 0) then
534               write(*,*) "calchim: Error, no HCO2+ tracer !!!"
535               write(*,*) "HCO2+ is needed if O2+ is in traceur.def"
536               stop
537            else
538               nbq = nbq + 1
539               niq(nbq) = i_hco2plus
540            end if
[2158]541         else
[2170]542            if (i_hco2plus /= 0) then
543               write(*,*) "calchim: Error: HCO2+ is present, but O2+ is not!!!"
544               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
545            endif
546         endif
[2158]547         i_elec = igcm_elec
[2170]548         if(ionchem) then
549            if (i_elec == 0) then
550               write(*,*) "calchim: Error, no e- tracer !!!"
551               write(*,*) "e- is needed if O2+ is in traceur.def"
552               stop
553            else
554               nbq = nbq + 1
555               niq(nbq) = i_elec
556            end if
[2158]557         else
[2170]558            if (i_elec /= 0) then
559               write(*,*) "calchim: Error: e- is present, but O2+ is not!!!"
560               write(*,*) "Both must be in traceur.def if ion chemistry wanted"
561            endif
562         endif
[1495]563         write(*,*) 'calchim: found nbq    = ',nbq,' tracers'
564               
565         firstcall = .false.
566      end if ! if (firstcall)
567
568! Initializations
569
570      zycol(:,:)    = 0.
571      dqchim(:,:,:) = 0.
572      dqschim(:,:)  = 0.
573
[2030]574      kb = 1.3806e-23
575
[1495]576!     latvl1= 22.27
577!     lonvl1= -47.94
[2029]578!     ig_vl1= 1+ int( (1.5-(latvl1-90.)*48./180.)  -2 )*64. +    &
579!             int(1.5+(lonvl1+180)*64./360.)
[1495]580
581!=======================================================================
582!     loop over grid
583!=======================================================================
584
585      do ig = 1,ngrid
586         
587         foundswitch = 0
588         do l = 1,nlayer
589            do i = 1,nbq
590               iq = niq(i) ! get tracer index
591               zq(ig,l,iq) = pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep
592               zycol(l,iq) = zq(ig,l,iq)*mmean(ig,l)/mmol(iq)
593            end do
594            zt(ig,l)  = pt(ig,l) + pdt(ig,l)*ptimestep
595            zu(ig,l)  = pu(ig,l) + pdu(ig,l)*ptimestep
596            zv(ig,l)  = pv(ig,l) + pdv(ig,l)*ptimestep
597            zpress(l) = pplay(ig,l)/100.
598            ztemp(l)  = zt(ig,l)
599            zdens(l)  = zpress(l)/(kb*1.e4*ztemp(l))
600            zlocal(l) = zzlay(ig,l)/1000.
601            zmmean(l) = mmean(ig,l)
[2158]602            ztelec(l) = temp_elect(zlocal(l),ztemp(l),1)
603            !Electronic temperature. Index 1 -> Viking; Index 2-> MAVEN
[1495]604
605!           surfdust1d and surfice1d: conversion from m2/m3 to cm2/cm3
606
607            surfdust1d(l) = surfdust(ig,l)*1.e-2
608            surfice1d(l)  = surfice(ig,l)*1.e-2
609
[2158]610!           search for switch index between regions 
[2170]611
612            if (unichim) then
613               lswitch = nlayer + 1
[2158]614            else
[1495]615               if (foundswitch == 0 .and. pplay(ig,l) < 1.e-1) then
616                  lswitch = l
617                  foundswitch = 1
618               end if
[2158]619            endif
[1495]620         end do ! of do l=1,nlayer
621
622         szacol = acos(mu0(ig))*180./pi
[2031]623         taucol = tau(ig)
[1495]624
625!=======================================================================
626!     call chemical subroutines
627!=======================================================================
628
629         if (photochem) then
[2170]630            ! set number of reactions, depending on ion chemistry or not
631            if (ionchem) then
632               nb_reaction_4_max = 60   ! set number of bimolecular reactions
633               nb_reaction_3_max = 6    ! set number of quadratic reactions
634               nquench           = 9    ! set number of quenching + heterogeneous reactions
635               nphotion          = 18   ! set number of photoionizations
636            else
637               nb_reaction_4_max = 31   ! set number of bimolecular reactions
638               nb_reaction_3_max = 6    ! set number of quadratic reactions
639               nquench           = 9    ! set number of quenching + heterogeneous reactions
640               nphotion          = 0    ! set number of photoionizations
641            end if
642
643!        nb_phot_max is the total number of processes that are treated
644!        numerically as a photolysis:
645
646            nb_phot_max = nphot + nphotion + nquench
647
648!        call main photochemistry routine
649
650            call photochemistry(nlayer,nq,nbq,ionchem,nb_reaction_3_max,  &
651                           nb_reaction_4_max, nb_phot_max, nphotion,      &
652                           jonline, ig,lswitch,zycol,szacol,ptimestep,    &
653                           zpress,zlocal,ztemp,ztelec,zdens,zmmean,       &
654                           dist_sol,zday,surfdust1d,surfice1d,            &
[2158]655                           jo3,jh2o,taucol,iter)
[1495]656
657!        ozone photolysis, for output
658
659            do l = 1,nlayer
660               jo3_3d(ig,l) = jo3(l)
[2030]661               jh2o_3d(ig,l) = jh2o(l)
[1495]662               iter_3d(ig,l) = iter(l)
663            end do
[2170]664
[1495]665!        condensation of h2o2
666
667            call perosat(ngrid, nlayer, nq,                        &
668                         ig,ptimestep,pplev,pplay,                 &
669                         ztemp,zycol,dqcloud,dqscloud)
670
[2170]671!        case of separate photochemical model in the thermosphere
672
673            if (.not.unichim) then
674               chemthermod = 3   !C/O/H/N/ions
[2158]675               call chemthermos(ig,nlayer,lswitch,chemthermod,zycol,ztemp,&
676                    zdens,zpress,zlocal,szacol,ptimestep,zday,&
677                    em_no,em_o2)
[2170]678               do l = 1,nlayer
679                  emission_no(ig,l) = em_no(l)
680                  emission_o2(ig,l) = em_o2(l)
681               end do
[2158]682            end if
[1495]683
[2170]684         end if  ! photochem
685
[1495]686!        dry deposition
687
688         if (depos) then
689            call deposition(ngrid, nlayer, nq,                      &
690                            ig, ig_vl1, pplay, pplev, zzlay, zzlev, &
691                            zu, zv, zt, zycol, ptimestep, co2ice)
692         end if
693
694!=======================================================================
695!     tendencies
696!=======================================================================
697
698!     index of the most abundant species at each level
699
700!         major(:) = maxloc(zycol, dim = 2)
701
702!     tendency for the most abundant species = - sum of others
703
704         do l = 1,nlayer
705            iloc=maxloc(zycol(l,:))
706            iqmax=iloc(1)
707            do i = 1,nbq
708               iq = niq(i) ! get tracer index
709               if (iq /= iqmax) then
710                  dqchim(ig,l,iq) = (zycol(l,iq)*mmol(iq)/mmean(ig,l)  &
711                                   - zq(ig,l,iq))/ptimestep
712                  dqchim(ig,l,iqmax) = dqchim(ig,l,iqmax)              &
713                                     - dqchim(ig,l,iq)
714               end if
715            end do
716         end do ! of do l = 1,nlayer
717
718!=======================================================================
719!     end of loop over grid
720!=======================================================================
721
722      end do ! of do ig=1,ngrid
723
724!=======================================================================
725!     write outputs
726!=======================================================================
727
728! value of parameter 'output' to trigger writting of outputs
729! is set above at the declaration of the variable.
730
731      if (photochem .and. output) then
732         if (ngrid > 1) then
733            call writediagfi(ngrid,'jo3','j o3->o1d',    &
734                             's-1',3,jo3_3d(1,1))
[2030]735            call writediagfi(ngrid,'jh2o','jh2o',    &
736                             's-1',3,jh2o_3d(1,1))
[1495]737            call writediagfi(ngrid,'iter','iterations',  &
738                             ' ',3,iter_3d(1,1))
[2170]739
740            if (.not. unichim) then
[2158]741               call writediagfi(ngrid,'emission_no',        &
742                    'NO nightglow emission rate','cm-3 s-1',3,emission_no)
743               call writediagfi(ngrid,'emission_o2',        &
744                    'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
745            endif
[2042]746           
[1495]747            if (callstats) then
748               call wstats(ngrid,'jo3','j o3->o1d',       &
749                           's-1',3,jo3_3d(1,1))
[1888]750               call wstats(ngrid,'emission_no',           &
751                   'NO nightglow emission rate','cm-3 s-1',3,emission_no)
752               call wstats(ngrid,'emission_o2',           &
753                   'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
[1495]754               call wstats(ngrid,'mmean','mean molecular mass',       &
755                           'g.mole-1',3,mmean(1,1))
756            endif
757         end if ! of if (ngrid.gt.1)
758      end if ! of if (output)
759
[2162]760      end subroutine calchim
761
762
763    subroutine ini_calchim_mod(ngrid,nlayer,nq)
764 
765      implicit none
766 
767      integer,intent(in) :: ngrid ! number of atmospheric columns
768      integer,intent(in) :: nlayer ! number of atmospheric layers
769      integer,intent(in) :: nq ! number of tracers
770
771      allocate(zdqchim(ngrid,nlayer,nq))
772      zdqchim(:,:,:)=0
773      allocate(zdqschim(ngrid,nq))
774      zdqschim(:,:)=0
775
776    end subroutine ini_calchim_mod
777
778
779    subroutine end_calchim_mod
780
781      implicit none
782
783      if (allocated(zdqchim))      deallocate(zdqchim)
784      if (allocated(zdqschim))      deallocate(zdqschim)
785
786    end subroutine end_calchim_mod
787
788END MODULE calchim_mod
789
Note: See TracBrowser for help on using the repository browser.