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

Last change on this file since 2505 was 2461, checked in by emillour, 4 years ago

Mars GCM:

  • Adding the deuterium chemistry now that the HDO cycle is included.
  • Chemistry still works as before if deuterium tracers are not present.
  • Added handling of hdo in molecular diffusion (moldiff_red).

FGG+JYC+EM

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