source: trunk/LMDZ.GENERIC/libf/aeronostd/calchim_asis.F90 @ 1800

Last change on this file since 1800 was 1796, checked in by bclmd, 8 years ago

Adding photochemistry to LMDZ Generic

  • Property svn:executable set to *
File size: 18.8 KB
Line 
1      subroutine calchim_asis(ngrid,nlayer,nq,                              &
2                         ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0,fract,   &
3                         zzlev,zzlay,zday,pq,pdq,dqchim,dqschim,            &
4                         tauref,co2ice,                                     &
5                         pu,pdu,pv,pdv,surfdust,surfice)
6
7      use tracer_h, 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_n, igcm_no, igcm_no2, igcm_n2d,   &
11                            mmol
12
13      use conc_mod, only: mmean ! mean molecular mass of the atmosphere
14      USE comcstfi_mod
15      use callkeys_mod
16
17      implicit none
18
19!=======================================================================
20!
21!   subject:
22!   --------
23!
24!  Prepare the call for the photochemical module, and send back the
25!  tendencies from photochemistry in the chemical species mass mixing ratios
26!
27!   Author:   Sebastien Lebonnois (08/11/2002)
28!   -------
29!    update 12/06/2003 for water ice clouds and compatibility with dust
30!    update 07/2003 for coupling with thermosphere (Monica Angelats-i-Coll)
31!    update 03/05/2005 cosmetic changes (Franck Lefevre)
32!    update sept. 2008 identify tracers by their names (Ehouarn Millour)
33!    update 05/12/2011 synchronize with latest version of chemistry (Franck Lefevre)
34!    update 16/03/2012 optimization (Franck Lefevre)
35!    update 11/12/2013 optimization (Franck Lefevre)
36!
37!   Arguments:
38!   ----------
39!
40!  Input:
41!
42!    ptimestep                  timestep (s)
43!    pplay(ngrid,nlayer)    Pressure at the middle of the layers (Pa)
44!    pplev(ngrid,nlayer+1)  Intermediate pressure levels (Pa)
45!    pt(ngrid,nlayer)       Temperature (K)
46!    pdt(ngrid,nlayer)      Temperature tendency (K)
47!    pu(ngrid,nlayer)       u component of the wind (ms-1)
48!    pdu(ngrid,nlayer)      u component tendency (K)
49!    pv(ngrid,nlayer)       v component of the wind (ms-1)
50!    pdv(ngrid,nlayer)      v component tendency (K)
51!    dist_sol                   distance of the sun (AU)
52!    mu0(ngrid)               cos of solar zenith angle (=1 when sun at zenith)
53!    fract(ngrid)           day fraction (for diurnal=false)
54!    pq(ngrid,nlayer,nq)  Advected fields, ie chemical species here
55!    pdq(ngrid,nlayer,nq) Previous tendencies on pq
56!    tauref(ngrid)            Optical depth at 7 hPa
57!    co2ice(ngrid)            co2 ice surface layer (kg.m-2)
58!    surfdust(ngrid,nlayer) dust surface area (m2/m3)
59!    surfice(ngrid,nlayer)  ice surface area (m2/m3)
60!
61!  Output:
62!
63!    dqchim(ngrid,nlayer,nq) ! tendencies on pq due to chemistry
64!    dqschim(ngrid,nq)         ! tendencies on qsurf
65!
66!=======================================================================
67
68#include "chimiedata.h"
69
70!     input:
71
72      integer,intent(in) :: ngrid    ! number of atmospheric columns
73      integer,intent(in) :: nlayer   ! number of atmospheric layers
74      integer,intent(in) :: nq       ! number of tracers
75      real :: ptimestep
76      real :: pplay(ngrid,nlayer)    ! pressure at the middle of the layers
77      real :: zzlay(ngrid,nlayer)    ! pressure at the middle of the layers
78      real :: pplev(ngrid,nlayer+1)  ! intermediate pressure levels
79      real :: zzlev(ngrid,nlayer+1)  ! altitude at layer boundaries
80      real :: pt(ngrid,nlayer)       ! temperature
81      real :: pdt(ngrid,nlayer)      ! temperature tendency
82      real :: pu(ngrid,nlayer)       ! u component of the wind (m.s-1)
83      real :: pdu(ngrid,nlayer)      ! u component tendency
84      real :: pv(ngrid,nlayer)       ! v component of the wind (m.s-1)
85      real :: pdv(ngrid,nlayer)      ! v component tendency
86      real :: dist_sol               ! distance of the sun (AU)
87      real :: mu0(ngrid)             ! cos of solar zenith angle (=1 when sun at zenith)
88      real :: pq(ngrid,nlayer,nq)    ! tracers mass mixing ratio
89      real :: pdq(ngrid,nlayer,nq)   ! previous tendencies
90      real :: zday                   ! date (time since Ls=0, in martian days)
91      real :: tauref(ngrid)          ! optical depth at 7 hPa
92      real :: co2ice(ngrid)          ! co2 ice surface layer (kg.m-2)
93      real :: surfdust(ngrid,nlayer) ! dust surface area (m2/m3)
94      real :: surfice(ngrid,nlayer)  !  ice surface area (m2/m3)
95
96!     output:
97
98      real :: dqchim(ngrid,nlayer,nq)   ! tendencies on pq due to chemistry
99      real :: dqschim(ngrid,nq)         ! tendencies on qsurf
100
101!     local variables:
102
103      integer,save :: nbq                   ! number of tracers used in the chemistry
104      integer,allocatable,save :: niq(:)    ! array storing the indexes of the tracers
105      integer :: iloc(1)                    ! index of major species
106      integer :: ig,l,i,iq,iqmax
107      integer :: foundswitch, lswitch
108      integer,save :: chemthermod
109
110      integer,save :: i_co2  = 0
111      integer,save :: i_co   = 0
112      integer,save :: i_o    = 0
113      integer,save :: i_o1d  = 0
114      integer,save :: i_o2   = 0
115      integer,save :: i_o3   = 0
116      integer,save :: i_h    = 0
117      integer,save :: i_h2   = 0
118      integer,save :: i_oh   = 0
119      integer,save :: i_ho2  = 0
120      integer,save :: i_h2o2 = 0
121      integer,save :: i_ch4  = 0
122      integer,save :: i_n2   = 0
123      integer,save :: i_h2o  = 0
124      integer,save :: i_n    = 0
125      integer,save :: i_no    = 0
126      integer,save :: i_no2    = 0
127      integer,save :: i_n2d    = 0
128
129      integer :: ig_vl1
130
131      real    :: latvl1, lonvl1
132      real    :: zq(ngrid,nlayer,nq)   ! pq+pdq*ptimestep before chemistry
133                                       ! new mole fraction after
134      real    :: zt(ngrid,nlayer)      ! temperature
135      real    :: zu(ngrid,nlayer)      ! u component of the wind
136      real    :: zv(ngrid,nlayer)      ! v component of the wind
137      real    :: taucol                ! optical depth at 7 hPa
138
139      logical,save :: firstcall = .true.
140      logical,save :: depos = .false.  ! switch for dry deposition
141
142!     for each column of atmosphere:
143
144      real :: zpress(nlayer)       ! Pressure (mbar)
145      real :: zdens(nlayer)        ! Density  (cm-3)
146      real :: ztemp(nlayer)        ! Temperature (K)
147      real :: zlocal(nlayer)       ! Altitude (km)
148      real :: zycol(nlayer,nq)     ! Composition (mole fractions)
149      real :: zmmean(nlayer)       ! Mean molecular mass (g.mole-1)
150      real :: szacol               ! Solar zenith angle
151      real :: fract(ngrid)         ! day fraction
152      real :: fractcol             ! day fraction
153      real :: surfice1d(nlayer)    ! Ice surface area (cm2/cm3)
154      real :: surfdust1d(nlayer)   ! Dust surface area (cm2/cm3)
155      real :: jo3(nlayer)          ! Photodissociation rate O3->O1D (s-1)
156
157      integer :: iter(nlayer)      !  Number of chemical iterations
158                                   !  within one physical timestep
159
160!     for output:
161
162      logical :: output             ! to issue calls to writediagfi and stats
163      parameter (output = .true.)
164      real :: jo3_3d(ngrid,nlayer)  ! Photodissociation rate O3->O1D (s-1)
165      real :: iter_3d(ngrid,nlayer) ! Number of chemical iterations
166                                    !  within one physical timestep
167
168!=======================================================================
169!     initialization of the chemistry (first call only)
170!=======================================================================
171
172      if (firstcall) then
173
174         if (photochem) then
175            print*,'calchim: Read photolysis lookup table'
176            call read_phototable
177         end if
178         ! find index of chemical tracers to use
179         allocate(niq(nq))
180         ! Listed here are all tracers that can go into photochemistry
181         nbq = 0 ! to count number of tracers
182         ! Species ALWAYS present if photochem=.T. or thermochem=.T.
183         i_co2 = igcm_co2
184         if (i_co2 == 0) then
185            write(*,*) "calchim: Error; no CO2 tracer !!!"
186            stop
187         else
188            nbq = nbq + 1
189            niq(nbq) = i_co2
190         end if
191         i_co = igcm_co
192         if (i_co == 0) then
193            write(*,*) "calchim: Error; no CO tracer !!!"
194            stop
195         else
196            nbq = nbq + 1
197            niq(nbq) = i_co
198         end if
199         i_o = igcm_o
200         if (i_o == 0) then
201            write(*,*) "calchim: Error; no O tracer !!!"
202            stop
203         else
204            nbq = nbq + 1
205            niq(nbq) = i_o
206         end if
207         i_o1d = igcm_o1d
208         if (i_o1d == 0) then
209            write(*,*) "calchim: Error; no O1D tracer !!!"
210            stop
211         else
212            nbq = nbq + 1
213            niq(nbq) = i_o1d
214         end if
215         i_o2 = igcm_o2
216         if (i_o2 == 0) then
217            write(*,*) "calchim: Error; no O2 tracer !!!"
218            stop
219         else
220            nbq = nbq + 1
221            niq(nbq) = i_o2
222         end if
223         i_o3 = igcm_o3
224         if (i_o3 == 0) then
225            write(*,*) "calchim: Error; no O3 tracer !!!"
226            stop
227         else
228            nbq = nbq + 1
229            niq(nbq) = i_o3
230         end if
231         i_h = igcm_h
232         if (i_h == 0) then
233            write(*,*) "calchim: Error; no H tracer !!!"
234            stop
235         else
236            nbq = nbq + 1
237            niq(nbq) = i_h
238         end if
239         i_h2 = igcm_h2
240         if (i_h2 == 0) then
241            write(*,*) "calchim: Error; no H2 tracer !!!"
242            stop
243         else
244            nbq = nbq + 1
245            niq(nbq) = i_h2
246         end if
247         i_oh = igcm_oh
248         if (i_oh == 0) then
249            write(*,*) "calchim: Error; no OH tracer !!!"
250            stop
251         else
252            nbq = nbq + 1
253            niq(nbq) = i_oh
254         end if
255         i_ho2 = igcm_ho2
256         if (i_ho2 == 0) then
257            write(*,*) "calchim: Error; no HO2 tracer !!!"
258            stop
259         else
260            nbq = nbq + 1
261            niq(nbq) = i_ho2
262         end if
263         i_h2o2 = igcm_h2o2
264         if (i_h2o2 == 0) then
265            write(*,*) "calchim: Error; no H2O2 tracer !!!"
266            stop
267         else
268            nbq = nbq + 1
269            niq(nbq) = i_h2o2
270         end if
271         i_ch4 = igcm_ch4
272         if (i_ch4 == 0) then
273            write(*,*) "calchim: Error; no CH4 tracer !!!"
274            write(*,*) "CH4 will be ignored in the chemistry"
275         else
276            nbq = nbq + 1
277            niq(nbq) = i_ch4
278         end if
279         i_n2 = igcm_n2
280         if (i_n2 == 0) then
281            write(*,*) "calchim: Error; no N2 tracer !!!"
282            stop
283         else
284            nbq = nbq + 1
285            niq(nbq) = i_n2
286         end if
287         i_n = igcm_n
288         if (i_n == 0) then
289            if (photochem) then
290               write(*,*) "calchim: Error; no N tracer !!!"
291               write(*,*) "N will be ignored in the chemistry"
292            end if
293         else
294            nbq = nbq + 1
295            niq(nbq) = i_n
296         end if
297         i_n2d = igcm_n2d
298         if (i_n2d == 0) then
299            if (photochem) then
300               write(*,*) "calchim: Error; no N2D tracer !!!"
301               write(*,*) "N2d will be ignored in the chemistry"
302            end if
303         else
304            nbq = nbq + 1
305            niq(nbq) = i_n2d
306         end if
307         i_no = igcm_no
308         if (i_no == 0) then
309            if (photochem) then
310               write(*,*) "calchim: Error; no NO tracer !!!"
311               write(*,*) "NO will be ignored in the chemistry"
312            end if
313         else
314            nbq = nbq + 1
315            niq(nbq) = i_no
316         end if
317         i_no2 = igcm_no2
318         if (i_no2 == 0) then
319            if (photochem) then
320               write(*,*) "calchim: Error; no NO2 tracer !!!"
321               write(*,*) "NO2 will be ignored in the chemistry"
322            end if
323         else
324            nbq = nbq + 1
325            niq(nbq) = i_no2
326         end if
327         i_h2o = igcm_h2o_vap
328         if (i_h2o == 0) then
329            write(*,*) "calchim: Error; no water vapor tracer !!!"
330            stop
331         else
332            nbq = nbq + 1
333            niq(nbq) = i_h2o
334         end if
335         !Check tracers needed for thermospheric chemistry
336!         if(thermochem) then
337!            chemthermod=0  !Default: C/O/H chemistry
338!            !Nitrogen chemistry
339!            !NO is used to determine if N chemistry is wanted
340!            !chemthermod=2 -> N chemistry
341!            if (i_no == 0) then
342!               write(*,*) "calchim: no NO tracer"
343!               write(*,*) "C/O/H themosp chemistry only "
344!            else
345!               chemthermod=2
346!               write(*,*) "calchim: NO in traceur.def"
347!               write(*,*) "Nitrogen chemistry included"
348!            end if
349!            ! N
350!            if(chemthermod == 2) then
351!               if (i_n == 0) then
352!                  write(*,*) "calchim: Error; no N tracer !!!"
353!                  write(*,*) "N is needed if NO is in traceur.def"
354!                  stop
355!               end if
356!            ! NO2
357!               if (i_no2 == 0) then
358!                  write(*,*) "calchim: Error; no NO2 tracer !!!"
359!                  write(*,*) "NO2 is needed if NO is in traceur.def"
360!                  stop
361!               end if
362!            ! N(2D)
363!               if (i_n2d == 0) then
364!                  write(*,*) "calchim: Error; no N2D !!!"
365!                  write(*,*) "N2D is needed if NO is in traceur.def"
366!                  stop
367!               end if
368!            endif    !Of if(chemthermod == 2)
369!         endif      !Of thermochem
370
371         write(*,*) 'calchim: found nbq    = ',nbq,' tracers'
372               
373         firstcall = .false.
374      end if ! if (firstcall)
375
376! Initializations
377
378      zycol(:,:)    = 0.
379      dqchim(:,:,:) = 0.
380      dqschim(:,:)  = 0.
381
382!     latvl1= 22.27
383!     lonvl1= -47.94
384!     ig_vl1= 1+ int( (1.5-(latvl1-90.)*jjm/180.)  -2 )*iim +    &
385!             int(1.5+(lonvl1+180)*iim/360.)
386
387!=======================================================================
388!     loop over grid
389!=======================================================================
390
391      do ig = 1,ngrid
392         
393         foundswitch = 0
394         do l = 1,nlayer
395            do i = 1,nbq
396               iq = niq(i) ! get tracer index
397               zq(ig,l,iq) = pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep
398               zycol(l,iq) = zq(ig,l,iq)*mmean(ig,l)/mmol(iq)
399            end do
400            zt(ig,l)  = pt(ig,l) + pdt(ig,l)*ptimestep
401            zu(ig,l)  = pu(ig,l) + pdu(ig,l)*ptimestep
402            zv(ig,l)  = pv(ig,l) + pdv(ig,l)*ptimestep
403            zpress(l) = pplay(ig,l)/100.
404            ztemp(l)  = zt(ig,l)
405            zdens(l)  = zpress(l)/(kb*1.e4*ztemp(l))
406            zlocal(l) = zzlay(ig,l)/1000.
407            zmmean(l) = mmean(ig,l)
408
409!           surfdust1d and surfice1d: conversion from m2/m3 to cm2/cm3
410
411            surfdust1d(l) = surfdust(ig,l)*1.e-2
412            surfice1d(l)  = surfice(ig,l)*1.e-2
413
414!           search for switch index between regions
415
416!            if (photochem .and. thermochem) then
417!               if (foundswitch == 0 .and. pplay(ig,l) < 1.e-1) then
418!                  lswitch = l
419!                  foundswitch = 1
420!               end if
421!            end if
422            if (.not. photochem) then
423               lswitch = 22
424            end if
425!            if (.not. thermochem) then
426               lswitch = min(50,nlayer+1)
427!            end if
428
429         end do ! of do l=1,nlayer
430
431         szacol = acos(mu0(ig))*180./pi
432         taucol = tauref(ig)*(700./610.)  ! provisoire en attente de nouveau jmars
433         fractcol=fract(ig)
434
435!=======================================================================
436!     call chemical subroutines
437!=======================================================================
438
439!        chemistry in lower atmosphere
440
441         if (photochem) then
442
443            call photochemistry_asis(nlayer,nq,ngrid,                          &
444                                ig,lswitch,zycol,szacol,fractcol,ptimestep,    &
445                                zpress,ztemp,zdens,zmmean,dist_sol,            &
446                                surfdust1d,surfice1d,jo3,taucol,iter)
447
448!        ozone photolysis, for output
449
450            do l = 1,nlayer
451               jo3_3d(ig,l) = jo3(l)
452               iter_3d(ig,l) = iter(l)
453            end do
454
455!        condensation of h2o2
456
457!            call perosat(ngrid, nlayer, nq,                        &
458!                         ig,ptimestep,pplev,pplay,                 &
459!                         ztemp,zycol,dqcloud,dqscloud)
460         end if
461
462!        chemistry in upper atmosphere
463       
464!         if (thermochem) then
465!            call chemthermos(ig,nlayer,lswitch,chemthermod,zycol,ztemp,&
466!                             zdens,zpress,zlocal,szacol,ptimestep,zday)
467!         end if
468
469!        dry deposition
470
471!         if (depos) then
472!            call deposition(ngrid, nlayer, nq,                      &
473!                            ig, ig_vl1, pplay, pplev, zzlay, zzlev, &
474!                            zu, zv, zt, zycol, ptimestep, co2ice)
475!         end if
476
477!=======================================================================
478!     tendencies
479!=======================================================================
480
481!     index of the most abundant species at each level
482
483!         major(:) = maxloc(zycol, dim = 2)
484
485!     tendency for the most abundant species = - sum of others
486
487         do l = 1,nlayer
488            iloc=maxloc(zycol(l,:))
489            iqmax=iloc(1)
490            do i = 1,nbq
491               iq = niq(i) ! get tracer index
492               if (iq /= iqmax) then
493                  dqchim(ig,l,iq) = (zycol(l,iq)*mmol(iq)/mmean(ig,l) - zq(ig,l,iq))/ptimestep
494                  dqchim(ig,l,iq) = max(dqchim(ig,l,iq),-zq(ig,l,iq)/ptimestep)
495                  dqchim(ig,l,iqmax) = dqchim(ig,l,iqmax) - dqchim(ig,l,iq)
496               end if
497            end do
498         end do ! of do l = 1,nlayer
499
500
501
502
503
504
505!=======================================================================
506!     end of loop over grid
507!=======================================================================
508
509      end do ! of do ig=1,ngrid
510
511!=======================================================================
512!     write outputs
513!=======================================================================
514
515! value of parameter 'output' to trigger writting of outputs
516! is set above at the declaration of the variable.
517
518      if (photochem .and. output) then
519            call writediagfi(ngrid,'jo3','j o3->o1d',    &
520                             's-1',3,jo3_3d(1,1))
521            call writediagfi(ngrid,'iter','iterations',  &
522                             ' ',3,iter_3d(1,1))
523            if (callstats) then
524               call wstats(ngrid,'jo3','j o3->o1d',       &
525                           's-1',3,jo3_3d(1,1))
526               call wstats(ngrid,'mmean','mean molecular mass',       &
527                           'g.mole-1',3,mmean(1,1))
528            endif
529      end if ! of if (output)
530
531      return
532      end
Note: See TracBrowser for help on using the repository browser.