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

Last change on this file since 168 was 164, checked in by emillour, 14 years ago

Mars GCM:

Updates and corrections (to enable compiling/running in debug mode with

ifort)

  • removed option "-free-form" from makegcm_ifort and set mod_loc_dir="." so that module files (produced in local directory by ifort) are moved to LIBO
  • updated initracer.F, physdem1.F, physiq.F, inichim_newstart.F to avoid referencing out-of-bound array indexes (even if unused)
  • cosmetic updates on inwrite.F, datareadnc.F
  • updated newstart.F to initialize and use 'datadir' when looking for files
  • corrected bug on interpolation of sub-surface temperatures in lect_start_archive.F

EM

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