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

Last change on this file since 2537 was 1811, checked in by bclmd, 7 years ago

adding comment

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