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

Last change on this file since 455 was 455, checked in by tnavarro, 13 years ago

tiny change : nuice_sed initialisation is now done in inifis. Also changed initracer and improvedclouds.

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