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

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

Mars GCM:

  • added modifications (from JYC & FGG) to tracer.h & initracer.F for ions
  • minor improvement to newstart.F (q=x option, check that tracer index provided by user is valid).

EM

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