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

Last change on this file since 283 was 283, checked in by aslmd, 13 years ago

LMDZ.MARS:

AS: J'ai fait ce commit et fait a la main un merge avec les modifications entre temps

24/08/11 == TN

Attempts to tune the water cycle by adding outliers

+ A few structural changes !!

  • watercap.h is now obsolete and removed -- all is in surfdat.h
  • watercaptag initialized in surfini.F (up to 7 areas defined) instead of initracer.F
    • settings proposed by AS commented
    • experiments by TN decommented. use with caution.
  • water ice albedo and thermal inertia in callphys.def and inifis.F
  • water ice albedo in surfini.F
  • water ice albedo computation in albedocaps.F90
  • alb_surfice is now obsolete in physiq.F, albedo_h2o_ice is used instead
  • frost_albedo_threshold defined in surfdat.h
  • water ice thermal inertia in soil.F

TODO: * calibrate thermal inertia and ice albedo

  • have a look at subgrid-scale ice with dryness ?
File size: 23.4 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."ar") then
345          igcm_ar=iq
346          mmol(igcm_ar)=40.
347          count=count+1
348        endif
349        if (noms(iq).eq."n") then
350          igcm_n=iq
351          mmol(igcm_n)=14.
352          count=count+1
353        endif
354        if (noms(iq).eq."no") then
355          igcm_no=iq
356          mmol(igcm_no)=30.
357          count=count+1
358        endif
359        if (noms(iq).eq."no2") then
360          igcm_no2=iq
361          mmol(igcm_no2)=46.
362          count=count+1
363        endif
364        if (noms(iq).eq."n2d") then
365          igcm_n2d=iq
366          mmol(igcm_n2d)=28.
367          count=count+1
368        endif
369        if (noms(iq).eq."co2plus") then
370          igcm_co2plus=iq
371          mmol(igcm_co2plus)=44.
372          count=count+1
373        endif
374        if (noms(iq).eq."oplus") then
375          igcm_oplus=iq
376          mmol(igcm_oplus)=16.
377          count=count+1
378        endif
379        if (noms(iq).eq."o2plus") then
380          igcm_o2plus=iq
381          mmol(igcm_o2plus)=32.
382          count=count+1
383        endif
384        if (noms(iq).eq."coplus") then
385          igcm_coplus=iq
386          mmol(igcm_coplus)=28.
387          count=count+1
388        endif
389        if (noms(iq).eq."cplus") then
390          igcm_cplus=iq
391          mmol(igcm_cplus)=12.
392          count=count+1
393        endif
394        if (noms(iq).eq."nplus") then
395          igcm_nplus=iq
396          mmol(igcm_nplus)=14.
397          count=count+1
398        endif
399        if (noms(iq).eq."noplus") then
400          igcm_noplus=iq
401          mmol(igcm_noplus)=30.
402          count=count+1
403        endif
404        if (noms(iq).eq."n2plus") then
405          igcm_n2plus=iq
406          mmol(igcm_n2plus)=28.
407          count=count+1
408        endif
409        if (noms(iq).eq."hplus") then
410          igcm_hplus=iq
411          mmol(igcm_hplus)=1.
412          count=count+1
413        endif
414        if (noms(iq).eq."elec") then
415          igcm_elec=iq
416          mmol(igcm_elec)=1./1822.89
417          count=count+1
418        endif
419        if (noms(iq).eq."h2o_vap") then
420          igcm_h2o_vap=iq
421          mmol(igcm_h2o_vap)=18.
422          count=count+1
423        endif
424        if (noms(iq).eq."h2o_ice") then
425          igcm_h2o_ice=iq
426          mmol(igcm_h2o_ice)=18.
427          count=count+1
428        endif
429        ! Other stuff: e.g. for simulations using co2 + neutral gaz
430        if (noms(iq).eq."Ar_N2") then
431          igcm_ar_n2=iq
432          mmol(igcm_ar_n2)=30.
433          count=count+1
434        endif
435       
436
437      enddo ! of do iq=1,nqmx
438!      count=count+nbqchem
439     
440      ! check that we identified all tracers:
441      if (count.ne.nqmx) then
442        write(*,*) "initracer: found only ",count," tracers"
443        write(*,*) "               expected ",nqmx
444        do iq=1,count
445          write(*,*)'      ',iq,' ',trim(noms(iq))
446        enddo
447        stop
448      else
449        write(*,*) "initracer: found all expected tracers, namely:"
450        do iq=1,nqmx
451          write(*,*)'      ',iq,' ',trim(noms(iq))
452        enddo
453      endif
454
455      ! if water cycle but iceparty=.false., there will nevertheless be
456      ! water ice at the surface (iceparty is not used anymore, but this
457      ! part is still relevant, as we want to stay compatible with the
458      ! older versions).
459      if (water.and.(igcm_h2o_ice.eq.0)) then
460        igcm_h2o_ice=igcm_h2o_vap ! so that qsurf(i_h2o_ice) is identified
461                                  ! even though there is no q(i_h2o_ice)
462      else
463       ! surface ice qsurf(i_h2o_ice) was loaded twice by phyetat0,
464       ! as qsurf(i_h2o_vap) & as qsurf(i_h2o_ice), so to be clean:
465       if (igcm_h2o_vap.ne.0) then
466         qsurf(1:ngridmx,igcm_h2o_vap)=0
467       endif
468      endif
469
470c------------------------------------------------------------
471c     Initialisation tracers ....
472c------------------------------------------------------------
473      call zerophys(nqmx,rho_q)
474
475      rho_dust=2500.  ! Mars dust density (kg.m-3)
476      rho_ice=920.    ! Water ice density (kg.m-3)
477      nuice_ref=0.1   ! Effective variance nueff of the
478                      ! water-ice size distributions
479
480      if (doubleq) then
481c       "doubleq" technique
482c       -------------------
483c      (transport of mass and number mixing ratio)
484c       iq=1: Q mass mixing ratio, iq=2: N number mixing ratio
485
486        if( (nqmx.lt.2).or.(water.and.(nqmx.lt.4)) ) then
487          write(*,*)'initracer: nqmx is too low : nqmx=', nqmx
488          write(*,*)'water= ',water,' doubleq= ',doubleq   
489        end if
490
491        nueff_lift = 0.5
492        varian=sqrt(log(1.+nueff_lift))
493
494        rho_q(igcm_dust_mass)=rho_dust
495        rho_q(igcm_dust_number)=rho_dust
496
497c       Intermediate calcul for computing geometric mean radius r0
498c       as a function of mass and number mixing ratio Q and N
499c       (r0 = (r3n_q * Q/ N)^(1/3))
500        r3n_q = exp(-4.5*varian**2)*(3./4.)/(pi*rho_dust)
501
502c       Intermediate calcul for computing effective radius reff
503c       from geometric mean radius r0
504c       (reff = ref_r0 * r0)
505        ref_r0 = exp(2.5*varian**2)
506       
507c       lifted dust :
508c       '''''''''''
509        reff_lift = 3.0e-6 !3.e-6 !Effective radius of lifted dust (m)
510        alpha_devil(igcm_dust_mass)=9.e-9   !  dust devil lift mass coeff
511c       alpha_lift(igcm_dust_mass)=3.0e-15  !  Lifted mass coeff
512        alpha_lift(igcm_dust_mass)=1.e-6 !1.e-6 !Lifted mass coeff
513
514        r0_lift = reff_lift/ref_r0
515        alpha_devil(igcm_dust_number)=r3n_q*
516     &                        alpha_devil(igcm_dust_mass)/r0_lift**3
517        alpha_lift(igcm_dust_number)=r3n_q*
518     &                        alpha_lift(igcm_dust_mass)/r0_lift**3
519
520        radius(igcm_dust_mass) = reff_lift
521        radius(igcm_dust_number) = reff_lift
522
523        write(*,*) "initracer: doubleq_param reff_lift:", reff_lift
524        write(*,*) "initracer: doubleq_param nueff_lift:", nueff_lift
525        write(*,*) "initracer: doubleq_param alpha_lift:",
526     &    alpha_lift(igcm_dust_mass)
527
528      else
529
530       if (dustbin.gt.1) then
531        print*,'initracer: STOP!',
532     $   ' properties of dust need to be set in initracer !!!'
533        stop
534
535       else if (dustbin.eq.1) then
536
537c       This will be used for 1 dust particle size:
538c       ------------------------------------------
539        radius(igcm_dustbin(1))=3.e-6
540        alpha_lift(igcm_dustbin(1))=0.0e-6
541        alpha_devil(igcm_dustbin(1))=7.65e-9
542        rho_q(igcm_dustbin(1))=rho_dust
543
544       endif
545      end if    ! (doubleq)
546
547c     Submicron dust mode:
548c     --------------------
549
550      if (submicron) then
551        radius(igcm_dust_submicron)=0.1e-6
552        rho_q(igcm_dust_submicron)=rho_dust
553        if (doubleq) then
554c         If doubleq is also active, we use the population ratio:
555          alpha_lift(igcm_dust_submicron) =
556     &      alpha_lift(igcm_dust_number)*popratio*
557     &      rho_q(igcm_dust_submicron)*4./3.*pi*
558     &      radius(igcm_dust_submicron)**3.
559          alpha_devil(igcm_dust_submicron)=1.e-30
560        else
561          alpha_lift(igcm_dust_submicron)=1e-6
562          alpha_devil(igcm_dust_submicron)=1.e-30
563        endif ! (doubleq)
564      end if  ! (submicron)
565
566c     Initialization for photochemistry:
567c     ---------------------------------
568      if (photochem) then
569      ! initialize chemistry+water (water will be correctly initialized below)
570      ! by initializing everything which is not dust ...
571        do iq=1,nqmx
572          txt=noms(iq)
573          if (txt(1:4).ne."dust") then
574            radius(iq)=0.
575            alpha_lift(iq) =0.
576            alpha_devil(iq)=0.
577          endif
578        enddo ! do iq=1,nqmx
579      endif
580
581c     Initialization for water vapor
582c     ------------------------------
583      if(water) then
584         radius(igcm_h2o_vap)=0.
585         alpha_lift(igcm_h2o_vap) =0.
586         alpha_devil(igcm_h2o_vap)=0.
587         if(water.and.(nqmx.ge.2)) then
588           radius(igcm_h2o_ice)=3.e-6
589           rho_q(igcm_h2o_ice)=rho_ice
590           alpha_lift(igcm_h2o_ice) =0.
591           alpha_devil(igcm_h2o_ice)=0.
592         elseif(water.and.(nqmx.lt.2)) then
593            write(*,*) 'nqmx is too low : nqmx=', nqmx
594            write(*,*) 'water= ',water
595         endif
596
597      end if  ! (water)
598
599c     Output for records:
600c     ~~~~~~~~~~~~~~~~~~
601      write(*,*)
602      Write(*,*) '******** initracer : dust transport parameters :'
603      write(*,*) 'alpha_lift = ', alpha_lift
604      write(*,*) 'alpha_devil = ', alpha_devil
605      write(*,*) 'radius  = ', radius
606      if(doubleq) then
607        write(*,*) 'reff_lift (um) =  ', reff_lift
608        write(*,*) 'size distribution variance  = ', varian
609        write(*,*) 'r3n_q , ref_r0 : ', r3n_q , ref_r0
610      end if
611
612!
613!     some extra (possibly redundant) sanity checks for tracers:
614!     ---------------------------------------------------------
615
616       if (doubleq) then
617       ! verify that we indeed have dust_mass and dust_number tracers
618         if (igcm_dust_mass.eq.0) then
619           write(*,*) "initracer: error !!"
620           write(*,*) "  cannot use doubleq option without ",
621     &                "a dust_mass tracer !"
622           stop
623         endif
624         if (igcm_dust_number.eq.0) then
625           write(*,*) "initracer: error !!"
626           write(*,*) "  cannot use doubleq option without ",
627     &                "a dust_number tracer !"
628           stop
629         endif
630       endif
631
632       if ((.not.doubleq).and.(dustbin.gt.0)) then
633       ! verify that we indeed have 'dustbin' dust tracers
634         count=0
635         do iq=1,dustbin
636           if (igcm_dustbin(iq).ne.0) then
637             count=count+1
638           endif
639         enddo
640         if (count.ne.dustbin) then
641           write(*,*) "initracer: error !!"
642           write(*,*) "  dusbin is set to ",dustbin,
643     &                " but we only have the following dust tracers:"
644           do iq=1,count
645             write(*,*)"   ",trim(noms(igcm_dustbin(iq)))
646           enddo
647           stop
648         endif
649       endif
650
651       if (water) then
652       ! verify that we indeed have h2o_vap and h2o_ice tracers
653         if (igcm_h2o_vap.eq.0) then
654           write(*,*) "initracer: error !!"
655           write(*,*) "  cannot use water option without ",
656     &                "an h2o_vap tracer !"
657           stop
658         endif
659         if (igcm_h2o_ice.eq.0) then
660           write(*,*) "initracer: error !!"
661           write(*,*) "  cannot use water option without ",
662     &                "an h2o_ice tracer !"
663           stop
664         endif
665       endif
666
667       if (photochem .or. callthermos) then
668       ! verify that we indeed have the chemistry tracers
669         if (igcm_co2.eq.0) then
670           write(*,*) "initracer: error !!"
671           write(*,*) "  cannot use chemistry option without ",
672     &                "a co2 tracer !"
673         stop
674         endif
675         if (igcm_co.eq.0) then
676           write(*,*) "initracer: error !!"
677           write(*,*) "  cannot use chemistry option without ",
678     &                "a co tracer !"
679         stop
680         endif
681         if (igcm_o.eq.0) then
682           write(*,*) "initracer: error !!"
683           write(*,*) "  cannot use chemistry option without ",
684     &                "a o tracer !"
685         stop
686         endif
687         if (igcm_o1d.eq.0) then
688           write(*,*) "initracer: error !!"
689           write(*,*) "  cannot use chemistry option without ",
690     &                "a o1d tracer !"
691         stop
692         endif
693         if (igcm_o2.eq.0) then
694           write(*,*) "initracer: error !!"
695           write(*,*) "  cannot use chemistry option without ",
696     &                "an o2 tracer !"
697         stop
698         endif
699         if (igcm_o3.eq.0) then
700           write(*,*) "initracer: error !!"
701           write(*,*) "  cannot use chemistry option without ",
702     &                "an o3 tracer !"
703         stop
704         endif
705         if (igcm_h.eq.0) then
706           write(*,*) "initracer: error !!"
707           write(*,*) "  cannot use chemistry option without ",
708     &                "a h tracer !"
709         stop
710         endif
711         if (igcm_h2.eq.0) then
712           write(*,*) "initracer: error !!"
713           write(*,*) "  cannot use chemistry option without ",
714     &                "a h2 tracer !"
715         stop
716         endif
717         if (igcm_oh.eq.0) then
718           write(*,*) "initracer: error !!"
719           write(*,*) "  cannot use chemistry option without ",
720     &                "an oh tracer !"
721         stop
722         endif
723         if (igcm_ho2.eq.0) then
724           write(*,*) "initracer: error !!"
725           write(*,*) "  cannot use chemistry option without ",
726     &                "a ho2 tracer !"
727         stop
728         endif
729         if (igcm_h2o2.eq.0) then
730           write(*,*) "initracer: error !!"
731           write(*,*) "  cannot use chemistry option without ",
732     &                "a h2o2 tracer !"
733         stop
734         endif
735         if (igcm_n2.eq.0) then
736           write(*,*) "initracer: error !!"
737           write(*,*) "  cannot use chemistry option without ",
738     &                "a n2 tracer !"
739         stop
740         endif
741         if (igcm_ar.eq.0) then
742           write(*,*) "initracer: error !!"
743           write(*,*) "  cannot use chemistry option without ",
744     &                "an ar tracer !"
745         stop
746         endif
747       endif ! of if (photochem .or. callthermos)
748
749      end
Note: See TracBrowser for help on using the repository browser.