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

Last change on this file since 2566 was 2562, checked in by cmathe, 4 years ago

GCM MARS: CO2 clouds microphysics improvements

File size: 33.0 KB
RevLine 
[1224]1      SUBROUTINE initracer(ngrid,nq,qsurf)
[38]2
[1036]3       use tracer_mod
[2321]4       use comcstfi_h, only: pi
[2409]5       use dust_param_mod, only: doubleq, submicron, dustbin
[38]6       IMPLICIT NONE
7c=======================================================================
8c   subject:
9c   --------
10c   Initialization related to tracer
11c   (transported dust, water, chemical species, ice...)
12c
13c   Name of the tracer
14c
15c   Test of dimension :
[1036]16c   Initialize tracer related data in tracer_mod, using tracer names provided
17c   by the dynamics in "infotrac"
[38]18c
19c
20c   author: F.Forget
21c   ------
22c    Modifs: Franck Montmessin, Sebastien Lebonnois (june 2003)
23c            Ehouarn Millour (oct. 2008) identify tracers by their names
24c=======================================================================
25
26
[1974]27      include "callkeys.h"
[38]28
[1036]29      integer,intent(in) :: ngrid ! number of atmospheric columns
30      integer,intent(in) :: nq ! number of tracers
31      real,intent(out) :: qsurf(ngrid,nq) ! tracer on surface (e.g.  kg.m-2)
32
[283]33      integer iq,ig,count
[38]34      real r0_lift , reff_lift, nueff_lift
[1974]35      real r0_storm,reff_storm
[38]36c     Ratio of small over large dust particles (used when both
37c       doubleq and the submicron mode are active); In Montmessin
38c       et al. (2002), a value of 25 has been deduced;
39      real, parameter :: popratio = 25.
[1974]40      character(len=30) :: txt ! to store some text
[38]41
42c-----------------------------------------------------------------------
[1036]43c  radius(nq)      ! aerosol particle radius (m)
44c  rho_q(nq)       ! tracer densities (kg.m-3)
45c  alpha_lift(nq)  ! saltation vertical flux/horiz flux ratio (m-1)
46c  alpha_devil(nq) ! lifting coeeficient by dust devil
[38]47c  rho_dust          ! Mars dust density
48c  rho_ice           ! Water ice density
49c  nuice_ref         ! Effective variance nueff of the
50c                    !   water-ice size distributions
51c  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
52c  varian            ! Characteristic variance of log-normal distribution
53c-----------------------------------------------------------------------
54
[1036]55
[38]56c------------------------------------------------------------
57c         NAME and molar mass of the tracer
58c------------------------------------------------------------
59   
60! Identify tracers by their names: (and set corresponding values of mmol)
61      ! 0. initialize tracer indexes to zero:
[1036]62      igcm_dustbin(1:nq)=0
[1617]63      igcm_co2_ice=0
64      igcm_ccnco2_mass=0
65      igcm_ccnco2_number=0
[2562]66      igcm_ccnco2_h2o_mass_ice=0
67      igcm_ccnco2_h2o_mass_ccn=0
68      igcm_ccnco2_h2o_number=0
[38]69      igcm_dust_mass=0
70      igcm_dust_number=0
[358]71      igcm_ccn_mass=0
72      igcm_ccn_number=0
[38]73      igcm_dust_submicron=0
74      igcm_h2o_vap=0
75      igcm_h2o_ice=0
[2312]76      igcm_hdo_vap=0
77      igcm_hdo_ice=0
[1974]78      igcm_stormdust_mass=0
79      igcm_stormdust_number=0
[2199]80      igcm_topdust_mass=0
81      igcm_topdust_number=0
[38]82      igcm_co2=0
83      igcm_co=0
84      igcm_o=0
85      igcm_o1d=0
86      igcm_o2=0
87      igcm_o3=0
88      igcm_h=0
[2461]89      igcm_d=0
90      igcm_hd=0
[38]91      igcm_h2=0
[2461]92      igcm_od=0
93      igcm_do2=0
94      igcm_hdo2=0
[38]95      igcm_oh=0
96      igcm_ho2=0
97      igcm_h2o2=0
[459]98      igcm_ch4=0
[38]99      igcm_n2=0
100      igcm_ar=0
101      igcm_ar_n2=0
[171]102      igcm_n=0
103      igcm_no=0
104      igcm_no2=0
105      igcm_n2d=0
[1660]106      igcm_he=0
[171]107      igcm_co2plus=0
108      igcm_oplus=0
109      igcm_o2plus=0
110      igcm_coplus=0
111      igcm_cplus=0
112      igcm_nplus=0
113      igcm_noplus=0
114      igcm_n2plus=0
115      igcm_hplus=0
[635]116      igcm_hco2plus=0
[2284]117      igcm_hcoplus=0
[2302]118      igcm_h2oplus=0
[2321]119      igcm_h3oplus=0
120      igcm_ohplus=0
[171]121      igcm_elec=0
[38]122
123      ! 1. find dust tracers
124      count=0
125      if (dustbin.gt.0) then
[1036]126        do iq=1,nq
[38]127          txt=" "
128          write(txt,'(a4,i2.2)')'dust',count+1
129          if (noms(iq).eq.txt) then
130            count=count+1
131            igcm_dustbin(count)=iq
132            mmol(iq)=100.
133          endif
[1036]134        enddo !do iq=1,nq
[38]135      endif ! of if (dustbin.gt.0)
136      if (doubleq) then
[1036]137        do iq=1,nq
[38]138          if (noms(iq).eq."dust_mass") then
139            igcm_dust_mass=iq
140            count=count+1
141          endif
142          if (noms(iq).eq."dust_number") then
143            igcm_dust_number=iq
144            count=count+1
145          endif
146        enddo
147      endif ! of if (doubleq)
[740]148      if (microphys) then
[1036]149        do iq=1,nq
[358]150          if (noms(iq).eq."ccn_mass") then
151            igcm_ccn_mass=iq
152            count=count+1
153          endif
154          if (noms(iq).eq."ccn_number") then
155            igcm_ccn_number=iq
156            count=count+1
157          endif
158        enddo
[740]159      endif ! of if (microphys)
[38]160      if (submicron) then
[1036]161        do iq=1,nq
[38]162          if (noms(iq).eq."dust_submicron") then
163            igcm_dust_submicron=iq
164            mmol(iq)=100.
165            count=count+1
166          endif
167        enddo
168      endif ! of if (submicron)
[1974]169       if (rdstorm) then
170        do iq=1,nq
171          if (noms(iq).eq."stormdust_mass") then
172            igcm_stormdust_mass=iq
173            count=count+1
174          endif
175          if (noms(iq).eq."stormdust_number") then
176            igcm_stormdust_number=iq
177            count=count+1
178          endif
179        enddo
[2199]180      endif ! of if (rdstorm)
181       if (slpwind) then
182        do iq=1,nq
183          if (noms(iq).eq."topdust_mass") then
184            igcm_topdust_mass=iq
185            count=count+1
186          endif
187          if (noms(iq).eq."topdust_number") then
188            igcm_topdust_number=iq
189            count=count+1
190          endif
191        enddo
192      endif ! of if (slpwind)   
[38]193      ! 2. find chemistry and water tracers
[1036]194      do iq=1,nq
[38]195        if (noms(iq).eq."co2") then
196          igcm_co2=iq
197          mmol(igcm_co2)=44.
198          count=count+1
199        endif
200        if (noms(iq).eq."co") then
201          igcm_co=iq
202          mmol(igcm_co)=28.
203          count=count+1
204        endif
205        if (noms(iq).eq."o") then
206          igcm_o=iq
207          mmol(igcm_o)=16.
208          count=count+1
209        endif
210        if (noms(iq).eq."o1d") then
211          igcm_o1d=iq
212          mmol(igcm_o1d)=16.
213          count=count+1
214        endif
215        if (noms(iq).eq."o2") then
216          igcm_o2=iq
217          mmol(igcm_o2)=32.
218          count=count+1
219        endif
220        if (noms(iq).eq."o3") then
221          igcm_o3=iq
222          mmol(igcm_o3)=48.
223          count=count+1
224        endif
225        if (noms(iq).eq."h") then
226          igcm_h=iq
227          mmol(igcm_h)=1.
228          count=count+1
229        endif
230        if (noms(iq).eq."h2") then
231          igcm_h2=iq
232          mmol(igcm_h2)=2.
233          count=count+1
234        endif
235        if (noms(iq).eq."oh") then
236          igcm_oh=iq
237          mmol(igcm_oh)=17.
238          count=count+1
239        endif
240        if (noms(iq).eq."ho2") then
241          igcm_ho2=iq
242          mmol(igcm_ho2)=33.
243          count=count+1
244        endif
245        if (noms(iq).eq."h2o2") then
246          igcm_h2o2=iq
247          mmol(igcm_h2o2)=34.
248          count=count+1
249        endif
250        if (noms(iq).eq."n2") then
251          igcm_n2=iq
252          mmol(igcm_n2)=28.
253          count=count+1
254        endif
[324]255        if (noms(iq).eq."ch4") then
256          igcm_ch4=iq
257          mmol(igcm_ch4)=16.
258          count=count+1
259        endif
[38]260        if (noms(iq).eq."ar") then
261          igcm_ar=iq
262          mmol(igcm_ar)=40.
263          count=count+1
264        endif
[324]265        if (noms(iq).eq."n") then
[171]266          igcm_n=iq
267          mmol(igcm_n)=14.
268          count=count+1
269        endif
270        if (noms(iq).eq."no") then
271          igcm_no=iq
272          mmol(igcm_no)=30.
273          count=count+1
274        endif
275        if (noms(iq).eq."no2") then
276          igcm_no2=iq
277          mmol(igcm_no2)=46.
278          count=count+1
279        endif
280        if (noms(iq).eq."n2d") then
281          igcm_n2d=iq
282          mmol(igcm_n2d)=28.
283          count=count+1
284        endif
[1660]285        if (noms(iq).eq."he") then
286          igcm_he=iq
287          mmol(igcm_he)=4.
288          count=count+1
289        endif
[171]290        if (noms(iq).eq."co2plus") then
291          igcm_co2plus=iq
292          mmol(igcm_co2plus)=44.
293          count=count+1
294        endif
295        if (noms(iq).eq."oplus") then
296          igcm_oplus=iq
297          mmol(igcm_oplus)=16.
298          count=count+1
299        endif
300        if (noms(iq).eq."o2plus") then
301          igcm_o2plus=iq
302          mmol(igcm_o2plus)=32.
303          count=count+1
304        endif
305        if (noms(iq).eq."coplus") then
306          igcm_coplus=iq
307          mmol(igcm_coplus)=28.
308          count=count+1
309        endif
310        if (noms(iq).eq."cplus") then
311          igcm_cplus=iq
312          mmol(igcm_cplus)=12.
313          count=count+1
314        endif
315        if (noms(iq).eq."nplus") then
316          igcm_nplus=iq
317          mmol(igcm_nplus)=14.
318          count=count+1
319        endif
320        if (noms(iq).eq."noplus") then
321          igcm_noplus=iq
322          mmol(igcm_noplus)=30.
323          count=count+1
324        endif
[324]325        if (noms(iq).eq."n2plus") then
[171]326          igcm_n2plus=iq
327          mmol(igcm_n2plus)=28.
328          count=count+1
329        endif
330        if (noms(iq).eq."hplus") then
331          igcm_hplus=iq
332          mmol(igcm_hplus)=1.
333          count=count+1
334        endif
[635]335        if (noms(iq).eq."hco2plus") then
336          igcm_hco2plus=iq
337          mmol(igcm_hco2plus)=45.
338          count=count+1
339        endif
[2284]340        if (noms(iq).eq."hcoplus") then
341          igcm_hcoplus=iq
342          mmol(igcm_hcoplus)=29.
343          count=count+1
344        endif
[2302]345        if (noms(iq).eq."h2oplus") then
346          igcm_h2oplus=iq
347          mmol(igcm_h2oplus)=18.
348          count=count+1
349        endif
[2321]350        if (noms(iq).eq."h3oplus") then
351          igcm_h3oplus=iq
352          mmol(igcm_h3oplus)=19.
353          count=count+1
354        endif
355        if (noms(iq).eq."ohplus") then
356          igcm_ohplus=iq
357          mmol(igcm_ohplus)=17.
358          count=count+1
359        endif
[171]360        if (noms(iq).eq."elec") then
361          igcm_elec=iq
362          mmol(igcm_elec)=1./1822.89
363          count=count+1
364        endif
[38]365        if (noms(iq).eq."h2o_vap") then
366          igcm_h2o_vap=iq
367          mmol(igcm_h2o_vap)=18.
368          count=count+1
369        endif
[2312]370        if (noms(iq).eq."hdo_vap") then
371          igcm_hdo_vap=iq
[2462]372          mmol(igcm_hdo_vap)=19.
[2312]373          count=count+1
374        endif
[2461]375        if (noms(iq).eq."od") then
376          igcm_od=iq
377          mmol(igcm_od)=18.
378          count=count+1
379        endif
380        if (noms(iq).eq."d") then
381           igcm_d=iq
382           mmol(igcm_d)=2.
383           count=count+1
384        endif
385        if (noms(iq).eq."hd") then
386           igcm_hd=iq
387           mmol(igcm_hd)=3.
388           count=count+1
389        endif
390        if (noms(iq).eq."do2") then
391           igcm_do2=iq
392           mmol(igcm_do2)=34.
393           count=count+1
394        endif
395        if (noms(iq).eq."hdo2") then
396           igcm_hdo2=iq
397           mmol(igcm_hdo2)=35.
398           count=count+1
399        endif
[1617]400        if (noms(iq).eq."co2_ice") then
401          igcm_co2_ice=iq
402          mmol(igcm_co2_ice)=44.
403          count=count+1
404        endif
[38]405        if (noms(iq).eq."h2o_ice") then
406          igcm_h2o_ice=iq
407          mmol(igcm_h2o_ice)=18.
408          count=count+1
409        endif
[2312]410        if (noms(iq).eq."hdo_ice") then
411          igcm_hdo_ice=iq
[2463]412          mmol(igcm_hdo_ice)=19.
[2312]413          count=count+1
414        endif
[38]415        ! Other stuff: e.g. for simulations using co2 + neutral gaz
416        if (noms(iq).eq."Ar_N2") then
417          igcm_ar_n2=iq
418          mmol(igcm_ar_n2)=30.
419          count=count+1
420        endif
[1720]421        if (co2clouds) then
[1617]422           if (noms(iq).eq."ccnco2_mass") then
423              igcm_ccnco2_mass=iq
424              count=count+1
425           endif
426           if (noms(iq).eq."ccnco2_number") then
427              igcm_ccnco2_number=iq
428              count=count+1
429           endif
[2562]430           if (co2useh2o) then
431           if (noms(iq).eq."ccnco2_h2o_number") then
432              igcm_ccnco2_h2o_number=iq
433              count=count+1
434           endif
435           if (noms(iq).eq."ccnco2_h2o_mass_ice") then
436              igcm_ccnco2_h2o_mass_ice=iq
437              count=count+1
438           endif
439           if (noms(iq).eq."ccnco2_h2o_mass_ccn") then
440              igcm_ccnco2_h2o_mass_ccn=iq
441              count=count+1
442           endif
443           end if
[1617]444        endif
445      enddo                     ! of do iq=1,nq
446     
[38]447      ! check that we identified all tracers:
[1036]448      if (count.ne.nq) then
[38]449        write(*,*) "initracer: found only ",count," tracers"
[1036]450        write(*,*) "               expected ",nq
[38]451        do iq=1,count
452          write(*,*)'      ',iq,' ',trim(noms(iq))
453        enddo
[2302]454        call abort_physic("initracer","tracer mismatch",1)
[38]455      else
456        write(*,*) "initracer: found all expected tracers, namely:"
[1036]457        do iq=1,nq
[38]458          write(*,*)'      ',iq,' ',trim(noms(iq))
459        enddo
460      endif
461
462      ! if water cycle but iceparty=.false., there will nevertheless be
463      ! water ice at the surface (iceparty is not used anymore, but this
464      ! part is still relevant, as we want to stay compatible with the
465      ! older versions).
466      if (water.and.(igcm_h2o_ice.eq.0)) then
467        igcm_h2o_ice=igcm_h2o_vap ! so that qsurf(i_h2o_ice) is identified
468                                  ! even though there is no q(i_h2o_ice)
469      else
470       ! surface ice qsurf(i_h2o_ice) was loaded twice by phyetat0,
471       ! as qsurf(i_h2o_vap) & as qsurf(i_h2o_ice), so to be clean:
472       if (igcm_h2o_vap.ne.0) then
[1036]473         qsurf(1:ngrid,igcm_h2o_vap)=0
[38]474       endif
475      endif
476
[2312]477      ! Additional test required for HDO
478      ! We need to compute some things for H2O before HDO
479      if (hdo) then
480        if (igcm_h2o_vap.gt.igcm_hdo_vap) then
[2398]481           call abort_physic("initracer",
482     &          "Tracer H2O must be initialized before HDO",1)
[2324]483        else if ((nqfils(igcm_h2o_ice).lt.1)
484     &             .or. (nqfils(igcm_h2o_vap).lt.1)) then
485           write(*,*) "HDO must be transported as a son",
486     &                " of H2O: complete traceur.def"
[2398]487           call abort_physic("initracer","adapt your tracer.def",1)
[2324]488        else if ((igcm_hdo_ice.lt.nq-2)
489     &             .or. (igcm_hdo_vap.lt.nq-2)) then
490           write(*,*) "The isotopes (HDO) must be placed at",
491     &                " the end of the list in traceur.def"
[2398]492           call abort_physic("initracer","adapt your tracer.def",1)
[2312]493        endif
494      endif
495
[38]496c------------------------------------------------------------
[1036]497c     Initialize tracers .... (in tracer_mod)
[38]498c------------------------------------------------------------
[1005]499      ! start by setting everything to (default) zero
[1036]500      rho_q(1:nq)=0     ! tracer density (kg.m-3)
501      radius(1:nq)=0.   ! tracer particle radius (m)
502      alpha_lift(1:nq) =0.  ! tracer saltation vertical flux/horiz flux ratio (m-1)
503      alpha_devil(1:nq)=0.  ! tracer lifting coefficient by dust devils
[38]504
[1005]505
506      ! some reference values
[38]507      rho_dust=2500.  ! Mars dust density (kg.m-3)
508      rho_ice=920.    ! Water ice density (kg.m-3)
[1685]509      rho_ice_co2=1650.
510      !Mangan et al., Icarus 2017 :CO2 density = 1.72391-2.53×10−4T – 2.87×10−6T^2
[38]511      nuice_ref=0.1   ! Effective variance nueff of the
[358]512                      ! water-ice size distribution
[455]513      !!!nuice_sed=0.45   ! Sedimentation effective variance
[358]514                      ! of the water-ice size distribution
[38]515
[2494]516      nuiceco2_ref=0.2    !C.M. Effective variance nueff of the
517                          ! co2-ice size distribution
[38]518      if (doubleq) then
519c       "doubleq" technique
520c       -------------------
521c      (transport of mass and number mixing ratio)
522c       iq=1: Q mass mixing ratio, iq=2: N number mixing ratio
523
[2312]524        if( (nq.lt.2).or.(water.and.(nq.lt.4))
525     *       .or.(hdo.and.(nq.lt.6) )) then
[1036]526          write(*,*)'initracer: nq is too low : nq=', nq
[38]527          write(*,*)'water= ',water,' doubleq= ',doubleq   
528        end if
529
530        nueff_lift = 0.5
531        varian=sqrt(log(1.+nueff_lift))
532
533        rho_q(igcm_dust_mass)=rho_dust
534        rho_q(igcm_dust_number)=rho_dust
535
536c       Intermediate calcul for computing geometric mean radius r0
537c       as a function of mass and number mixing ratio Q and N
538c       (r0 = (r3n_q * Q/ N)^(1/3))
539        r3n_q = exp(-4.5*varian**2)*(3./4.)/(pi*rho_dust)
540
541c       Intermediate calcul for computing effective radius reff
542c       from geometric mean radius r0
543c       (reff = ref_r0 * r0)
544        ref_r0 = exp(2.5*varian**2)
545       
546c       lifted dust :
547c       '''''''''''
548        reff_lift = 3.0e-6 !3.e-6 !Effective radius of lifted dust (m)
549        alpha_devil(igcm_dust_mass)=9.e-9   !  dust devil lift mass coeff
550c       alpha_lift(igcm_dust_mass)=3.0e-15  !  Lifted mass coeff
[1455]551
552!! default lifting settings
553!! -- GCM: alpha_lift not zero because large-scale lifting by default
554!! -- MESOSCALE: alpha_lift zero because no lifting at all in mesoscale by default
555#ifdef MESOSCALE
556        alpha_lift(igcm_dust_mass)=0.0
557#else
[38]558        alpha_lift(igcm_dust_mass)=1.e-6 !1.e-6 !Lifted mass coeff
[1974]559        IF (dustinjection.ge.1) THEN
560                reff_lift = 3.0e-6 ! Effective radius of lifted dust (m)
561                alpha_lift(igcm_dust_mass)=(4/3.)*reff_lift*rho_dust
562     &                                          /2.4
563        ENDIF
[1455]564#endif
[38]565
566        r0_lift = reff_lift/ref_r0
567        alpha_devil(igcm_dust_number)=r3n_q*
568     &                        alpha_devil(igcm_dust_mass)/r0_lift**3
569        alpha_lift(igcm_dust_number)=r3n_q*
570     &                        alpha_lift(igcm_dust_mass)/r0_lift**3
571
572        radius(igcm_dust_mass) = reff_lift
573        radius(igcm_dust_number) = reff_lift
574
575        write(*,*) "initracer: doubleq_param reff_lift:", reff_lift
576        write(*,*) "initracer: doubleq_param nueff_lift:", nueff_lift
577        write(*,*) "initracer: doubleq_param alpha_lift:",
578     &    alpha_lift(igcm_dust_mass)
[1974]579!c ----------------------------------------------------------------------
580!c rocket dust storm scheme
581!c lifting tracer stormdust using same distribution than
582!c normal dust
583        if (rdstorm) then
584          reff_storm=3.e-6 ! reff_lift !3.e-6
585          r0_storm=reff_storm/ref_r0
586          rho_q(igcm_stormdust_mass)=rho_dust
587          rho_q(igcm_stormdust_number)=rho_dust
588
589          alpha_devil(igcm_stormdust_mass)=9.e-9   ! dust devil lift mass coeff
590          alpha_lift(igcm_stormdust_mass)=4./3./2.4*reff_storm*rho_dust
591
592          write(*,*) 'alpha_lift(rds):',alpha_lift(igcm_stormdust_mass)
593 
594          alpha_devil(igcm_stormdust_number)=r3n_q*
595     &                      alpha_devil(igcm_stormdust_mass)/r0_storm**3
596          alpha_lift(igcm_stormdust_number)=r3n_q*
597     &                       alpha_lift(igcm_stormdust_mass)/r0_storm**3
598 
599          radius(igcm_stormdust_mass) = reff_storm
600          radius(igcm_stormdust_number) = reff_storm
601        end if !(rdstorm)
602!c ----------------------------------------------------------------------
[2199]603!c slope wind scheme
604!c you need a radius value for topdust to active its sedimentation
605!c we take the same value as for the normal dust
606        if (slpwind) then
607          rho_q(igcm_topdust_mass)=rho_dust
608          rho_q(igcm_topdust_number)=rho_dust
609          radius(igcm_topdust_mass) = 3.e-6
610          radius(igcm_topdust_number) = 3.e-6
611        end if !(slpwind)
612!c ----------------------------------------------------------------------
[1974]613     
[38]614      else
615
[648]616       ! initialize varian, which may be used (e.g. by surfacearea)
617       ! even with conrath dust
618       nueff_lift = 0.5
619       varian=sqrt(log(1.+nueff_lift))
620
[38]621       if (dustbin.gt.1) then
622        print*,'initracer: STOP!',
623     $   ' properties of dust need to be set in initracer !!!'
[2302]624        call abort_physic("initracer","dustbin properties issue",1)
[38]625
626       else if (dustbin.eq.1) then
627
628c       This will be used for 1 dust particle size:
629c       ------------------------------------------
630        radius(igcm_dustbin(1))=3.e-6
631        alpha_lift(igcm_dustbin(1))=0.0e-6
632        alpha_devil(igcm_dustbin(1))=7.65e-9
633        rho_q(igcm_dustbin(1))=rho_dust
634
635       endif
636      end if    ! (doubleq)
637
[358]638
639c     Scavenging of dust particles by H2O clouds:
640c     ------------------------------------------
641c     Initialize the two tracers used for the CCNs
642      if (water.AND.doubleq.AND.scavenging) then
643        radius(igcm_ccn_mass) = radius(igcm_dust_mass)
644        alpha_lift(igcm_ccn_mass) = 1e-30
645        alpha_devil(igcm_ccn_mass) = 1e-30
646        rho_q(igcm_ccn_mass) = rho_dust
647
648        radius(igcm_ccn_number) = radius(igcm_ccn_mass)
649        alpha_lift(igcm_ccn_number) = alpha_lift(igcm_ccn_mass)
650        alpha_devil(igcm_ccn_number) = alpha_devil(igcm_ccn_mass)
651        rho_q(igcm_ccn_number) = rho_q(igcm_ccn_mass)
652      endif ! of if (water.AND.doubleq.AND.scavenging)
653
[38]654c     Submicron dust mode:
655c     --------------------
656
657      if (submicron) then
658        radius(igcm_dust_submicron)=0.1e-6
659        rho_q(igcm_dust_submicron)=rho_dust
660        if (doubleq) then
661c         If doubleq is also active, we use the population ratio:
662          alpha_lift(igcm_dust_submicron) =
663     &      alpha_lift(igcm_dust_number)*popratio*
664     &      rho_q(igcm_dust_submicron)*4./3.*pi*
665     &      radius(igcm_dust_submicron)**3.
666          alpha_devil(igcm_dust_submicron)=1.e-30
667        else
668          alpha_lift(igcm_dust_submicron)=1e-6
669          alpha_devil(igcm_dust_submicron)=1.e-30
670        endif ! (doubleq)
671      end if  ! (submicron)
672
673c     Initialization for water vapor
674c     ------------------------------
675      if(water) then
676         radius(igcm_h2o_vap)=0.
677         alpha_lift(igcm_h2o_vap) =0.
678         alpha_devil(igcm_h2o_vap)=0.
[1036]679         if(water.and.(nq.ge.2)) then
[38]680           radius(igcm_h2o_ice)=3.e-6
681           rho_q(igcm_h2o_ice)=rho_ice
682           alpha_lift(igcm_h2o_ice) =0.
683           alpha_devil(igcm_h2o_ice)=0.
[1036]684         elseif(water.and.(nq.lt.2)) then
685            write(*,*) 'nq is too low : nq=', nq
[38]686            write(*,*) 'water= ',water
687         endif
688
689      end if  ! (water)
[2312]690
691c     Initialization for hdo vapor
692c     ------------------------------
693      if (hdo) then
694         radius(igcm_hdo_vap)=0.
695         alpha_lift(igcm_hdo_vap) =0.
696         alpha_devil(igcm_hdo_vap)=0.
697         if(water.and.(nq.ge.2)) then
698           radius(igcm_hdo_ice)=3.e-6
699           rho_q(igcm_hdo_ice)=rho_ice
700           alpha_lift(igcm_hdo_ice) =0.
701           alpha_devil(igcm_hdo_ice)=0.
702         elseif(hdo.and.(nq.lt.6)) then
703            write(*,*) 'nq is too low : nq=', nq
704            write(*,*) 'hdo= ',hdo
705         endif
706
707      end if  ! (hdo)
708
[1617]709     
710! Initialisation for CO2 clouds
711      if (co2clouds ) then
712        radius(igcm_ccnco2_mass) = radius(igcm_dust_mass)
713        alpha_lift(igcm_ccnco2_mass) = 1e-30
714        alpha_devil(igcm_ccnco2_mass) = 1e-30
715        rho_q(igcm_ccnco2_mass) = rho_dust
716        radius(igcm_ccnco2_number) = radius(igcm_ccnco2_mass)
717        alpha_lift(igcm_ccnco2_number) = alpha_lift(igcm_ccnco2_mass)
718        alpha_devil(igcm_ccnco2_number) = alpha_devil(igcm_ccnco2_mass)
719        rho_q(igcm_ccnco2_number) = rho_q(igcm_ccnco2_mass)
720     
721        radius(igcm_co2)=0.
722        alpha_lift(igcm_co2) =0.
723        alpha_devil(igcm_co2)=0.
[1685]724        radius(igcm_co2_ice)=1.e-8
[1617]725        rho_q(igcm_co2_ice)=rho_ice_co2
726        alpha_lift(igcm_co2_ice) =0.
727        alpha_devil(igcm_co2_ice)=0.
[38]728
[1617]729      endif
730     
[38]731c     Output for records:
732c     ~~~~~~~~~~~~~~~~~~
733      write(*,*)
734      Write(*,*) '******** initracer : dust transport parameters :'
735      write(*,*) 'alpha_lift = ', alpha_lift
736      write(*,*) 'alpha_devil = ', alpha_devil
737      write(*,*) 'radius  = ', radius
738      if(doubleq) then
739        write(*,*) 'reff_lift (um) =  ', reff_lift
740        write(*,*) 'size distribution variance  = ', varian
741        write(*,*) 'r3n_q , ref_r0 : ', r3n_q , ref_r0
742      end if
743
744!
745!     some extra (possibly redundant) sanity checks for tracers:
746!     ---------------------------------------------------------
747
748       if (doubleq) then
749       ! verify that we indeed have dust_mass and dust_number tracers
750         if (igcm_dust_mass.eq.0) then
751           write(*,*) "initracer: error !!"
752           write(*,*) "  cannot use doubleq option without ",
753     &                "a dust_mass tracer !"
[2302]754           call abort_physic("initracer","doubleq issue",1)
[38]755         endif
756         if (igcm_dust_number.eq.0) then
757           write(*,*) "initracer: error !!"
758           write(*,*) "  cannot use doubleq option without ",
759     &                "a dust_number tracer !"
[2302]760           call abort_physic("initracer","doubleq issue",1)
[38]761         endif
762       endif
763
764       if ((.not.doubleq).and.(dustbin.gt.0)) then
765       ! verify that we indeed have 'dustbin' dust tracers
766         count=0
767         do iq=1,dustbin
768           if (igcm_dustbin(iq).ne.0) then
769             count=count+1
770           endif
771         enddo
772         if (count.ne.dustbin) then
773           write(*,*) "initracer: error !!"
[2302]774           write(*,*) "  dustbin is set to ",dustbin,
[38]775     &                " but we only have the following dust tracers:"
776           do iq=1,count
777             write(*,*)"   ",trim(noms(igcm_dustbin(iq)))
778           enddo
[2302]779           call abort_physic("initracer","dustbin issue",1)
[38]780         endif
781       endif
782
783       if (water) then
784       ! verify that we indeed have h2o_vap and h2o_ice tracers
785         if (igcm_h2o_vap.eq.0) then
786           write(*,*) "initracer: error !!"
787           write(*,*) "  cannot use water option without ",
788     &                "an h2o_vap tracer !"
[2302]789           call abort_physic("initracer","water cycle issue",1)
[38]790         endif
791         if (igcm_h2o_ice.eq.0) then
792           write(*,*) "initracer: error !!"
793           write(*,*) "  cannot use water option without ",
794     &                "an h2o_ice tracer !"
[2302]795           call abort_physic("initracer","water cycle issue",1)
[38]796         endif
797       endif
798
[2312]799       if (hdo) then
800       ! verify that we indeed have hdo_vap and hdo_ice tracers
801         if (igcm_hdo_vap.eq.0) then
802           write(*,*) "initracer: error !!"
803           write(*,*) "  cannot use hdo option without ",
804     &                "an hdo_vap tracer !"
[2398]805           call abort_physic("initracer","hdo cycle issue",1)
[2312]806         endif
807         if (igcm_hdo_ice.eq.0) then
808           write(*,*) "initracer: error !!"
809           write(*,*) "  cannot use hdo option without ",
810     &                "an hdo_ice tracer !"
[2398]811           call abort_physic("initracer","hdo cycle issue",1)
[2312]812         endif
813       endif
814
815
[1617]816       if (co2clouds) then
817          !verify that we have co2_ice and co2 tracers
818          if (igcm_co2 .eq. 0) then
819             write(*,*) "initracer: error !!"
820             write(*,*) "  cannot use co2 clouds option without ",
821     &            "a co2 tracer !"
[2302]822             call abort_physic("initracer","co2 clouds issue",1)
[1617]823          endif
824          if (igcm_co2_ice .eq. 0) then
825             write(*,*) "initracer: error !!"
826             write(*,*) "  cannot use co2 clouds option without ",
827     &            "a co2_ice tracer !"
[2302]828             call abort_physic("initracer","co2 clouds issue",1)
[1617]829          endif
830       endif
[1974]831 
832       if (rdstorm) then
833       ! verify that we indeed have stormdust_mass and stormdust_number tracers
834         if (igcm_stormdust_mass.eq.0) then
835           write(*,*) "initracer: error !!"
836           write(*,*) "  cannot use rdstorm option without ",
837     &                "a stormdust_mass tracer !"
[2302]838           call abort_physic("initracer","rdstorm issue",1)
[1974]839         endif
840         if (igcm_stormdust_number.eq.0) then
841           write(*,*) "initracer: error !!"
842           write(*,*) "  cannot use rdstorm option without ",
843     &                "a stormdust_number tracer !"
[2302]844           call abort_physic("initracer","rdstorm issue",1)
[1974]845         endif
846       endif
[2199]847
848       if (slpwind) then
849       ! verify that we indeed have topdust_mass and topdust_number tracers
850         if (igcm_topdust_mass.eq.0) then
851           write(*,*) "initracer: error !!"
852           write(*,*) "  cannot use slpwind option without ",
853     &                "a topdust_mass tracer !"
[2302]854           call abort_physic("initracer","slpwind issue",1)
[2199]855         endif
856         if (igcm_topdust_number.eq.0) then
857           write(*,*) "initracer: error !!"
858           write(*,*) "  cannot use slpwind option without ",
859     &                "a topdust_number tracer !"
[2302]860           call abort_physic("initracer","slpwind issue",1)
[2199]861         endif
862       endif
[1974]863     
[1380]864       if (callnlte) then ! NLTE requirements
865         if (nltemodel.ge.1) then
866           ! check that co2, co, o and n2 tracers are available
867           if (igcm_co2.eq.0) then
868             write(*,*) "initracer: error !!"
869             write(*,*) "  with nltemodel>0, we need the co2 tracer!"
[2302]870             call abort_physic("initracer","missing co2 tracer",1)
[1380]871           endif
872           if (igcm_co.eq.0) then
873             write(*,*) "initracer: error !!"
874             write(*,*) "  with nltemodel>0, we need the co tracer!"
[2302]875             call abort_physic("initracer","missing co tracer",1)
[1380]876           endif
877           if (igcm_o.eq.0) then
878             write(*,*) "initracer: error !!"
879             write(*,*) "  with nltemodel>0, we need the o tracer!"
[2302]880             call abort_physic("initracer","missing o tracer",1)
[1380]881           endif
882           if (igcm_n2.eq.0) then
883             write(*,*) "initracer: error !!"
884             write(*,*) "  with nltemodel>0, we need the n2 tracer!"
[2302]885             call abort_physic("initracer","missing n2 tracer",1)
[1380]886           endif
887         endif
888       endif
889
[358]890       if (scavenging) then
891       ! verify that we indeed have ccn_mass and ccn_number tracers
[1617]892         if (igcm_ccn_mass.eq.0 .and. igcm_ccnco2_mass.eq.0) then
[358]893           write(*,*) "initracer: error !!"
894           write(*,*) "  cannot use scavenging option without ",
[1617]895     &                "a ccn_mass or ccnco2_mass tracer !"
[2302]896             call abort_physic("initracer","scavenging issue",1)
[358]897         endif
[1617]898         if (igcm_ccn_number.eq.0 .and. igcm_ccnco2_number.eq.0 ) then
[358]899           write(*,*) "initracer: error !!"
900           write(*,*) "  cannot use scavenging option without ",
[1617]901     &                "a ccn_number or ccnco2_number tracer !"
[2302]902             call abort_physic("initracer","scavenging issue",1)
[358]903         endif
904       endif ! of if (scavenging)
905
[38]906       if (photochem .or. callthermos) then
907       ! verify that we indeed have the chemistry tracers
908         if (igcm_co2.eq.0) then
909           write(*,*) "initracer: error !!"
910           write(*,*) "  cannot use chemistry option without ",
911     &                "a co2 tracer !"
[2302]912           call abort_physic("initracer","missing co2 tracer",1)
[38]913         endif
914         if (igcm_co.eq.0) then
915           write(*,*) "initracer: error !!"
916           write(*,*) "  cannot use chemistry option without ",
917     &                "a co tracer !"
[2302]918           call abort_physic("initracer","missing co tracer",1)
[38]919         endif
920         if (igcm_o.eq.0) then
921           write(*,*) "initracer: error !!"
922           write(*,*) "  cannot use chemistry option without ",
923     &                "a o tracer !"
[2302]924           call abort_physic("initracer","missing o tracer",1)
[38]925         endif
926         if (igcm_o1d.eq.0) then
927           write(*,*) "initracer: error !!"
928           write(*,*) "  cannot use chemistry option without ",
929     &                "a o1d tracer !"
[2302]930           call abort_physic("initracer","missing o1d tracer",1)
[38]931         endif
932         if (igcm_o2.eq.0) then
933           write(*,*) "initracer: error !!"
934           write(*,*) "  cannot use chemistry option without ",
935     &                "an o2 tracer !"
[2302]936           call abort_physic("initracer","missing o2 tracer",1)
[38]937         endif
938         if (igcm_o3.eq.0) then
939           write(*,*) "initracer: error !!"
940           write(*,*) "  cannot use chemistry option without ",
941     &                "an o3 tracer !"
[2302]942           call abort_physic("initracer","missing o3 tracer",1)
[38]943         endif
944         if (igcm_h.eq.0) then
945           write(*,*) "initracer: error !!"
946           write(*,*) "  cannot use chemistry option without ",
947     &                "a h tracer !"
[2302]948           call abort_physic("initracer","missing h tracer",1)
[38]949         endif
950         if (igcm_h2.eq.0) then
951           write(*,*) "initracer: error !!"
952           write(*,*) "  cannot use chemistry option without ",
953     &                "a h2 tracer !"
[2302]954           call abort_physic("initracer","missing h2 tracer",1)
[38]955         endif
956         if (igcm_oh.eq.0) then
957           write(*,*) "initracer: error !!"
958           write(*,*) "  cannot use chemistry option without ",
959     &                "an oh tracer !"
[2302]960           call abort_physic("initracer","missing oh tracer",1)
[38]961         endif
962         if (igcm_ho2.eq.0) then
963           write(*,*) "initracer: error !!"
964           write(*,*) "  cannot use chemistry option without ",
965     &                "a ho2 tracer !"
[2302]966           call abort_physic("initracer","missing ho2 tracer",1)
[1720]967      endif
[38]968         if (igcm_h2o2.eq.0) then
969           write(*,*) "initracer: error !!"
970           write(*,*) "  cannot use chemistry option without ",
971     &                "a h2o2 tracer !"
[2302]972           call abort_physic("initracer","missing h2o2 tracer",1)
[38]973         endif
974         if (igcm_n2.eq.0) then
975           write(*,*) "initracer: error !!"
976           write(*,*) "  cannot use chemistry option without ",
977     &                "a n2 tracer !"
[2302]978           call abort_physic("initracer","missing n2 tracer",1)
[38]979         endif
980         if (igcm_ar.eq.0) then
981           write(*,*) "initracer: error !!"
982           write(*,*) "  cannot use chemistry option without ",
983     &                "an ar tracer !"
[2302]984           call abort_physic("initracer","missing ar tracer",1)
[38]985         endif
986       endif ! of if (photochem .or. callthermos)
987
988      end
Note: See TracBrowser for help on using the repository browser.