source: trunk/mars/libf/phymars/initracer.F @ 131

Last change on this file since 131 was 91, checked in by aslmd, 14 years ago

mars: test outliers [dans initracer.F, commente] LMD_MM_MARS: modifications mineures [retrocompatibilite ancienne physique, callphys.def test pour nouvelle physique] PLOT: generalisation de la routine map_uvt pour pouvoir tracer des figures en projection polaire complete

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