source: trunk/LMDZ.MARS/libf/aeronomars/calchim.F90 @ 2041

Last change on this file since 2041 was 2031, checked in by flefevre, 7 years ago

Optimisation photolyse on-line

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