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

Last change on this file since 1543 was 1528, checked in by emillour, 9 years ago

Mars GCM:

  • Got rid of references to "dimensions.h" from physics packages: use nbp_lon (=iim), nbp_lat (==jjp1) and nbp_lev from module mod_grid_phy_lmdz (in phy_common) instead.
  • Added "ioipsl_getin_p_mod.F90" (getin_p routine) in phy_common to correctly read in parameters from *.def files in a parallel environment.

EM

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