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

Last change on this file since 459 was 459, checked in by emillour, 13 years ago

Mars GCM: Include changes and updates for photochemistry by FL:

  • aeronomars/calchim.F : change in units of surface density.
  • aeronomars/surfacearea.F : new routine to compute ice and dust surface area

(m2/m3) available for heterogeneous reactions.

  • phymars/initracer.F : bug correction: initialize igcm_ch4 and change loop

bounds (when guessing tracer names/properties with
old input files).

  • phymars/watercloud.F : cleanup.
  • phymars/physiq.F : add call to surfacearea; photochemistry is now called

after sedimentation to take into acount updated rdust
and rice.

EM

File size: 25.2 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-2
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_ccn_mass=0
216      igcm_ccn_number=0
217      igcm_dust_submicron=0
218      igcm_h2o_vap=0
219      igcm_h2o_ice=0
220      igcm_co2=0
221      igcm_co=0
222      igcm_o=0
223      igcm_o1d=0
224      igcm_o2=0
225      igcm_o3=0
226      igcm_h=0
227      igcm_h2=0
228      igcm_oh=0
229      igcm_ho2=0
230      igcm_h2o2=0
231      igcm_ch4=0
232      igcm_n2=0
233      igcm_ar=0
234      igcm_ar_n2=0
235      igcm_n=0
236      igcm_no=0
237      igcm_no2=0
238      igcm_n2d=0
239      igcm_co2plus=0
240      igcm_oplus=0
241      igcm_o2plus=0
242      igcm_coplus=0
243      igcm_cplus=0
244      igcm_nplus=0
245      igcm_noplus=0
246      igcm_n2plus=0
247      igcm_hplus=0
248      igcm_elec=0
249
250
251      ! 1. find dust tracers
252      count=0
253      if (dustbin.gt.0) then
254        do iq=1,nqmx
255          txt=" "
256          write(txt,'(a4,i2.2)')'dust',count+1
257          if (noms(iq).eq.txt) then
258            count=count+1
259            igcm_dustbin(count)=iq
260            mmol(iq)=100.
261          endif
262        enddo !do iq=1,nqmx
263      endif ! of if (dustbin.gt.0)
264      if (doubleq) then
265        do iq=1,nqmx
266          if (noms(iq).eq."dust_mass") then
267            igcm_dust_mass=iq
268            count=count+1
269          endif
270          if (noms(iq).eq."dust_number") then
271            igcm_dust_number=iq
272            count=count+1
273          endif
274        enddo
275      endif ! of if (doubleq)
276      if (scavenging) then
277        do iq=1,nqmx
278          if (noms(iq).eq."ccn_mass") then
279            igcm_ccn_mass=iq
280            count=count+1
281          endif
282          if (noms(iq).eq."ccn_number") then
283            igcm_ccn_number=iq
284            count=count+1
285          endif
286        enddo
287      endif ! of if (scavenging)
288      if (submicron) then
289        do iq=1,nqmx
290          if (noms(iq).eq."dust_submicron") then
291            igcm_dust_submicron=iq
292            mmol(iq)=100.
293            count=count+1
294          endif
295        enddo
296      endif ! of if (submicron)
297      ! 2. find chemistry and water tracers
298      do iq=1,nqmx
299        if (noms(iq).eq."co2") then
300          igcm_co2=iq
301          mmol(igcm_co2)=44.
302          count=count+1
303        endif
304        if (noms(iq).eq."co") then
305          igcm_co=iq
306          mmol(igcm_co)=28.
307          count=count+1
308        endif
309        if (noms(iq).eq."o") then
310          igcm_o=iq
311          mmol(igcm_o)=16.
312          count=count+1
313        endif
314        if (noms(iq).eq."o1d") then
315          igcm_o1d=iq
316          mmol(igcm_o1d)=16.
317          count=count+1
318        endif
319        if (noms(iq).eq."o2") then
320          igcm_o2=iq
321          mmol(igcm_o2)=32.
322          count=count+1
323        endif
324        if (noms(iq).eq."o3") then
325          igcm_o3=iq
326          mmol(igcm_o3)=48.
327          count=count+1
328        endif
329        if (noms(iq).eq."h") then
330          igcm_h=iq
331          mmol(igcm_h)=1.
332          count=count+1
333        endif
334        if (noms(iq).eq."h2") then
335          igcm_h2=iq
336          mmol(igcm_h2)=2.
337          count=count+1
338        endif
339        if (noms(iq).eq."oh") then
340          igcm_oh=iq
341          mmol(igcm_oh)=17.
342          count=count+1
343        endif
344        if (noms(iq).eq."ho2") then
345          igcm_ho2=iq
346          mmol(igcm_ho2)=33.
347          count=count+1
348        endif
349        if (noms(iq).eq."h2o2") then
350          igcm_h2o2=iq
351          mmol(igcm_h2o2)=34.
352          count=count+1
353        endif
354        if (noms(iq).eq."n2") then
355          igcm_n2=iq
356          mmol(igcm_n2)=28.
357          count=count+1
358        endif
359        if (noms(iq).eq."ch4") then
360          igcm_ch4=iq
361          mmol(igcm_ch4)=16.
362          count=count+1
363        endif
364        if (noms(iq).eq."ar") then
365          igcm_ar=iq
366          mmol(igcm_ar)=40.
367          count=count+1
368        endif
369        if (noms(iq).eq."n") then
370          igcm_n=iq
371          mmol(igcm_n)=14.
372          count=count+1
373        endif
374        if (noms(iq).eq."no") then
375          igcm_no=iq
376          mmol(igcm_no)=30.
377          count=count+1
378        endif
379        if (noms(iq).eq."no2") then
380          igcm_no2=iq
381          mmol(igcm_no2)=46.
382          count=count+1
383        endif
384        if (noms(iq).eq."n2d") then
385          igcm_n2d=iq
386          mmol(igcm_n2d)=28.
387          count=count+1
388        endif
389        if (noms(iq).eq."co2plus") then
390          igcm_co2plus=iq
391          mmol(igcm_co2plus)=44.
392          count=count+1
393        endif
394        if (noms(iq).eq."oplus") then
395          igcm_oplus=iq
396          mmol(igcm_oplus)=16.
397          count=count+1
398        endif
399        if (noms(iq).eq."o2plus") then
400          igcm_o2plus=iq
401          mmol(igcm_o2plus)=32.
402          count=count+1
403        endif
404        if (noms(iq).eq."coplus") then
405          igcm_coplus=iq
406          mmol(igcm_coplus)=28.
407          count=count+1
408        endif
409        if (noms(iq).eq."cplus") then
410          igcm_cplus=iq
411          mmol(igcm_cplus)=12.
412          count=count+1
413        endif
414        if (noms(iq).eq."nplus") then
415          igcm_nplus=iq
416          mmol(igcm_nplus)=14.
417          count=count+1
418        endif
419        if (noms(iq).eq."noplus") then
420          igcm_noplus=iq
421          mmol(igcm_noplus)=30.
422          count=count+1
423        endif
424        if (noms(iq).eq."n2plus") then
425          igcm_n2plus=iq
426          mmol(igcm_n2plus)=28.
427          count=count+1
428        endif
429        if (noms(iq).eq."hplus") then
430          igcm_hplus=iq
431          mmol(igcm_hplus)=1.
432          count=count+1
433        endif
434        if (noms(iq).eq."elec") then
435          igcm_elec=iq
436          mmol(igcm_elec)=1./1822.89
437          count=count+1
438        endif
439        if (noms(iq).eq."h2o_vap") then
440          igcm_h2o_vap=iq
441          mmol(igcm_h2o_vap)=18.
442          count=count+1
443        endif
444        if (noms(iq).eq."h2o_ice") then
445          igcm_h2o_ice=iq
446          mmol(igcm_h2o_ice)=18.
447          count=count+1
448        endif
449        ! Other stuff: e.g. for simulations using co2 + neutral gaz
450        if (noms(iq).eq."Ar_N2") then
451          igcm_ar_n2=iq
452          mmol(igcm_ar_n2)=30.
453          count=count+1
454        endif
455
456
457      enddo ! of do iq=1,nqmx
458!      count=count+nbqchem
459     
460      ! check that we identified all tracers:
461      if (count.ne.nqmx) then
462        write(*,*) "initracer: found only ",count," tracers"
463        write(*,*) "               expected ",nqmx
464        do iq=1,count
465          write(*,*)'      ',iq,' ',trim(noms(iq))
466        enddo
467        stop
468      else
469        write(*,*) "initracer: found all expected tracers, namely:"
470        do iq=1,nqmx
471          write(*,*)'      ',iq,' ',trim(noms(iq))
472        enddo
473      endif
474
475      ! if water cycle but iceparty=.false., there will nevertheless be
476      ! water ice at the surface (iceparty is not used anymore, but this
477      ! part is still relevant, as we want to stay compatible with the
478      ! older versions).
479      if (water.and.(igcm_h2o_ice.eq.0)) then
480        igcm_h2o_ice=igcm_h2o_vap ! so that qsurf(i_h2o_ice) is identified
481                                  ! even though there is no q(i_h2o_ice)
482      else
483       ! surface ice qsurf(i_h2o_ice) was loaded twice by phyetat0,
484       ! as qsurf(i_h2o_vap) & as qsurf(i_h2o_ice), so to be clean:
485       if (igcm_h2o_vap.ne.0) then
486         qsurf(1:ngridmx,igcm_h2o_vap)=0
487       endif
488      endif
489
490c------------------------------------------------------------
491c     Initialisation tracers ....
492c------------------------------------------------------------
493      call zerophys(nqmx,rho_q)
494
495      rho_dust=2500.  ! Mars dust density (kg.m-3)
496      rho_ice=920.    ! Water ice density (kg.m-3)
497      nuice_ref=0.1   ! Effective variance nueff of the
498                      ! water-ice size distribution
499      !!!nuice_sed=0.45   ! Sedimentation effective variance
500                      ! of the water-ice size distribution
501
502      if (doubleq) then
503c       "doubleq" technique
504c       -------------------
505c      (transport of mass and number mixing ratio)
506c       iq=1: Q mass mixing ratio, iq=2: N number mixing ratio
507
508        if( (nqmx.lt.2).or.(water.and.(nqmx.lt.4)) ) then
509          write(*,*)'initracer: nqmx is too low : nqmx=', nqmx
510          write(*,*)'water= ',water,' doubleq= ',doubleq   
511        end if
512
513        nueff_lift = 0.5
514        varian=sqrt(log(1.+nueff_lift))
515
516        rho_q(igcm_dust_mass)=rho_dust
517        rho_q(igcm_dust_number)=rho_dust
518
519c       Intermediate calcul for computing geometric mean radius r0
520c       as a function of mass and number mixing ratio Q and N
521c       (r0 = (r3n_q * Q/ N)^(1/3))
522        r3n_q = exp(-4.5*varian**2)*(3./4.)/(pi*rho_dust)
523
524c       Intermediate calcul for computing effective radius reff
525c       from geometric mean radius r0
526c       (reff = ref_r0 * r0)
527        ref_r0 = exp(2.5*varian**2)
528       
529c       lifted dust :
530c       '''''''''''
531        reff_lift = 3.0e-6 !3.e-6 !Effective radius of lifted dust (m)
532        alpha_devil(igcm_dust_mass)=9.e-9   !  dust devil lift mass coeff
533c       alpha_lift(igcm_dust_mass)=3.0e-15  !  Lifted mass coeff
534        alpha_lift(igcm_dust_mass)=1.e-6 !1.e-6 !Lifted mass coeff
535
536        r0_lift = reff_lift/ref_r0
537        alpha_devil(igcm_dust_number)=r3n_q*
538     &                        alpha_devil(igcm_dust_mass)/r0_lift**3
539        alpha_lift(igcm_dust_number)=r3n_q*
540     &                        alpha_lift(igcm_dust_mass)/r0_lift**3
541
542        radius(igcm_dust_mass) = reff_lift
543        radius(igcm_dust_number) = reff_lift
544
545        write(*,*) "initracer: doubleq_param reff_lift:", reff_lift
546        write(*,*) "initracer: doubleq_param nueff_lift:", nueff_lift
547        write(*,*) "initracer: doubleq_param alpha_lift:",
548     &    alpha_lift(igcm_dust_mass)
549      else
550
551       if (dustbin.gt.1) then
552        print*,'initracer: STOP!',
553     $   ' properties of dust need to be set in initracer !!!'
554        stop
555
556       else if (dustbin.eq.1) then
557
558c       This will be used for 1 dust particle size:
559c       ------------------------------------------
560        radius(igcm_dustbin(1))=3.e-6
561        alpha_lift(igcm_dustbin(1))=0.0e-6
562        alpha_devil(igcm_dustbin(1))=7.65e-9
563        rho_q(igcm_dustbin(1))=rho_dust
564
565       endif
566      end if    ! (doubleq)
567
568
569c     Scavenging of dust particles by H2O clouds:
570c     ------------------------------------------
571c     Initialize the two tracers used for the CCNs
572      if (water.AND.doubleq.AND.scavenging) then
573        radius(igcm_ccn_mass) = radius(igcm_dust_mass)
574        alpha_lift(igcm_ccn_mass) = 1e-30
575        alpha_devil(igcm_ccn_mass) = 1e-30
576        rho_q(igcm_ccn_mass) = rho_dust
577
578        radius(igcm_ccn_number) = radius(igcm_ccn_mass)
579        alpha_lift(igcm_ccn_number) = alpha_lift(igcm_ccn_mass)
580        alpha_devil(igcm_ccn_number) = alpha_devil(igcm_ccn_mass)
581        rho_q(igcm_ccn_number) = rho_q(igcm_ccn_mass)
582      endif ! of if (water.AND.doubleq.AND.scavenging)
583
584c     Submicron dust mode:
585c     --------------------
586
587      if (submicron) then
588        radius(igcm_dust_submicron)=0.1e-6
589        rho_q(igcm_dust_submicron)=rho_dust
590        if (doubleq) then
591c         If doubleq is also active, we use the population ratio:
592          alpha_lift(igcm_dust_submicron) =
593     &      alpha_lift(igcm_dust_number)*popratio*
594     &      rho_q(igcm_dust_submicron)*4./3.*pi*
595     &      radius(igcm_dust_submicron)**3.
596          alpha_devil(igcm_dust_submicron)=1.e-30
597        else
598          alpha_lift(igcm_dust_submicron)=1e-6
599          alpha_devil(igcm_dust_submicron)=1.e-30
600        endif ! (doubleq)
601      end if  ! (submicron)
602
603c     Initialization for photochemistry:
604c     ---------------------------------
605      if (photochem) then
606      ! initialize chemistry+water (water will be correctly initialized below)
607      ! by initializing everything which is not dust ...
608        do iq=1,nqmx
609          txt=noms(iq)
610          if (txt(1:4).ne."dust") then
611            radius(iq)=0.
612            alpha_lift(iq) =0.
613            alpha_devil(iq)=0.
614          endif
615        enddo ! do iq=1,nqmx
616      endif
617
618c     Initialization for water vapor
619c     ------------------------------
620      if(water) then
621         radius(igcm_h2o_vap)=0.
622         alpha_lift(igcm_h2o_vap) =0.
623         alpha_devil(igcm_h2o_vap)=0.
624         if(water.and.(nqmx.ge.2)) then
625           radius(igcm_h2o_ice)=3.e-6
626           rho_q(igcm_h2o_ice)=rho_ice
627           alpha_lift(igcm_h2o_ice) =0.
628           alpha_devil(igcm_h2o_ice)=0.
629         elseif(water.and.(nqmx.lt.2)) then
630            write(*,*) 'nqmx is too low : nqmx=', nqmx
631            write(*,*) 'water= ',water
632         endif
633
634      end if  ! (water)
635
636c     Output for records:
637c     ~~~~~~~~~~~~~~~~~~
638      write(*,*)
639      Write(*,*) '******** initracer : dust transport parameters :'
640      write(*,*) 'alpha_lift = ', alpha_lift
641      write(*,*) 'alpha_devil = ', alpha_devil
642      write(*,*) 'radius  = ', radius
643      if(doubleq) then
644        write(*,*) 'reff_lift (um) =  ', reff_lift
645        write(*,*) 'size distribution variance  = ', varian
646        write(*,*) 'r3n_q , ref_r0 : ', r3n_q , ref_r0
647      end if
648
649!
650!     some extra (possibly redundant) sanity checks for tracers:
651!     ---------------------------------------------------------
652
653       if (doubleq) then
654       ! verify that we indeed have dust_mass and dust_number tracers
655         if (igcm_dust_mass.eq.0) then
656           write(*,*) "initracer: error !!"
657           write(*,*) "  cannot use doubleq option without ",
658     &                "a dust_mass tracer !"
659           stop
660         endif
661         if (igcm_dust_number.eq.0) then
662           write(*,*) "initracer: error !!"
663           write(*,*) "  cannot use doubleq option without ",
664     &                "a dust_number tracer !"
665           stop
666         endif
667       endif
668
669       if ((.not.doubleq).and.(dustbin.gt.0)) then
670       ! verify that we indeed have 'dustbin' dust tracers
671         count=0
672         do iq=1,dustbin
673           if (igcm_dustbin(iq).ne.0) then
674             count=count+1
675           endif
676         enddo
677         if (count.ne.dustbin) then
678           write(*,*) "initracer: error !!"
679           write(*,*) "  dusbin is set to ",dustbin,
680     &                " but we only have the following dust tracers:"
681           do iq=1,count
682             write(*,*)"   ",trim(noms(igcm_dustbin(iq)))
683           enddo
684           stop
685         endif
686       endif
687
688       if (water) then
689       ! verify that we indeed have h2o_vap and h2o_ice tracers
690         if (igcm_h2o_vap.eq.0) then
691           write(*,*) "initracer: error !!"
692           write(*,*) "  cannot use water option without ",
693     &                "an h2o_vap tracer !"
694           stop
695         endif
696         if (igcm_h2o_ice.eq.0) then
697           write(*,*) "initracer: error !!"
698           write(*,*) "  cannot use water option without ",
699     &                "an h2o_ice tracer !"
700           stop
701         endif
702       endif
703
704       if (scavenging) then
705       ! verify that we indeed have ccn_mass and ccn_number tracers
706         if (igcm_ccn_mass.eq.0) then
707           write(*,*) "initracer: error !!"
708           write(*,*) "  cannot use scavenging option without ",
709     &                "a ccn_mass tracer !"
710           stop
711         endif
712         if (igcm_ccn_number.eq.0) then
713           write(*,*) "initracer: error !!"
714           write(*,*) "  cannot use scavenging option without ",
715     &                "a ccn_number tracer !"
716           stop
717         endif
718       endif ! of if (scavenging)
719
720       if (photochem .or. callthermos) then
721       ! verify that we indeed have the chemistry tracers
722         if (igcm_co2.eq.0) then
723           write(*,*) "initracer: error !!"
724           write(*,*) "  cannot use chemistry option without ",
725     &                "a co2 tracer !"
726         stop
727         endif
728         if (igcm_co.eq.0) then
729           write(*,*) "initracer: error !!"
730           write(*,*) "  cannot use chemistry option without ",
731     &                "a co tracer !"
732         stop
733         endif
734         if (igcm_o.eq.0) then
735           write(*,*) "initracer: error !!"
736           write(*,*) "  cannot use chemistry option without ",
737     &                "a o tracer !"
738         stop
739         endif
740         if (igcm_o1d.eq.0) then
741           write(*,*) "initracer: error !!"
742           write(*,*) "  cannot use chemistry option without ",
743     &                "a o1d tracer !"
744         stop
745         endif
746         if (igcm_o2.eq.0) then
747           write(*,*) "initracer: error !!"
748           write(*,*) "  cannot use chemistry option without ",
749     &                "an o2 tracer !"
750         stop
751         endif
752         if (igcm_o3.eq.0) then
753           write(*,*) "initracer: error !!"
754           write(*,*) "  cannot use chemistry option without ",
755     &                "an o3 tracer !"
756         stop
757         endif
758         if (igcm_h.eq.0) then
759           write(*,*) "initracer: error !!"
760           write(*,*) "  cannot use chemistry option without ",
761     &                "a h tracer !"
762         stop
763         endif
764         if (igcm_h2.eq.0) then
765           write(*,*) "initracer: error !!"
766           write(*,*) "  cannot use chemistry option without ",
767     &                "a h2 tracer !"
768         stop
769         endif
770         if (igcm_oh.eq.0) then
771           write(*,*) "initracer: error !!"
772           write(*,*) "  cannot use chemistry option without ",
773     &                "an oh tracer !"
774         stop
775         endif
776         if (igcm_ho2.eq.0) then
777           write(*,*) "initracer: error !!"
778           write(*,*) "  cannot use chemistry option without ",
779     &                "a ho2 tracer !"
780         stop
781         endif
782         if (igcm_h2o2.eq.0) then
783           write(*,*) "initracer: error !!"
784           write(*,*) "  cannot use chemistry option without ",
785     &                "a h2o2 tracer !"
786         stop
787         endif
788         if (igcm_n2.eq.0) then
789           write(*,*) "initracer: error !!"
790           write(*,*) "  cannot use chemistry option without ",
791     &                "a n2 tracer !"
792         stop
793         endif
794         if (igcm_ar.eq.0) then
795           write(*,*) "initracer: error !!"
796           write(*,*) "  cannot use chemistry option without ",
797     &                "an ar tracer !"
798         stop
799         endif
800       endif ! of if (photochem .or. callthermos)
801
802      end
Note: See TracBrowser for help on using the repository browser.