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

Last change on this file since 2137 was 2044, checked in by flefevre, 7 years ago

Photolyse on-line : desormais par defaut

File size: 28.3 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            if (jonline) then
196               print*,'calchim: Read UV absorption cross-sections'
197               call init_photolysis     ! on-line photolysis
198            else
199               print*,'calchim: Read photolysis lookup table'
200               call read_phototable     ! off-line photolysis
201            end if
202         end if
203         ! find index of chemical tracers to use
204         allocate(niq(nq))
205         ! Listed here are all tracers that can go into photochemistry
206         nbq = 0 ! to count number of tracers
207         ! Species ALWAYS present if photochem=.T. or thermochem=.T.
208         i_co2 = igcm_co2
209         if (i_co2 == 0) then
210            write(*,*) "calchim: Error; no CO2 tracer !!!"
211            stop
212         else
213            nbq = nbq + 1
214            niq(nbq) = i_co2
215         end if
216         i_co = igcm_co
217         if (i_co == 0) then
218            write(*,*) "calchim: Error; no CO tracer !!!"
219            stop
220         else
221            nbq = nbq + 1
222            niq(nbq) = i_co
223         end if
224         i_o = igcm_o
225         if (i_o == 0) then
226            write(*,*) "calchim: Error; no O tracer !!!"
227            stop
228         else
229            nbq = nbq + 1
230            niq(nbq) = i_o
231         end if
232         i_o1d = igcm_o1d
233         if (i_o1d == 0) then
234            write(*,*) "calchim: Error; no O1D tracer !!!"
235            stop
236         else
237            nbq = nbq + 1
238            niq(nbq) = i_o1d
239         end if
240         i_o2 = igcm_o2
241         if (i_o2 == 0) then
242            write(*,*) "calchim: Error; no O2 tracer !!!"
243            stop
244         else
245            nbq = nbq + 1
246            niq(nbq) = i_o2
247         end if
248         i_o3 = igcm_o3
249         if (i_o3 == 0) then
250            write(*,*) "calchim: Error; no O3 tracer !!!"
251            stop
252         else
253            nbq = nbq + 1
254            niq(nbq) = i_o3
255         end if
256         i_h = igcm_h
257         if (i_h == 0) then
258            write(*,*) "calchim: Error; no H tracer !!!"
259            stop
260         else
261            nbq = nbq + 1
262            niq(nbq) = i_h
263         end if
264         i_h2 = igcm_h2
265         if (i_h2 == 0) then
266            write(*,*) "calchim: Error; no H2 tracer !!!"
267            stop
268         else
269            nbq = nbq + 1
270            niq(nbq) = i_h2
271         end if
272         i_oh = igcm_oh
273         if (i_oh == 0) then
274            write(*,*) "calchim: Error; no OH tracer !!!"
275            stop
276         else
277            nbq = nbq + 1
278            niq(nbq) = i_oh
279         end if
280         i_ho2 = igcm_ho2
281         if (i_ho2 == 0) then
282            write(*,*) "calchim: Error; no HO2 tracer !!!"
283            stop
284         else
285            nbq = nbq + 1
286            niq(nbq) = i_ho2
287         end if
288         i_h2o2 = igcm_h2o2
289         if (i_h2o2 == 0) then
290            write(*,*) "calchim: Error; no H2O2 tracer !!!"
291            stop
292         else
293            nbq = nbq + 1
294            niq(nbq) = i_h2o2
295         end if
296         i_ch4 = igcm_ch4
297         if (i_ch4 == 0) then
298            write(*,*) "calchim: Error; no CH4 tracer !!!"
299            write(*,*) "CH4 will be ignored in the chemistry"
300         else
301            nbq = nbq + 1
302            niq(nbq) = i_ch4
303         end if
304         i_n2 = igcm_n2
305         if (i_n2 == 0) then
306            write(*,*) "calchim: Error; no N2 tracer !!!"
307            stop
308         else
309            nbq = nbq + 1
310            niq(nbq) = i_n2
311         end if
312         i_n = igcm_n
313         if (i_n == 0) then
314            if (photochem) then
315               write(*,*) "calchim: Error; no N tracer !!!"
316               stop
317            end if
318         else
319            nbq = nbq + 1
320            niq(nbq) = i_n
321         end if
322         i_n2d = igcm_n2d
323         if (i_n2d == 0) then
324            if (photochem) then
325               write(*,*) "calchim: Error; no N2D tracer !!!"
326               stop
327            end if
328         else
329            nbq = nbq + 1
330            niq(nbq) = i_n2d
331         end if
332         i_no = igcm_no
333         if (i_no == 0) then
334            if (photochem) then
335               write(*,*) "calchim: Error; no NO tracer !!!"
336               stop
337            end if
338         else
339            nbq = nbq + 1
340            niq(nbq) = i_no
341         end if
342         i_no2 = igcm_no2
343         if (i_no2 == 0) then
344            if (photochem) then
345               write(*,*) "calchim: Error; no NO2 tracer !!!"
346               stop
347            end if
348         else
349            nbq = nbq + 1
350            niq(nbq) = i_no2
351         end if
352         i_h2o = igcm_h2o_vap
353         if (i_h2o == 0) then
354            write(*,*) "calchim: Error; no water vapor tracer !!!"
355            stop
356         else
357            nbq = nbq + 1
358            niq(nbq) = i_h2o
359         end if
360         !Check tracers needed for thermospheric chemistry
361         if(thermochem) then
362            chemthermod=0  !Default: C/O/H chemistry
363            !Nitrogen chemistry
364            !NO is used to determine if N chemistry is wanted
365            !chemthermod=2 -> N chemistry
366            if (i_no == 0) then
367               write(*,*) "calchim: no NO tracer"
368               write(*,*) "C/O/H themosp chemistry only "
369            else
370               chemthermod=2
371               write(*,*) "calchim: NO in traceur.def"
372               write(*,*) "Nitrogen chemistry included"
373            end if
374            ! N
375            if(chemthermod == 2) then
376               if (i_n == 0) then
377                  write(*,*) "calchim: Error; no N tracer !!!"
378                  write(*,*) "N is needed if NO is in traceur.def"
379                  stop
380               end if
381            ! NO2
382               if (i_no2 == 0) then
383                  write(*,*) "calchim: Error; no NO2 tracer !!!"
384                  write(*,*) "NO2 is needed if NO is in traceur.def"
385                  stop
386               end if
387            ! N(2D)
388               if (i_n2d == 0) then
389                  write(*,*) "calchim: Error; no N2D !!!"
390                  write(*,*) "N2D is needed if NO is in traceur.def"
391                  stop
392               end if
393            endif    !Of if(chemthermod == 2)
394            ! Ions
395            ! O2+ is used to determine if ion chemistry is needed
396            ! chemthermod=3 -> ion chemistry
397            i_o2plus = igcm_o2plus
398            if(chemthermod == 2) then
399               if (i_o2plus == 0) then
400                  write(*,*) "calchim: no O2+ tracer; no ion chemistry"
401               else
402                  nbq = nbq + 1
403                  niq(nbq) = i_o2plus
404                  chemthermod = 3
405                  write(*,*) "calchim: O2+ in traceur.def"
406                  write(*,*) "Ion chemistry included"
407               end if
408            else
409               if (i_o2plus /= 0) then
410                  write(*,*) "calchim: O2+ is present, but NO is not!!!"
411                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
412                  stop
413               endif
414            endif   !Of if(chemthermod == 2)
415            ! CO2+
416            i_co2plus = igcm_co2plus
417            if(chemthermod == 3) then
418               if (i_co2plus == 0) then
419                  write(*,*) "calchim: Error; no CO2+ tracer !!!"
420                  write(*,*) "CO2+ is needed if O2+ is in traceur.def"
421                  stop
422               else
423                  nbq = nbq + 1
424                  niq(nbq) = i_co2plus
425               end if
426            else
427               if (i_co2plus /= 0) then
428                  write(*,*) "calchim: Error: CO2+ is present, but O2+ is not!!!"
429                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
430                  stop
431               endif
432            endif    !Of if(chemthermod == 3)
433            ! O+
434            i_oplus = igcm_oplus
435            if(chemthermod == 3) then
436               if (i_oplus == 0) then
437                  write(*,*) "calchim: Error; no O+ tracer !!!"
438                  write(*,*) "O+ is needed if O2+ is in traceur.def"
439                  stop
440               else
441                  nbq = nbq + 1
442                  niq(nbq) = i_oplus
443               end if
444            else
445               if (i_oplus /= 0) then
446                  write(*,*) "calchim: Error: O+ is present, but O2+ is not!!!"
447                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
448                  stop
449               endif
450            endif   !Of if (chemthermod == 3)
451            ! CO+
452            i_coplus = igcm_coplus
453            if(chemthermod == 3) then
454               if (i_coplus == 0) then
455                  write(*,*) "calchim: Error; no CO+ tracer !!!"
456                  write(*,*) "CO+ is needed if O2+ is in traceur.def"
457                  stop
458               else
459                  nbq = nbq + 1
460                  niq(nbq) = i_coplus
461               end if
462            else
463               if (i_coplus /= 0) then
464                  write(*,*) "calchim: Error: CO+ is present, but O2+ is not!!!"
465                  write(*,*) " Both must be in traceur.def if ionosphere wanted"
466                  stop
467               endif
468            endif   ! Of if (chemthermod == 3)
469            ! C+
470            i_cplus = igcm_cplus
471            if(chemthermod == 3) then
472               if (i_cplus == 0) then
473                  write(*,*) "calchim: Error; no C+ tracer !!!"
474                  write(*,*) "C+ is needed if O2+ is in traceur.def"
475                  stop
476               else
477                  nbq = nbq + 1
478                  niq(nbq) = i_cplus
479               end if
480            else
481               if (i_cplus /= 0) then
482                  write(*,*) "calchim: Error; C+ is present, but O2+ is not!!!"
483                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
484                  stop
485               endif
486            endif   ! Of if (chemthermod == 3)
487            ! N+
488            i_nplus = igcm_nplus
489            if(chemthermod == 3) then
490               if (i_nplus == 0) then
491                  write(*,*) "calchim: Error; no N+ tracer !!!"
492                  write(*,*) "N+ is needed if O2+ is in traceur.def"
493                  stop
494               else
495                  nbq = nbq + 1
496                  niq(nbq) = i_nplus
497               end if
498            else
499               if (i_nplus /= 0) then
500                  write(*,*) "calchim: Error: N+ is present, but O2+ is not!!!"
501                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
502                  stop
503               endif
504            endif   !Of if (chemthermod == 3)
505            ! NO+
506            i_noplus = igcm_noplus
507            if(chemthermod == 3) then
508               if (i_noplus == 0) then
509                  write(*,*) "calchim: Error; no NO+ tracer !!!"
510                  write(*,*) "NO+ is needed if O2+ is in traceur.def"
511                  stop
512               else
513                  nbq = nbq + 1
514                  niq(nbq) = i_noplus
515               end if
516            else
517               if (i_noplus /= 0) then
518                  write(*,*) "calchim: Error: NO+ is present, but O2+ is not!!!"
519                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
520                  stop
521               endif
522            endif   !Of if (chemthermod == 3)
523            ! N2+
524            i_n2plus = igcm_n2plus
525            if (chemthermod == 3) then
526               if (i_n2plus == 0) then
527                  write(*,*) "calchim: Error; no N2+ tracer !!!"
528                  write(*,*) "N2+ is needed if O2+ is in traceur.def"
529                  stop
530               else
531                  nbq = nbq + 1
532                  niq(nbq) = i_n2plus
533               end if
534            else
535               if (i_n2plus /= 0) then
536                  write(*,*) "calchim: Error: N2+ is present, but O2+ is not!!!"
537                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
538                  stop
539               endif
540            endif   !Of if (chemthermod == 3)
541            !H+
542            i_hplus = igcm_hplus
543            if (chemthermod == 3) then
544               if (i_hplus == 0) then
545                  write(*,*) "calchim: Error; no H+ tracer !!!"
546                  write(*,*) "H+ is needed if O2+ is in traceur.def"
547                  stop
548               else
549                  nbq = nbq + 1
550                  niq(nbq) = i_hplus
551               end if
552            else
553               if (i_hplus /= 0) then
554                  write(*,*) "calchim: Error: H+ is present, but O2+ is not!!!"
555                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
556                  stop
557               endif
558            endif   !Of if (chemthermod == 3)
559            ! HCO2+
560            i_hco2plus = igcm_hco2plus
561            if(chemthermod == 3) then
562               if (i_hco2plus == 0) then
563                  write(*,*) "calchim: Error; no HCO2+ tracer !!!"
564                  write(*,*) "HCO2+ is needed if O2+ is in traceur.def"
565                  stop
566               else
567                  nbq = nbq + 1
568                  niq(nbq) = i_hco2plus
569               end if
570            else
571               if (i_hco2plus /= 0) then
572                  write(*,*) "calchim: Error: HCO2+ is present, but O2+ is not!!!"
573                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
574                  stop
575               endif
576            endif    !Of if(chemthermod == 3)
577            !e-
578            i_elec = igcm_elec
579            if(chemthermod == 3) then
580               if (i_elec == 0) then
581                  write(*,*) "calchim: Error; no e- tracer !!!"
582                  write(*,*) "e- is needed if O2+ is in traceur.def"
583                  stop
584               else
585                  nbq = nbq + 1
586                  niq(nbq) = i_elec
587               end if
588            else
589               if(i_elec /= 0) then
590                  write(*,*) "calchim: Error: e- is present, but O2+ is not!!!"
591                  write(*,*) "Both must be in traceur.def if ionosphere wanted"
592                  stop
593               endif
594            endif   !Of if (chemthermod == 3)
595         endif      !Of thermochem
596
597         write(*,*) 'calchim: found nbq    = ',nbq,' tracers'
598               
599         firstcall = .false.
600      end if ! if (firstcall)
601
602! Initializations
603
604      zycol(:,:)    = 0.
605      dqchim(:,:,:) = 0.
606      dqschim(:,:)  = 0.
607
608      kb = 1.3806e-23
609
610!     latvl1= 22.27
611!     lonvl1= -47.94
612!     ig_vl1= 1+ int( (1.5-(latvl1-90.)*48./180.)  -2 )*64. +    &
613!             int(1.5+(lonvl1+180)*64./360.)
614
615!=======================================================================
616!     loop over grid
617!=======================================================================
618
619      do ig = 1,ngrid
620         
621         foundswitch = 0
622         do l = 1,nlayer
623            do i = 1,nbq
624               iq = niq(i) ! get tracer index
625               zq(ig,l,iq) = pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep
626               zycol(l,iq) = zq(ig,l,iq)*mmean(ig,l)/mmol(iq)
627            end do
628            zt(ig,l)  = pt(ig,l) + pdt(ig,l)*ptimestep
629            zu(ig,l)  = pu(ig,l) + pdu(ig,l)*ptimestep
630            zv(ig,l)  = pv(ig,l) + pdv(ig,l)*ptimestep
631            zpress(l) = pplay(ig,l)/100.
632            ztemp(l)  = zt(ig,l)
633            zdens(l)  = zpress(l)/(kb*1.e4*ztemp(l))
634            zlocal(l) = zzlay(ig,l)/1000.
635            zmmean(l) = mmean(ig,l)
636
637!           surfdust1d and surfice1d: conversion from m2/m3 to cm2/cm3
638
639            surfdust1d(l) = surfdust(ig,l)*1.e-2
640            surfice1d(l)  = surfice(ig,l)*1.e-2
641
642!           search for switch index between regions
643
644            if (photochem .and. thermochem) then
645               if (foundswitch == 0 .and. pplay(ig,l) < 1.e-1) then
646                  lswitch = l
647                  foundswitch = 1
648               end if
649            end if
650            if (.not. photochem) then
651               lswitch = 22
652            end if
653            if (.not. thermochem) then
654               lswitch = min(50,nlayer+1)
655            end if
656
657         end do ! of do l=1,nlayer
658
659         szacol = acos(mu0(ig))*180./pi
660         taucol = tau(ig)
661
662!=======================================================================
663!     call chemical subroutines
664!=======================================================================
665
666!        chemistry in lower atmosphere
667
668         if (photochem) then
669            call photochemistry(nlayer,nq,                            &
670                           ig,lswitch,zycol,szacol,ptimestep,         &
671                           zpress,zlocal,ztemp,zdens,zmmean,dist_sol, &
672                           zday,surfdust1d,surfice1d,jo3,jh2o,taucol,iter)
673
674!        ozone photolysis, for output
675
676            do l = 1,nlayer
677               jo3_3d(ig,l) = jo3(l)
678               jh2o_3d(ig,l) = jh2o(l)
679               iter_3d(ig,l) = iter(l)
680            end do
681
682!        condensation of h2o2
683
684            call perosat(ngrid, nlayer, nq,                        &
685                         ig,ptimestep,pplev,pplay,                 &
686                         ztemp,zycol,dqcloud,dqscloud)
687         end if
688
689!        chemistry in upper atmosphere
690       
691         if (thermochem) then
692            call chemthermos(ig,nlayer,lswitch,chemthermod,zycol,ztemp,&
693                             zdens,zpress,zlocal,szacol,ptimestep,zday,&
694                             em_no,em_o2)
695            do l=1,nlayer
696               emission_no(ig,l)=em_no(l)
697               emission_o2(ig,l)=em_o2(l)
698            enddo
699         end if
700
701!        dry deposition
702
703         if (depos) then
704            call deposition(ngrid, nlayer, nq,                      &
705                            ig, ig_vl1, pplay, pplev, zzlay, zzlev, &
706                            zu, zv, zt, zycol, ptimestep, co2ice)
707         end if
708
709!=======================================================================
710!     tendencies
711!=======================================================================
712
713!     index of the most abundant species at each level
714
715!         major(:) = maxloc(zycol, dim = 2)
716
717!     tendency for the most abundant species = - sum of others
718
719         do l = 1,nlayer
720            iloc=maxloc(zycol(l,:))
721            iqmax=iloc(1)
722            do i = 1,nbq
723               iq = niq(i) ! get tracer index
724               if (iq /= iqmax) then
725                  dqchim(ig,l,iq) = (zycol(l,iq)*mmol(iq)/mmean(ig,l)  &
726                                   - zq(ig,l,iq))/ptimestep
727                  dqchim(ig,l,iqmax) = dqchim(ig,l,iqmax)              &
728                                     - dqchim(ig,l,iq)
729               end if
730            end do
731         end do ! of do l = 1,nlayer
732
733!=======================================================================
734!     end of loop over grid
735!=======================================================================
736
737      end do ! of do ig=1,ngrid
738
739!      stop
740!=======================================================================
741!     write outputs
742!=======================================================================
743
744! value of parameter 'output' to trigger writting of outputs
745! is set above at the declaration of the variable.
746
747      if (photochem .and. output) then
748         if (ngrid > 1) then
749            call writediagfi(ngrid,'jo3','j o3->o1d',    &
750                             's-1',3,jo3_3d(1,1))
751            call writediagfi(ngrid,'jh2o','jh2o',    &
752                             's-1',3,jh2o_3d(1,1))
753            call writediagfi(ngrid,'iter','iterations',  &
754                             ' ',3,iter_3d(1,1))
755            call writediagfi(ngrid,'emission_no',        &
756                 'NO nightglow emission rate','cm-3 s-1',3,emission_no)
757            call writediagfi(ngrid,'emission_o2',        &
758                 'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
759           
760            if (callstats) then
761               call wstats(ngrid,'jo3','j o3->o1d',       &
762                           's-1',3,jo3_3d(1,1))
763               call wstats(ngrid,'emission_no',           &
764                   'NO nightglow emission rate','cm-3 s-1',3,emission_no)
765               call wstats(ngrid,'emission_o2',           &
766                   'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
767               call wstats(ngrid,'mmean','mean molecular mass',       &
768                           'g.mole-1',3,mmean(1,1))
769            endif
770         end if ! of if (ngrid.gt.1)
771      end if ! of if (output)
772
773      return
774      end
Note: See TracBrowser for help on using the repository browser.