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

Last change on this file since 2437 was 2409, checked in by emillour, 5 years ago

Mars GCM:
Code tidying : make a "dust_param_mod" module to store dust cycle relevant flags
and variables (and remove them from callkeys.h)
EM

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