source: trunk/LMDZ.MARS/libf/aeronomars/inichim_newstart.F90 @ 627

Last change on this file since 627 was 618, checked in by emillour, 13 years ago

Mars GCM:

Update of the chemistry package:

  • calchim.F90 :
    • upgraded from calchim.F
    • can run with or without CH4
    • fixed the mass conservation scheme
  • photochemistry.F : -chemistry timestep is now independant from physics timestep -can run with or without a CH4 tracer
    • removed initial tests on species, since these are already done in calchim
  • removed inichim_readcallphys.F (not used any more)
  • concentrations.F: adaptated to better handle indexes and molar masses of

tracers. "ncomp" (in chimiedata.h) no longer needed.

  • inichim_newstart.F90 : rewritten and cleaned (won't be compatible with

old start files where tracer names are q01,...).
Now handles hybrid levels and automaticaly adapts
depending on which tracers are available.

  • newstart.F : adapted to followup changes in inchim_newstart.F90 and some

cleanup around the initialization of tracer names and surface
values.

FL+EM

File size: 15.2 KB
Line 
1      subroutine inichim_newstart(pq, qsurf, ps, flagh2o, flagthermo)
2
3      implicit none
4
5!=======================================================================
6!
7!  Purpose:
8!  --------
9!
10!  Initialization of the chemistry for use in the building of a new start file
11!  used by program newstart.F
12!  also used by program testphys1d.F
13!
14!  Authors:
15!  -------
16!  Initial version 11/2002 by Sebastien Lebonnois
17!  Revised 07/2003 by Monica Angelats-i-Coll to use input files
18!  Modified 10/2008 Identify tracers by their names Ehouarn Millour
19!  Modified 11/2011 Addition of methane Franck Lefevre
20!  Rewritten 04/2012 Franck Lefevre
21!
22!  Arguments:
23!  ----------
24!
25!  pq(iip1,jjp1,llm,nqmx)  Advected fields, ie chemical species here
26!  qsurf(ngridmx,nqmx)     Amount of tracer on the surface (kg/m2)
27!  ps(iip1,jjp1)           Surface pressure (Pa)
28!  flagh2o                 flag for initialisation of h2o (1: yes / 0: no)
29!  flagthermo              flag for initialisation of thermosphere only (1: yes / 0: no)
30!
31!=======================================================================
32
33#include "dimensions.h"
34#include "dimphys.h"
35#include "paramet.h"
36#include "tracer.h"
37#include "comvert.h"
38#include "callkeys.h"
39#include "datafile.h"
40
41! inputs :
42
43      real,intent(in) :: ps(iip1,jjp1)            ! surface pressure in the gcm (Pa)   
44      integer,intent(in) :: flagh2o               ! flag for h2o initialisation
45      integer,intent(in) :: flagthermo            ! flag for thermosphere initialisation only
46
47! outputs :
48
49      real,intent(out) :: pq(iip1,jjp1,llm,nqmx)  ! advected fields, ie chemical species
50      real,intent(out) :: qsurf(ngridmx,nqmx)     ! surface values (kg/m2) of tracers
51
52! local :
53
54      integer :: iq, i, j, l, n, nbqchem
55      integer :: count, ierr, dummy
56      real    :: mmean(iip1,jjp1,llm)             ! mean molecular mass (g)
57      real    :: pgcm                             ! pressure at each layer in the gcm (Pa)
58
59      integer, parameter         :: nalt = 252    ! number of altitudes in the initialization files
60      integer, parameter         :: nspe = 14     ! number of species in the initialization files
61      integer, dimension(nspe)   :: niq           ! local index of species in initialization files
62      real, dimension(nalt)      :: tinit, zzfile ! temperature in initialization files
63      real, dimension(nalt)      :: pinit         ! pressure in initialization files
64      real, dimension(nalt)      :: densinit      ! total number density in initialization files
65      real, dimension(nalt,nspe) :: vmrinit       ! mixing ratios in initialization files
66      real, dimension(nspe)      :: vmrint        ! mixing ratio interpolated onto the gcm vertical grid
67      real                       :: vmr
68
69      character(len=20)          :: txt           ! to store some text
70
71! 1. identify tracers by their names: (and set corresponding values of mmol)
72
73! 1.1 initialize tracer indexes to zero:
74
75      do iq = 1,nqmx
76        igcm_dustbin(iq) = 0
77      end do
78
79      igcm_dust_mass   = 0
80      igcm_dust_number = 0
81      igcm_ccn_mass    = 0
82      igcm_ccn_number  = 0
83      igcm_h2o_vap     = 0
84      igcm_h2o_ice     = 0
85      igcm_co2         = 0
86      igcm_co          = 0
87      igcm_o           = 0
88      igcm_o1d         = 0
89      igcm_o2          = 0
90      igcm_o3          = 0
91      igcm_h           = 0
92      igcm_h2          = 0
93      igcm_oh          = 0
94      igcm_ho2         = 0
95      igcm_h2o2        = 0
96      igcm_ch4         = 0
97      igcm_n2          = 0
98      igcm_ar          = 0
99      igcm_ar_n2       = 0
100      igcm_co2plus     = 0
101      igcm_oplus       = 0
102      igcm_o2plus      = 0
103      igcm_coplus      = 0
104      igcm_cplus       = 0
105      igcm_nplus       = 0
106      igcm_noplus      = 0
107      igcm_n2plus      = 0
108      igcm_hplus       = 0
109      igcm_elec        = 0
110
111! 1.2 find dust tracers
112
113      count = 0
114
115      if (dustbin > 0) then
116         do iq = 1,nqmx
117            txt = " "
118            write(txt,'(a4,i2.2)') 'dust', count + 1
119            if (noms(iq) == txt) then
120               count = count + 1
121               igcm_dustbin(count) = iq
122               mmol(iq) = 100.
123            end if
124         end do !do iq=1,nqmx
125      end if ! of if (dustbin.gt.0)
126
127      if (doubleq) then
128         do iq = 1,nqmx
129            if (noms(iq) == "dust_mass") then
130               igcm_dust_mass = iq
131               count = count + 1
132            end if
133            if (noms(iq) == "dust_number") then
134               igcm_dust_number = iq
135               count = count + 1
136            end if
137         end do
138      end if ! of if (doubleq)
139
140      if (scavenging) then
141         do iq = 1,nqmx
142            if (noms(iq) == "ccn_mass") then
143               igcm_ccn_mass = iq
144               count = count + 1
145            end if
146            if (noms(iq) == "ccn_number") then
147               igcm_ccn_number = iq
148               count = count + 1
149            end if
150         end do
151      end if ! of if (scavenging)
152
153      if (submicron) then
154         do iq=1,nqmx
155            if (noms(iq) == "dust_submicron") then
156               igcm_dust_submicron = iq
157               mmol(iq) = 100.
158               count = count + 1
159            end if
160         end do
161      end if ! of if (submicron)
162
163! 1.3 find chemistry and water tracers
164
165      nbqchem = 0
166      do iq = 1,nqmx
167         if (noms(iq) == "co2") then
168            igcm_co2 = iq
169            mmol(igcm_co2) = 44.
170            count = count + 1
171            nbqchem = nbqchem + 1
172        end if
173        if (noms(iq) == "co") then
174           igcm_co = iq
175           mmol(igcm_co) = 28.
176           count = count + 1
177           nbqchem = nbqchem + 1
178        end if
179        if (noms(iq) == "o") then
180           igcm_o = iq
181           mmol(igcm_o) = 16.
182           count = count + 1
183           nbqchem = nbqchem + 1
184        end if
185        if (noms(iq) == "o1d") then
186           igcm_o1d = iq
187           mmol(igcm_o1d) = 16.
188           count = count + 1
189           nbqchem = nbqchem + 1
190        end if
191        if (noms(iq) == "o2") then
192           igcm_o2 = iq
193           mmol(igcm_o2) = 32.
194           count = count + 1
195           nbqchem = nbqchem + 1
196        end if
197        if (noms(iq) == "o3") then
198           igcm_o3 = iq
199           mmol(igcm_o3) = 48.
200           count = count + 1
201           nbqchem = nbqchem + 1
202        end if
203        if (noms(iq) == "h") then
204           igcm_h = iq
205           mmol(igcm_h) = 1.
206           count = count + 1
207           nbqchem = nbqchem + 1
208        end if
209        if (noms(iq) == "h2") then
210           igcm_h2 = iq
211           mmol(igcm_h2) = 2.
212           count = count + 1
213           nbqchem = nbqchem + 1
214        end if
215        if (noms(iq) == "oh") then
216           igcm_oh = iq
217           mmol(igcm_oh) = 17.
218           count = count + 1
219           nbqchem = nbqchem + 1
220        end if
221        if (noms(iq) == "ho2") then
222           igcm_ho2 = iq
223           mmol(igcm_ho2) = 33.
224           count = count + 1
225           nbqchem = nbqchem + 1
226        end if
227        if (noms(iq) == "h2o2") then
228           igcm_h2o2 = iq
229           mmol(igcm_h2o2) = 34.
230           count = count + 1
231           nbqchem = nbqchem + 1
232        end if
233        if (noms(iq) == "ch4") then
234           igcm_ch4 = iq
235           mmol(igcm_ch4) = 16.
236           count = count + 1
237           nbqchem = nbqchem + 1
238        end if
239        if (noms(iq) == "n2") then
240           igcm_n2 = iq
241           mmol(igcm_n2) = 28.
242           count = count + 1
243           nbqchem = nbqchem + 1
244        end if
245        if (noms(iq) == "n") then
246           igcm_n = iq
247           mmol(igcm_n) = 14.
248           count = count + 1
249           nbqchem = nbqchem + 1
250        end if
251        if (noms(iq) == "n2d") then
252           igcm_n2d = iq
253           mmol(igcm_n2d) = 14.
254           count = count + 1
255           nbqchem = nbqchem + 1
256        end if
257        if (noms(iq) == "no") then
258           igcm_no = iq
259           mmol(igcm_no) = 30.
260           count = count + 1
261           nbqchem = nbqchem + 1
262        end if
263        if (noms(iq) == "no2") then
264           igcm_no2 = iq
265           mmol(igcm_no2) = 46.
266           count = count + 1
267           nbqchem = nbqchem + 1
268        end if
269        if (noms(iq) == "ar") then
270           igcm_ar = iq
271           mmol(igcm_ar) = 40.
272           count = count + 1
273           nbqchem = nbqchem + 1
274        end if
275        if (noms(iq) == "h2o_vap") then
276           igcm_h2o_vap = iq
277           mmol(igcm_h2o_vap) = 18.
278           count = count + 1
279           nbqchem = nbqchem + 1
280        end if
281        if (noms(iq) == "h2o_ice") then
282           igcm_h2o_ice = iq
283           mmol(igcm_h2o_ice) = 18.
284           count = count + 1
285           nbqchem = nbqchem + 1
286        end if
287
288! 1.4 find ions
289
290        if (noms(iq) == "co2plus") then
291           igcm_co2plus = iq
292           mmol(igcm_co2plus) = 44.
293           count = count + 1
294           nbqchem = nbqchem + 1
295        end if
296        if (noms(iq) == "oplus") then
297           igcm_oplus = iq
298           mmol(igcm_oplus) = 16.
299           count = count + 1
300           nbqchem = nbqchem + 1
301        end if
302        if (noms(iq) == "o2plus") then
303           igcm_o2plus = iq
304           mmol(igcm_o2plus) = 32.
305           count = count + 1
306           nbqchem = nbqchem + 1
307        end if
308        if (noms(iq) == "coplus") then
309           igcm_coplus = iq
310           mmol(igcm_coplus) = 28.
311           count = count + 1
312           nbqchem = nbqchem + 1
313        end if
314        if (noms(iq) == "cplus") then
315           igcm_cplus = iq
316           mmol(igcm_cplus) = 12.
317           count = count + 1
318           nbqchem = nbqchem + 1
319        end if
320        if (noms(iq) == "nplus") then
321           igcm_nplus = iq
322           mmol(igcm_nplus) = 14.
323           count = count + 1
324           nbqchem = nbqchem + 1
325        end if
326        if (noms(iq) == "noplus") then
327           igcm_noplus = iq
328           mmol(igcm_noplus) = 30.
329           count = count + 1
330           nbqchem = nbqchem + 1
331        end if
332        if (noms(iq) == "n2plus") then
333           igcm_n2plus = iq
334           mmol(igcm_n2plus) = 28.
335           count = count + 1
336           nbqchem = nbqchem + 1
337        end if
338        if (noms(iq) == "hplus") then
339           igcm_hplus = iq
340           mmol(igcm_hplus) = 1.
341           count = count + 1
342           nbqchem = nbqchem + 1
343        end if
344        if (noms(iq) == "elec") then
345           igcm_elec = iq
346           mmol(igcm_elec) = 1./1822.89
347           count = count + 1
348        end if
349
350! 1.5 find idealized non-condensible tracer
351
352        if (noms(iq) == "Ar_N2") then
353           igcm_ar_n2 = iq
354           mmol(igcm_ar_n2) = 30.
355           count = count + 1
356        end if
357
358      end do ! of do iq=1,nqmx
359     
360! 1.6 check that we identified all tracers:
361
362      if (count /= nqmx) then
363         write(*,*) "inichim_newstart: found only ",count," tracers"
364         write(*,*) "                  expected ",nqmx
365         do iq = 1,count
366            write(*,*) '      ', iq, ' ', trim(noms(iq))
367         end do
368         stop
369      else
370         write(*,*) "inichim_newstart: found all expected tracers"
371         do iq = 1,nqmx
372            write(*,*) '      ', iq, ' ', trim(noms(iq))
373         end do
374      end if
375
376! 2. load in chemistry data for initialization:
377
378! order of major species in initialization file:
379!
380!    1: co2
381!    2: ar
382!    3: n2 
383!    4: o2 
384!    5: co 
385!    6: o   
386!    7: h2
387!
388! order of minor species in initialization file:
389!
390!    1: h 
391!    2: oh
392!    3: ho2
393!    4: h2o
394!    5: h2o2
395!    6: o1d
396!    7: o3
397
398! major species:
399
400      niq(1) = igcm_co2
401      niq(2) = igcm_ar
402      niq(3) = igcm_n2
403      niq(4) = igcm_o2
404      niq(5) = igcm_co
405      niq(6) = igcm_o
406      niq(7) = igcm_h2
407
408! minor species:
409
410      niq(8)  = igcm_h
411      niq(9)  = igcm_oh
412      niq(10) = igcm_ho2
413      niq(11) = igcm_h2o_vap
414      niq(12) = igcm_h2o2
415      niq(13) = igcm_o1d
416      niq(14) = igcm_o3
417
418! 2.1 open initialization files
419
420      open(210, iostat=ierr,file=trim(datafile)// &
421                            '/atmosfera_LMD_may.dat')
422      if (ierr /= 0) then
423         write(*,*)'Error : cannot open file atmosfera_LMD_may.dat '
424         write(*,*)'(in aeronomars/inichim.F)'
425         write(*,*)'It should be in :', trim(datafile),'/'
426         write(*,*)'1) You can change this path in callphys.def with'
427         write(*,*)'   datadir=/path/to/datafiles/'
428         write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)'
429         write(*,*)'   can be obtained online on:'
430         write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
431         stop
432      end if
433      open(220, iostat=ierr,file=trim(datafile)// &
434                            '/atmosfera_LMD_min.dat')
435      if (ierr /= 0) then
436         write(*,*)'Error : cannot open file atmosfera_LMD_min.dat '
437         write(*,*)'(in aeronomars/inichim.F)'
438         write(*,*)'It should be in :', trim(datafile),'/'
439         write(*,*)'1) You can change this path in callphys.def with'
440         write(*,*)'   datadir=/path/to/datafiles/'
441         write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)'
442         write(*,*)'   can be obtained online on:'
443         write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
444         stop
445      end if
446
447! 2.2 read initialization files
448
449! major species
450
451      read(210,*)
452      do l = 1,nalt
453         read(210,*) dummy, tinit(l), pinit(l), densinit(l), &
454                     (vmrinit(l,n), n = 1,7)
455         pinit(l) = pinit(l)*100.              ! conversion in Pa
456         pinit(l) = log(pinit(l))              ! for the vertical interpolation
457      end do
458      close(210)
459
460! minor species
461
462      read(220,*)
463      do l = 1,nalt
464         read(220,*) dummy, (vmrinit(l,n), n = 8,14)
465      end do
466      close(220)
467
468! 3. initialization of tracers
469
470      do i = 1,iip1
471         do j = 1,jjp1
472            do l = 1,llm
473
474               pgcm = aps(l) + bps(l)*ps(i,j)  ! gcm pressure
475               pgcm = log(pgcm)                ! for the vertical interpolation
476               mmean(i,j,l) = 0.
477
478! 3.1 vertical interpolation
479
480               do n = 1,nspe
481                  call intrplf(pgcm,vmr,pinit,vmrinit(:,n),nalt)
482                  vmrint(n) = vmr
483                  iq = niq(n)
484                  mmean(i,j,l) = mmean(i,j,l) + vmrint(n)*mmol(iq)
485               end do
486
487! 3.2 attribute mixing ratio: - all layers or only thermosphere
488!                             - with our without h2o
489
490               if (flagthermo == 0 .or. (flagthermo == 1 .and. exp(pgcm) < 1.e-3)) then
491                  do n = 1,nspe
492                     iq = niq(n)
493                     if (iq /= igcm_h2o_vap .or. flagh2o == 1) then
494                        pq(i,j,l,iq) = vmrint(n)*mmol(iq)/mmean(i,j,l)
495                     end if
496                  end do
497               end if
498
499            end do
500         end do
501      end do
502
503      ! set surface values of chemistry tracers to zero
504      if (flagthermo == 0) then
505        ! NB: no problem for "surface water vapour" tracer which is always 0
506        do n=1,nspe
507          iq=niq(n)
508          qsurf(1:ngridmx,iq)=0
509        enddo
510      endif
511
512
513! 3.3 initialization of tracers not contained in the initialization files
514
515! methane : 10 ppbv
516
517      if (igcm_ch4 /= 0) then
518         vmr = 10.e-9       
519         do i = 1,iip1
520            do j = 1,jjp1
521               do l = 1,llm
522                  pq(i,j,l,igcm_ch4) = vmr*mmol(igcm_ch4)/mmean(i,j,l)
523               end do
524            end do
525         end do
526         ! set surface value to zero
527         qsurf(1:ngridmx,igcm_ch4)=0
528      end if
529
530      end
Note: See TracBrowser for help on using the repository browser.