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

Last change on this file since 2236 was 1918, checked in by emillour, 7 years ago

Mars GCM:
Code cleanup:

  • remove "comorbit.h" since it is no longer used.
  • turn "datafile.h" into module datafile_mod.F90 (and rename variable "datafile" as "datadir" since it stores the path to the datafile directory).

EM

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