source: trunk/LMDZ.MARS/libf/phymars/initracer.F @ 342

Last change on this file since 342 was 324, checked in by aslmd, 14 years ago

MESOSCALE : Preparatory commit for the ultimate option mars=42 which

would allow mesoscale modeling with photochemistry.

[see r76 method to add 'mars' options]
Modified module_lmd_driver.F and Registry.EM and runmeso
Modified readmeteo.F90 and introduced an option -DPHOTOCHEM
...
Transparent to the casual user
Option mars=42 not yet finished though -- so do not use!

File size: 23.5 KB
Line 
1      SUBROUTINE initracer(qsurf,co2ice)
2
3       IMPLICIT NONE
4c=======================================================================
5c   subject:
6c   --------
7c   Initialization related to tracer
8c   (transported dust, water, chemical species, ice...)
9c
10c   Name of the tracer
11c
12c   Test of dimension :
13c   Initialize COMMON tracer in tracer.h, using tracer names provided
14c   by the dynamics in "advtrac.h"
15c
16c   Old conventions: (not used any more)
17c
18c   If water=T : q(iq=nqmx) is the water mass mixing ratio
19c     and q(iq=nqmx-1) is the ice mass mixing ratio
20
21c   If there is transported dust, it uses iq=1 to iq=dustbin
22c   If there is no transported dust : dustbin=0
23c   If doubleq=T : q(iq=1) is the dust mass mixing ratio
24c                  q(iq=2) is the dust number mixing ratio
25
26c   If (photochem.or.thermochem) there is "ncomp" chemical species (ncomp
27c   is set in aeronomars/chimiedata.h) using the ncomp iq values starting at
28c      iq=nqchem_min = dustbin+1   (nqchem_min is defined in inifis.F)
29c   
30c
31c   author: F.Forget
32c   ------
33c    Modifs: Franck Montmessin, Sebastien Lebonnois (june 2003)
34c            Ehouarn Millour (oct. 2008) identify tracers by their names
35c=======================================================================
36
37
38#include "dimensions.h"
39#include "dimphys.h"
40#include "comcstfi.h"
41#include "callkeys.h"
42#include "tracer.h"
43#include "advtrac.h"
44#include "comgeomfi.h"
45#include "chimiedata.h"
46
47#include "surfdat.h"
48
49      real qsurf(ngridmx,nqmx)       ! tracer on surface (e.g.  kg.m-2)
50      real co2ice(ngridmx)           ! co2 ice mass on surface (e.g.  kg.m-2)
51      integer iq,ig,count
52      real r0_lift , reff_lift, nueff_lift
53c     Ratio of small over large dust particles (used when both
54c       doubleq and the submicron mode are active); In Montmessin
55c       et al. (2002), a value of 25 has been deduced;
56      real, parameter :: popratio = 25.
57      logical :: oldnames ! =.true. if old tracer naming convention (q01,...)
58      character(len=20) :: txt ! to store some text
59
60c-----------------------------------------------------------------------
61c  radius(nqmx)      ! aerosol particle radius (m)
62c  rho_q(nqmx)       ! tracer densities (kg.m-3)
63c  alpha_lift(nqmx)  ! saltation vertical flux/horiz flux ratio (m-1)
64c  alpha_devil(nqmx) ! lifting coeeficient by dust devil
65c  rho_dust          ! Mars dust density
66c  rho_ice           ! Water ice density
67c  nuice_ref         ! Effective variance nueff of the
68c                    !   water-ice size distributions
69c  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
70c  varian            ! Characteristic variance of log-normal distribution
71c-----------------------------------------------------------------------
72
73      integer :: nqchem_start
74
75! Initialization: get tracer names from the dynamics and check if we are
76!                 using 'old' tracer convention ('q01',q02',...)
77!                 or new convention (full tracer names)
78      ! check if tracers have 'old' names
79      oldnames=.false.
80      count=0
81      do iq=1,nqmx
82        txt=" "
83        write(txt,'(a1,i2.2)') 'q',iq
84        if (txt.eq.tnom(iq)) then
85          count=count+1
86        endif
87      enddo ! of do iq=1,nqmx
88     
89      if (count.eq.nqmx) then
90        write(*,*) "initracer: tracers seem to follow old naming ",
91     &             "convention (q01,q02,...)"
92        write(*,*) "   => will work for now ... "
93        write(*,*) "      but you should run newstart to rename them"
94        oldnames=.true.
95      endif
96
97      ! copy/set tracer names
98      if (oldnames) then ! old names (derived from iq & option values)
99        ! 1. dust:
100        if (dustbin.ne.0) then ! transported dust
101          do iq=1,dustbin
102            txt=" "
103            write(txt,'(a4,i2.2)') 'dust',iq
104            noms(iq)=txt
105            mmol(iq)=100.
106          enddo
107          ! special case if doubleq
108          if (doubleq) then
109            if (dustbin.ne.2) then
110              write(*,*) 'initracer: error expected dustbin=2'
111            else
112!              noms(1)='dust_mass'   ! dust mass mixing ratio
113!              noms(2)='dust_number' ! dust number mixing ratio
114! same as above, but avoid explict possible out-of-bounds on array
115               noms(1)='dust_mass'   ! dust mass mixing ratio
116               do iq=2,2
117               noms(iq)='dust_number' ! dust number mixing ratio
118               enddo
119            endif
120          endif
121        endif
122        ! 2. water & ice
123        if (water) then
124          noms(nqmx)='h2o_vap'
125          mmol(nqmx)=18.
126!            noms(nqmx-1)='h2o_ice'
127!            mmol(nqmx-1)=18.
128          !"loop" to avoid potential out-of-bounds in array
129          do iq=nqmx-1,nqmx-1
130            noms(iq)='h2o_ice'
131            mmol(iq)=18.
132          enddo
133        endif
134        ! 3. Chemistry
135        if (photochem .or. callthermos) then
136          if (doubleq) then
137            nqchem_start=3
138          else
139            nqchem_start=dustbin+1
140          end if
141          do iq=nqchem_start,nqchem_start+ncomp-1
142            noms(iq)=nomchem(iq-nqchem_start+1)
143            mmol(iq)=mmolchem(iq-nqchem_start+1)
144          enddo
145        endif ! of if (photochem .or. callthermos)
146        ! 4. Other tracers ????
147        if ((dustbin.eq.0).and.(.not.water)) then
148          noms(1)='co2'
149          mmol(1)=44
150          if (nqmx.eq.2) then
151            noms(nqmx)='Ar_N2'
152            mmol(nqmx)=30
153          endif
154        endif
155      else
156        ! copy tracer names from dynamics
157        do iq=1,nqmx
158          noms(iq)=tnom(iq)
159        enddo
160        ! mmol(:) array is initialized later (see below)
161      endif ! of if (oldnames)
162
163      ! special modification when using 'old' tracers:
164      ! qsurf(nqmx) was h2o ice whereas q(nqmx) was water vapour
165      ! and (if iceparty) q(nqmx-1) was null whereas q(nqmx-1) was water ice
166      if (oldnames.and.water) then
167        write(*,*)'initracer: moving surface water ice to index ',nqmx-1
168        ! "loop" to avoid potential out-of-bounds on arrays
169        do iq=nqmx-1,nqmx-1
170          qsurf(1:ngridmx,iq)=qsurf(1:ngridmx,iq+1)
171        enddo
172        qsurf(1:ngridmx,nqmx)=0
173      endif
174     
175c------------------------------------------------------------
176c     Test Dimensions tracers
177c------------------------------------------------------------
178
179! Ehouarn: testing number of tracers is obsolete...
180!      if(photochem.or.thermochem) then
181!          if (water) then
182!              if ((dustbin+ncomp+2).ne.nqmx) then
183!                 print*,'initracer: tracer dimension problem:'
184!                 print*,"(dustbin+ncomp+2).ne.nqmx"
185!                 print*,"ncomp: ",ncomp
186!                 print*,"dustbin: ",dustbin
187!                 print*,"nqmx: ",nqmx
188!                 print*,'Change ncomp in chimiedata.h'
189!               endif
190!          else
191!              if ((dustbin+ncomp+1).ne.nqmx) then
192!                 print*,'initracer: tracer dimension problem:'
193!                 print*,"(dustbin+ncomp+1).ne.nqmx"
194!                 print*,"ncomp: ",ncomp
195!                 print*,"dustbin: ",dustbin
196!                 print*,"nqmx: ",nqmx
197!                 print*,'Change ncomp in chimiedata.h'
198!                 STOP
199!               endif
200!            endif
201!      endif
202
203c------------------------------------------------------------
204c         NAME and molar mass of the tracer
205c------------------------------------------------------------
206
207   
208! Identify tracers by their names: (and set corresponding values of mmol)
209      ! 0. initialize tracer indexes to zero:
210      do iq=1,nqmx
211        igcm_dustbin(iq)=0
212      enddo
213      igcm_dust_mass=0
214      igcm_dust_number=0
215      igcm_dust_submicron=0
216      igcm_h2o_vap=0
217      igcm_h2o_ice=0
218      igcm_co2=0
219      igcm_co=0
220      igcm_o=0
221      igcm_o1d=0
222      igcm_o2=0
223      igcm_o3=0
224      igcm_h=0
225      igcm_h2=0
226      igcm_oh=0
227      igcm_ho2=0
228      igcm_h2o2=0
229      igcm_n2=0
230      igcm_ar=0
231      igcm_ar_n2=0
232      igcm_n=0
233      igcm_no=0
234      igcm_no2=0
235      igcm_n2d=0
236      igcm_co2plus=0
237      igcm_oplus=0
238      igcm_o2plus=0
239      igcm_coplus=0
240      igcm_cplus=0
241      igcm_nplus=0
242      igcm_noplus=0
243      igcm_n2plus=0
244      igcm_hplus=0
245      igcm_elec=0
246
247
248      ! 1. find dust tracers
249      count=0
250      if (dustbin.gt.0) then
251        do iq=1,nqmx
252          txt=" "
253          write(txt,'(a4,i2.2)')'dust',count+1
254          if (noms(iq).eq.txt) then
255            count=count+1
256            igcm_dustbin(count)=iq
257            mmol(iq)=100.
258          endif
259        enddo !do iq=1,nqmx
260      endif ! of if (dustbin.gt.0)
261      if (doubleq) then
262        do iq=1,nqmx
263          if (noms(iq).eq."dust_mass") then
264            igcm_dust_mass=iq
265            count=count+1
266          endif
267          if (noms(iq).eq."dust_number") then
268            igcm_dust_number=iq
269            count=count+1
270          endif
271        enddo
272      endif ! of if (doubleq)
273      if (submicron) then
274        do iq=1,nqmx
275          if (noms(iq).eq."dust_submicron") then
276            igcm_dust_submicron=iq
277            mmol(iq)=100.
278            count=count+1
279          endif
280        enddo
281      endif ! of if (submicron)
282      ! 2. find chemistry and water tracers
283      do iq=1,nqmx
284        if (noms(iq).eq."co2") then
285          igcm_co2=iq
286          mmol(igcm_co2)=44.
287          count=count+1
288        endif
289        if (noms(iq).eq."co") then
290          igcm_co=iq
291          mmol(igcm_co)=28.
292          count=count+1
293        endif
294        if (noms(iq).eq."o") then
295          igcm_o=iq
296          mmol(igcm_o)=16.
297          count=count+1
298        endif
299        if (noms(iq).eq."o1d") then
300          igcm_o1d=iq
301          mmol(igcm_o1d)=16.
302          count=count+1
303        endif
304        if (noms(iq).eq."o2") then
305          igcm_o2=iq
306          mmol(igcm_o2)=32.
307          count=count+1
308        endif
309        if (noms(iq).eq."o3") then
310          igcm_o3=iq
311          mmol(igcm_o3)=48.
312          count=count+1
313        endif
314        if (noms(iq).eq."h") then
315          igcm_h=iq
316          mmol(igcm_h)=1.
317          count=count+1
318        endif
319        if (noms(iq).eq."h2") then
320          igcm_h2=iq
321          mmol(igcm_h2)=2.
322          count=count+1
323        endif
324        if (noms(iq).eq."oh") then
325          igcm_oh=iq
326          mmol(igcm_oh)=17.
327          count=count+1
328        endif
329        if (noms(iq).eq."ho2") then
330          igcm_ho2=iq
331          mmol(igcm_ho2)=33.
332          count=count+1
333        endif
334        if (noms(iq).eq."h2o2") then
335          igcm_h2o2=iq
336          mmol(igcm_h2o2)=34.
337          count=count+1
338        endif
339        if (noms(iq).eq."n2") then
340          igcm_n2=iq
341          mmol(igcm_n2)=28.
342          count=count+1
343        endif
344        if (noms(iq).eq."ch4") then
345          igcm_ch4=iq
346          mmol(igcm_ch4)=16.
347          count=count+1
348        endif
349        if (noms(iq).eq."ar") then
350          igcm_ar=iq
351          mmol(igcm_ar)=40.
352          count=count+1
353        endif
354        if (noms(iq).eq."n") then
355          igcm_n=iq
356          mmol(igcm_n)=14.
357          count=count+1
358        endif
359        if (noms(iq).eq."no") then
360          igcm_no=iq
361          mmol(igcm_no)=30.
362          count=count+1
363        endif
364        if (noms(iq).eq."no2") then
365          igcm_no2=iq
366          mmol(igcm_no2)=46.
367          count=count+1
368        endif
369        if (noms(iq).eq."n2d") then
370          igcm_n2d=iq
371          mmol(igcm_n2d)=28.
372          count=count+1
373        endif
374        if (noms(iq).eq."co2plus") then
375          igcm_co2plus=iq
376          mmol(igcm_co2plus)=44.
377          count=count+1
378        endif
379        if (noms(iq).eq."oplus") then
380          igcm_oplus=iq
381          mmol(igcm_oplus)=16.
382          count=count+1
383        endif
384        if (noms(iq).eq."o2plus") then
385          igcm_o2plus=iq
386          mmol(igcm_o2plus)=32.
387          count=count+1
388        endif
389        if (noms(iq).eq."coplus") then
390          igcm_coplus=iq
391          mmol(igcm_coplus)=28.
392          count=count+1
393        endif
394        if (noms(iq).eq."cplus") then
395          igcm_cplus=iq
396          mmol(igcm_cplus)=12.
397          count=count+1
398        endif
399        if (noms(iq).eq."nplus") then
400          igcm_nplus=iq
401          mmol(igcm_nplus)=14.
402          count=count+1
403        endif
404        if (noms(iq).eq."noplus") then
405          igcm_noplus=iq
406          mmol(igcm_noplus)=30.
407          count=count+1
408        endif
409        if (noms(iq).eq."n2plus") then
410          igcm_n2plus=iq
411          mmol(igcm_n2plus)=28.
412          count=count+1
413        endif
414        if (noms(iq).eq."hplus") then
415          igcm_hplus=iq
416          mmol(igcm_hplus)=1.
417          count=count+1
418        endif
419        if (noms(iq).eq."elec") then
420          igcm_elec=iq
421          mmol(igcm_elec)=1./1822.89
422          count=count+1
423        endif
424        if (noms(iq).eq."h2o_vap") then
425          igcm_h2o_vap=iq
426          mmol(igcm_h2o_vap)=18.
427          count=count+1
428        endif
429        if (noms(iq).eq."h2o_ice") then
430          igcm_h2o_ice=iq
431          mmol(igcm_h2o_ice)=18.
432          count=count+1
433        endif
434        ! Other stuff: e.g. for simulations using co2 + neutral gaz
435        if (noms(iq).eq."Ar_N2") then
436          igcm_ar_n2=iq
437          mmol(igcm_ar_n2)=30.
438          count=count+1
439        endif
440
441
442      enddo ! of do iq=1,nqmx
443!      count=count+nbqchem
444     
445      ! check that we identified all tracers:
446      if (count.ne.nqmx) then
447        write(*,*) "initracer: found only ",count," tracers"
448        write(*,*) "               expected ",nqmx
449        do iq=1,count
450          write(*,*)'      ',iq,' ',trim(noms(iq))
451        enddo
452        stop
453      else
454        write(*,*) "initracer: found all expected tracers, namely:"
455        do iq=1,nqmx
456          write(*,*)'      ',iq,' ',trim(noms(iq))
457        enddo
458      endif
459
460      ! if water cycle but iceparty=.false., there will nevertheless be
461      ! water ice at the surface (iceparty is not used anymore, but this
462      ! part is still relevant, as we want to stay compatible with the
463      ! older versions).
464      if (water.and.(igcm_h2o_ice.eq.0)) then
465        igcm_h2o_ice=igcm_h2o_vap ! so that qsurf(i_h2o_ice) is identified
466                                  ! even though there is no q(i_h2o_ice)
467      else
468       ! surface ice qsurf(i_h2o_ice) was loaded twice by phyetat0,
469       ! as qsurf(i_h2o_vap) & as qsurf(i_h2o_ice), so to be clean:
470       if (igcm_h2o_vap.ne.0) then
471         qsurf(1:ngridmx,igcm_h2o_vap)=0
472       endif
473      endif
474
475c------------------------------------------------------------
476c     Initialisation tracers ....
477c------------------------------------------------------------
478      call zerophys(nqmx,rho_q)
479
480      rho_dust=2500.  ! Mars dust density (kg.m-3)
481      rho_ice=920.    ! Water ice density (kg.m-3)
482      nuice_ref=0.1   ! Effective variance nueff of the
483                      ! water-ice size distributions
484
485      if (doubleq) then
486c       "doubleq" technique
487c       -------------------
488c      (transport of mass and number mixing ratio)
489c       iq=1: Q mass mixing ratio, iq=2: N number mixing ratio
490
491        if( (nqmx.lt.2).or.(water.and.(nqmx.lt.4)) ) then
492          write(*,*)'initracer: nqmx is too low : nqmx=', nqmx
493          write(*,*)'water= ',water,' doubleq= ',doubleq   
494        end if
495
496        nueff_lift = 0.5
497        varian=sqrt(log(1.+nueff_lift))
498
499        rho_q(igcm_dust_mass)=rho_dust
500        rho_q(igcm_dust_number)=rho_dust
501
502c       Intermediate calcul for computing geometric mean radius r0
503c       as a function of mass and number mixing ratio Q and N
504c       (r0 = (r3n_q * Q/ N)^(1/3))
505        r3n_q = exp(-4.5*varian**2)*(3./4.)/(pi*rho_dust)
506
507c       Intermediate calcul for computing effective radius reff
508c       from geometric mean radius r0
509c       (reff = ref_r0 * r0)
510        ref_r0 = exp(2.5*varian**2)
511       
512c       lifted dust :
513c       '''''''''''
514        reff_lift = 3.0e-6 !3.e-6 !Effective radius of lifted dust (m)
515        alpha_devil(igcm_dust_mass)=9.e-9   !  dust devil lift mass coeff
516c       alpha_lift(igcm_dust_mass)=3.0e-15  !  Lifted mass coeff
517        alpha_lift(igcm_dust_mass)=1.e-6 !1.e-6 !Lifted mass coeff
518
519        r0_lift = reff_lift/ref_r0
520        alpha_devil(igcm_dust_number)=r3n_q*
521     &                        alpha_devil(igcm_dust_mass)/r0_lift**3
522        alpha_lift(igcm_dust_number)=r3n_q*
523     &                        alpha_lift(igcm_dust_mass)/r0_lift**3
524
525        radius(igcm_dust_mass) = reff_lift
526        radius(igcm_dust_number) = reff_lift
527
528        write(*,*) "initracer: doubleq_param reff_lift:", reff_lift
529        write(*,*) "initracer: doubleq_param nueff_lift:", nueff_lift
530        write(*,*) "initracer: doubleq_param alpha_lift:",
531     &    alpha_lift(igcm_dust_mass)
532
533      else
534
535       if (dustbin.gt.1) then
536        print*,'initracer: STOP!',
537     $   ' properties of dust need to be set in initracer !!!'
538        stop
539
540       else if (dustbin.eq.1) then
541
542c       This will be used for 1 dust particle size:
543c       ------------------------------------------
544        radius(igcm_dustbin(1))=3.e-6
545        alpha_lift(igcm_dustbin(1))=0.0e-6
546        alpha_devil(igcm_dustbin(1))=7.65e-9
547        rho_q(igcm_dustbin(1))=rho_dust
548
549       endif
550      end if    ! (doubleq)
551
552c     Submicron dust mode:
553c     --------------------
554
555      if (submicron) then
556        radius(igcm_dust_submicron)=0.1e-6
557        rho_q(igcm_dust_submicron)=rho_dust
558        if (doubleq) then
559c         If doubleq is also active, we use the population ratio:
560          alpha_lift(igcm_dust_submicron) =
561     &      alpha_lift(igcm_dust_number)*popratio*
562     &      rho_q(igcm_dust_submicron)*4./3.*pi*
563     &      radius(igcm_dust_submicron)**3.
564          alpha_devil(igcm_dust_submicron)=1.e-30
565        else
566          alpha_lift(igcm_dust_submicron)=1e-6
567          alpha_devil(igcm_dust_submicron)=1.e-30
568        endif ! (doubleq)
569      end if  ! (submicron)
570
571c     Initialization for photochemistry:
572c     ---------------------------------
573      if (photochem) then
574      ! initialize chemistry+water (water will be correctly initialized below)
575      ! by initializing everything which is not dust ...
576        do iq=1,nqmx
577          txt=noms(iq)
578          if (txt(1:4).ne."dust") then
579            radius(iq)=0.
580            alpha_lift(iq) =0.
581            alpha_devil(iq)=0.
582          endif
583        enddo ! do iq=1,nqmx
584      endif
585
586c     Initialization for water vapor
587c     ------------------------------
588      if(water) then
589         radius(igcm_h2o_vap)=0.
590         alpha_lift(igcm_h2o_vap) =0.
591         alpha_devil(igcm_h2o_vap)=0.
592         if(water.and.(nqmx.ge.2)) then
593           radius(igcm_h2o_ice)=3.e-6
594           rho_q(igcm_h2o_ice)=rho_ice
595           alpha_lift(igcm_h2o_ice) =0.
596           alpha_devil(igcm_h2o_ice)=0.
597         elseif(water.and.(nqmx.lt.2)) then
598            write(*,*) 'nqmx is too low : nqmx=', nqmx
599            write(*,*) 'water= ',water
600         endif
601
602      end if  ! (water)
603
604c     Output for records:
605c     ~~~~~~~~~~~~~~~~~~
606      write(*,*)
607      Write(*,*) '******** initracer : dust transport parameters :'
608      write(*,*) 'alpha_lift = ', alpha_lift
609      write(*,*) 'alpha_devil = ', alpha_devil
610      write(*,*) 'radius  = ', radius
611      if(doubleq) then
612        write(*,*) 'reff_lift (um) =  ', reff_lift
613        write(*,*) 'size distribution variance  = ', varian
614        write(*,*) 'r3n_q , ref_r0 : ', r3n_q , ref_r0
615      end if
616
617!
618!     some extra (possibly redundant) sanity checks for tracers:
619!     ---------------------------------------------------------
620
621       if (doubleq) then
622       ! verify that we indeed have dust_mass and dust_number tracers
623         if (igcm_dust_mass.eq.0) then
624           write(*,*) "initracer: error !!"
625           write(*,*) "  cannot use doubleq option without ",
626     &                "a dust_mass tracer !"
627           stop
628         endif
629         if (igcm_dust_number.eq.0) then
630           write(*,*) "initracer: error !!"
631           write(*,*) "  cannot use doubleq option without ",
632     &                "a dust_number tracer !"
633           stop
634         endif
635       endif
636
637       if ((.not.doubleq).and.(dustbin.gt.0)) then
638       ! verify that we indeed have 'dustbin' dust tracers
639         count=0
640         do iq=1,dustbin
641           if (igcm_dustbin(iq).ne.0) then
642             count=count+1
643           endif
644         enddo
645         if (count.ne.dustbin) then
646           write(*,*) "initracer: error !!"
647           write(*,*) "  dusbin is set to ",dustbin,
648     &                " but we only have the following dust tracers:"
649           do iq=1,count
650             write(*,*)"   ",trim(noms(igcm_dustbin(iq)))
651           enddo
652           stop
653         endif
654       endif
655
656       if (water) then
657       ! verify that we indeed have h2o_vap and h2o_ice tracers
658         if (igcm_h2o_vap.eq.0) then
659           write(*,*) "initracer: error !!"
660           write(*,*) "  cannot use water option without ",
661     &                "an h2o_vap tracer !"
662           stop
663         endif
664         if (igcm_h2o_ice.eq.0) then
665           write(*,*) "initracer: error !!"
666           write(*,*) "  cannot use water option without ",
667     &                "an h2o_ice tracer !"
668           stop
669         endif
670       endif
671
672       if (photochem .or. callthermos) then
673       ! verify that we indeed have the chemistry tracers
674         if (igcm_co2.eq.0) then
675           write(*,*) "initracer: error !!"
676           write(*,*) "  cannot use chemistry option without ",
677     &                "a co2 tracer !"
678         stop
679         endif
680         if (igcm_co.eq.0) then
681           write(*,*) "initracer: error !!"
682           write(*,*) "  cannot use chemistry option without ",
683     &                "a co tracer !"
684         stop
685         endif
686         if (igcm_o.eq.0) then
687           write(*,*) "initracer: error !!"
688           write(*,*) "  cannot use chemistry option without ",
689     &                "a o tracer !"
690         stop
691         endif
692         if (igcm_o1d.eq.0) then
693           write(*,*) "initracer: error !!"
694           write(*,*) "  cannot use chemistry option without ",
695     &                "a o1d tracer !"
696         stop
697         endif
698         if (igcm_o2.eq.0) then
699           write(*,*) "initracer: error !!"
700           write(*,*) "  cannot use chemistry option without ",
701     &                "an o2 tracer !"
702         stop
703         endif
704         if (igcm_o3.eq.0) then
705           write(*,*) "initracer: error !!"
706           write(*,*) "  cannot use chemistry option without ",
707     &                "an o3 tracer !"
708         stop
709         endif
710         if (igcm_h.eq.0) then
711           write(*,*) "initracer: error !!"
712           write(*,*) "  cannot use chemistry option without ",
713     &                "a h tracer !"
714         stop
715         endif
716         if (igcm_h2.eq.0) then
717           write(*,*) "initracer: error !!"
718           write(*,*) "  cannot use chemistry option without ",
719     &                "a h2 tracer !"
720         stop
721         endif
722         if (igcm_oh.eq.0) then
723           write(*,*) "initracer: error !!"
724           write(*,*) "  cannot use chemistry option without ",
725     &                "an oh tracer !"
726         stop
727         endif
728         if (igcm_ho2.eq.0) then
729           write(*,*) "initracer: error !!"
730           write(*,*) "  cannot use chemistry option without ",
731     &                "a ho2 tracer !"
732         stop
733         endif
734         if (igcm_h2o2.eq.0) then
735           write(*,*) "initracer: error !!"
736           write(*,*) "  cannot use chemistry option without ",
737     &                "a h2o2 tracer !"
738         stop
739         endif
740         if (igcm_n2.eq.0) then
741           write(*,*) "initracer: error !!"
742           write(*,*) "  cannot use chemistry option without ",
743     &                "a n2 tracer !"
744         stop
745         endif
746         if (igcm_ar.eq.0) then
747           write(*,*) "initracer: error !!"
748           write(*,*) "  cannot use chemistry option without ",
749     &                "an ar tracer !"
750         stop
751         endif
752       endif ! of if (photochem .or. callthermos)
753
754      end
Note: See TracBrowser for help on using the repository browser.