source: trunk/LMDZ.GENERIC/libf/phystd/initracer.F90 @ 3058

Last change on this file since 3058 was 2804, checked in by aslmd, 2 years ago

This is the RGCS scheme

Now we can use this scheme to take into account radiative effect of a Generic Condensable Specie. One needs to add the keyword is_rgcs=1 to the ModernTraceur?.def to activate the scheme + specify the number of aerogeneric specie in callphys.def. Also, one needs to go into suaer_corrk.F90 to add the Mie scattering file into the code for the specie you want. Currently, MnS, Fe, Cr, Mg2SiO4 are available. This has been tested on a Hot Jupiter case, might have issue in different configuration. Also, the numerical stability of the scheme is not guaranted when using strong scatterers such as silicates. Physical oscillations in the evaporation/condensation can occur, and need fixing. This will be for another commit.

File size: 24.2 KB
Line 
1      SUBROUTINE initracer(ngrid,nq)
2
3      use surfdat_h, ONLY: dryness
4      USE tracer_h
5      USE callkeys_mod, only: water
6      USE recombin_corrk_mod, ONLY: ini_recombin
7      USE mod_phys_lmdz_para, only: is_master, bcast
8      use generic_cloud_common_h
9      IMPLICIT NONE
10!=======================================================================
11!   subject:
12!   --------
13!   Initialization related to tracer
14!   (transported dust, water, chemical species, ice...)
15!
16!   Name of the tracer
17!
18!   Test of dimension :
19!   Initialize COMMON tracer in tracer.h, using tracer names provided
20!   by the argument nametrac
21!
22!   author: F.Forget
23!   ------
24!            Ehouarn Millour (oct. 2008) identify tracers by their names
25!            Y Jaziri & J. Vatant d'Ollone (2020) : Modern traceur.def
26!            L Teinturier (2022): Tracer names are now read here instead of
27!                                  inside interfaces
28!=======================================================================
29
30      integer,intent(in) :: ngrid,nq
31
32      character(len=500) :: tracline   ! to read traceur.def lines
33      ! character(len=30)  :: txt        ! to store some text
34      integer :: blank      !to store the index of 1st blank when reading tracers names
35      integer iq,ig,count,ierr
36
37!-----------------------------------------------------------------------
38!  radius(nq)      ! aerosol particle radius (m)
39!  rho_q(nq)       ! tracer densities (kg.m-3)
40!  qext(nq)        ! Single Scat. Extinction coeff at 0.67 um
41!  alpha_lift(nq)  ! saltation vertical flux/horiz flux ratio (m-1)
42!  alpha_devil(nq) ! lifting coeeficient by dust devil
43!  rho_dust          ! Mars dust density
44!  rho_ice           ! Water ice density
45!  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
46!  varian            ! Characteristic variance of log-normal distribution
47!-----------------------------------------------------------------------
48
49      if (is_master) then ! only the master proc/thread needs do this
50
51        moderntracdef=.false. ! For modern traceur.def (default false, old type)
52
53        open(407, form = 'formatted', status = 'old', &
54             file = 'traceur.def', iostat=ierr)
55        if (ierr /=0) then
56          ! call abort_physic('initracer',&
57          ! 'Problem in opening traceur.def',1)
58          print*,'Problem in opening traceur.def'
59          stop
60        end if
61!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
62        READ(407,'(A)') tracline
63        IF (trim(tracline).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
64          READ(tracline,*) nqtot ! Try standard traceur.def
65        ELSE
66         moderntracdef = .true.
67         DO
68           READ(407,'(A)',iostat=ierr) tracline
69           IF (ierr==0) THEN
70             IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
71               READ(tracline,*) nqtot
72               ! Temporary not implemented solution
73               if (nqtot/=nq) then
74      !            call abort_physic('initracer','Different number of '// &
75      ! 'tracers in dynamics and physics not managed yet',1)
76                print*,'!= nbr oftracers in dynamics and physics not managed yet'
77                stop
78               endif
79               EXIT
80             ENDIF
81           ELSE ! If pb, or if reached EOF without having found nqtot
82      !        call abort_physic('initracer','Unable to read numbers '// &
83      !  'of tracers in traceur.def',1)
84            print*,"unable to read numbers of tracer in tracer.def"
85            stop
86           ENDIF
87         ENDDO
88        ENDIF ! if modern or standard traceur.def
89       
90       endif ! of if (is_master)
91       
92       ! share the information with other procs/threads (if any)
93       CALL bcast(nqtot)
94       CALL bcast(moderntracdef)
95       
96!! -----------------------------------------------------------------------
97       !! For the moment number of tracers in dynamics and physics are equal
98       nqtot=nq
99       !! we allocate once for all arrays in common in tracer_h.F90
100       !! (supposedly those are not used before call to initracer)
101       IF (.NOT.ALLOCATED(noms))           ALLOCATE(noms(nq))
102       IF (.NOT.ALLOCATED(mmol))           ALLOCATE(mmol(nq))
103       IF (.NOT.ALLOCATED(aki))            ALLOCATE(aki(nqtot))
104       IF (.NOT.ALLOCATED(cpi))            ALLOCATE(cpi(nqtot))
105       IF (.NOT.ALLOCATED(radius))         ALLOCATE(radius(nq))
106       IF (.NOT.ALLOCATED(rho_q))          ALLOCATE(rho_q(nq))
107       IF (.NOT.ALLOCATED(qext))           ALLOCATE(qext(nq))
108       IF (.NOT.ALLOCATED(alpha_lift))     ALLOCATE(alpha_lift(nq))
109       IF (.NOT.ALLOCATED(alpha_devil))    ALLOCATE(alpha_devil(nq))
110       IF (.NOT.ALLOCATED(qextrhor))       ALLOCATE(qextrhor(nq))
111       IF (.NOT.ALLOCATED(igcm_dustbin))   ALLOCATE(igcm_dustbin(nq))
112       IF (.NOT.ALLOCATED(is_chim))        ALLOCATE(is_chim(nqtot))
113       IF (.NOT.ALLOCATED(is_rad))         ALLOCATE(is_rad(nqtot))
114       IF (.NOT.ALLOCATED(is_recomb))      ALLOCATE(is_recomb(nqtot))
115       IF (.NOT.ALLOCATED(is_recomb_qset)) THEN
116         ALLOCATE(is_recomb_qset(nqtot))
117       ENDIF
118       IF (.NOT.ALLOCATED(is_recomb_qotf)) THEN
119         ALLOCATE(is_recomb_qotf(nqtot))
120       ENDIF
121       IF (.NOT. allocated(is_condensable)) allocate(is_condensable(nq)) !LT
122       IF (.NOT. allocated(is_rgcs)) allocate(is_rgcs(nq)) !LT
123       IF (.NOT. allocated(constants_mass)) allocate(constants_mass(nq))
124       IF (.NOT. allocated(constants_delta_vapH)) allocate(constants_delta_vapH(nq))
125       IF (.NOT. allocated(constants_Tref)) allocate(constants_Tref(nq))
126       IF (.NOT. allocated(constants_Pref)) allocate(constants_Pref(nq))
127       IF (.NOT. allocated(constants_epsi_generic)) allocate(constants_epsi_generic(nq))
128       IF (.NOT. allocated(constants_RLVTT_generic)) allocate(constants_RLVTT_generic(nq))
129       IF (.NOT. allocated(constants_metallicity_coeff)) allocate(constants_metallicity_coeff(nq))
130
131       !! initialization
132       alpha_lift(:)     = 0.
133       alpha_devil(:)    = 0.
134       mmol(:)           = 0.
135       aki(:)            = 0.
136       cpi(:)            = 0.
137       is_chim(:)        = 0
138       is_rad(:)         = 0
139       is_recomb(:)      = 0
140       is_recomb_qset(:) = 0
141       is_recomb_qotf(:) = 0
142       
143       ! Added by JVO 2017 : these arrays are handled later
144       ! -> initialization is the least we can do, please !!!
145       radius(:)=0.
146       qext(:)=0.
147
148       ! For condensable tracers, by Lucas Teinturier and Noé Clément (2022)
149
150       is_condensable(:)= 0
151       is_rgcs(:) = 0
152       constants_mass(:)=0
153       constants_delta_vapH(:)=0
154       constants_Tref(:)=0
155       constants_Pref(:)=0
156       constants_epsi_generic(:)=0
157       constants_RLVTT_generic(:)=0
158       constants_metallicity_coeff(:)=0
159
160       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
161
162
163! Initialization: Read tracers names from traceur.def
164        do iq=1,nq
165          if (is_master) read(407,'(A)') tracline
166          call bcast(tracline)
167          blank = index(tracline,' ') ! Find position of 1st blank in tracline
168          noms(iq) = tracline(1:blank) !ensure that in Modern-trac case, noms is actually the name and not all properties
169        enddo
170
171! Identify tracers by their names: (and set corresponding values of mmol)
172      ! 0. initialize tracer indexes to zero:
173      ! NB: igcm_* indexes are commons in 'tracer.h'
174      do iq=1,nq
175        igcm_dustbin(iq)=0
176      enddo
177      igcm_dust_mass=0
178      igcm_dust_number=0
179      igcm_h2o_vap=0
180      igcm_h2o_ice=0
181      igcm_co2=0
182      igcm_co=0
183      igcm_o=0
184      igcm_o1d=0
185      igcm_o2=0
186      igcm_o3=0
187      igcm_h=0
188      igcm_h2=0
189      igcm_oh=0
190      igcm_ho2=0
191      igcm_h2o2=0
192      igcm_n2=0
193      igcm_n=0
194      igcm_n2d=0
195      igcm_no=0
196      igcm_no2=0
197      igcm_ar=0
198      igcm_ar_n2=0
199      igcm_co2_ice=0
200
201      igcm_ch4=0
202      igcm_ch3=0
203      igcm_ch=0
204      igcm_3ch2=0
205      igcm_1ch2=0
206      igcm_cho=0
207      igcm_ch2o=0
208      igcm_ch3o=0
209      igcm_c=0
210      igcm_c2=0
211      igcm_c2h=0
212      igcm_c2h2=0
213      igcm_c2h3=0
214      igcm_c2h4=0
215      igcm_c2h6=0
216      igcm_ch2co=0
217      igcm_ch3co=0
218      igcm_hcaer=0
219
220      ! 1. find dust tracers
221      count=0
222
223      ! 2. find chemistry and water tracers
224      do iq=1,nq
225        if (noms(iq).eq."co2") then
226          igcm_co2=iq
227          mmol(igcm_co2)=44.
228          count=count+1
229!          write(*,*) 'co2: count=',count
230        endif
231        if (noms(iq).eq."co2_ice") then
232          igcm_co2_ice=iq
233          mmol(igcm_co2_ice)=44.
234          count=count+1
235!          write(*,*) 'co2_ice: count=',count
236        endif
237        if (noms(iq).eq."h2o_vap") then
238          igcm_h2o_vap=iq
239          mmol(igcm_h2o_vap)=18.
240          count=count+1
241!          write(*,*) 'h2o_vap: count=',count
242        endif
243        if (noms(iq).eq."h2o_ice") then
244          igcm_h2o_ice=iq
245          mmol(igcm_h2o_ice)=18.
246          count=count+1
247!          write(*,*) 'h2o_ice: count=',count
248        endif
249        if (noms(iq).eq."co") then
250          igcm_co=iq
251          mmol(igcm_co)=28.
252          count=count+1
253        endif
254        if (noms(iq).eq."o") then
255          igcm_o=iq
256          mmol(igcm_o)=16.
257          count=count+1
258        endif
259        if (noms(iq).eq."o1d") then
260          igcm_o1d=iq
261          mmol(igcm_o1d)=16.
262          count=count+1
263        endif
264        if (noms(iq).eq."o2") then
265          igcm_o2=iq
266          mmol(igcm_o2)=32.
267          count=count+1
268        endif
269        if (noms(iq).eq."o3") then
270          igcm_o3=iq
271          mmol(igcm_o3)=48.
272          count=count+1
273        endif
274        if (noms(iq).eq."h") then
275          igcm_h=iq
276          mmol(igcm_h)=1.
277          count=count+1
278        endif
279        if (noms(iq).eq."h2") then
280          igcm_h2=iq
281          mmol(igcm_h2)=2.
282          count=count+1
283        endif
284        if (noms(iq).eq."oh") then
285          igcm_oh=iq
286          mmol(igcm_oh)=17.
287          count=count+1
288        endif
289        if (noms(iq).eq."ho2") then
290          igcm_ho2=iq
291          mmol(igcm_ho2)=33.
292          count=count+1
293        endif
294        if (noms(iq).eq."h2o2") then
295          igcm_h2o2=iq
296          mmol(igcm_h2o2)=34.
297          count=count+1
298        endif
299        if (noms(iq).eq."n2") then
300          igcm_n2=iq
301          mmol(igcm_n2)=28.
302          count=count+1
303        endif
304        if (noms(iq).eq."ch4") then
305          igcm_ch4=iq
306          mmol(igcm_ch4)=16.
307          count=count+1
308        endif
309        if (noms(iq).eq."ar") then
310          igcm_ar=iq
311          mmol(igcm_ar)=40.
312          count=count+1
313        endif
314        if (noms(iq).eq."n") then
315          igcm_n=iq
316          mmol(igcm_n)=14.
317          count=count+1
318        endif
319        if (noms(iq).eq."no") then
320          igcm_no=iq
321          mmol(igcm_no)=30.
322          count=count+1
323        endif
324        if (noms(iq).eq."no2") then
325          igcm_no2=iq
326          mmol(igcm_no2)=46.
327          count=count+1
328        endif
329        if (noms(iq).eq."n2d") then
330          igcm_n2d=iq
331          mmol(igcm_n2d)=28.
332          count=count+1
333        endif
334        if (noms(iq).eq."ch3") then
335          igcm_ch3=iq
336          mmol(igcm_ch3)=15.
337          count=count+1
338        endif
339        if (noms(iq).eq."ch") then
340          igcm_ch=iq
341          mmol(igcm_ch)=13.
342          count=count+1
343        endif
344        if (noms(iq).eq."3ch2") then
345          igcm_3ch2=iq
346          mmol(igcm_3ch2)=14.
347          count=count+1
348        endif
349        if (noms(iq).eq."1ch2") then
350          igcm_1ch2=iq
351          mmol(igcm_1ch2)=14.
352          count=count+1
353        endif
354        if (noms(iq).eq."cho") then
355          igcm_cho=iq
356          mmol(igcm_cho)=29.
357          count=count+1
358        endif
359        if (noms(iq).eq."ch2o") then
360          igcm_ch2o=iq
361          mmol(igcm_ch2o)=30.
362          count=count+1
363        endif
364        if (noms(iq).eq."ch3o") then
365          igcm_ch3o=iq
366          mmol(igcm_ch3o)=31.
367          count=count+1
368        endif
369        if (noms(iq).eq."c") then
370          igcm_c=iq
371          mmol(igcm_c)=12.
372          count=count+1
373        endif
374        if (noms(iq).eq."c2") then
375          igcm_c2=iq
376          mmol(igcm_c2)=24.
377          count=count+1
378        endif
379        if (noms(iq).eq."c2h") then
380          igcm_c2h=iq
381          mmol(igcm_c2h)=25.
382          count=count+1
383        endif
384        if (noms(iq).eq."c2h2") then
385          igcm_c2h2=iq
386          mmol(igcm_c2h2)=26.
387          count=count+1
388        endif
389        if (noms(iq).eq."c2h3") then
390          igcm_c2h3=iq
391          mmol(igcm_c2h3)=27.
392          count=count+1
393        endif
394        if (noms(iq).eq."c2h4") then
395          igcm_c2h4=iq
396          mmol(igcm_c2h4)=28.
397          count=count+1
398        endif
399        if (noms(iq).eq."c2h6") then
400          igcm_c2h6=iq
401          mmol(igcm_c2h6)=30.
402          count=count+1
403        endif
404        if (noms(iq).eq."ch2co") then
405          igcm_ch2co=iq
406          mmol(igcm_ch2co)=42.
407          count=count+1
408        endif
409        if (noms(iq).eq."ch3co") then
410          igcm_ch3co=iq
411          mmol(igcm_ch3co)=43.
412          count=count+1
413        endif
414        if (noms(iq).eq."hcaer") then
415          igcm_hcaer=iq
416          mmol(igcm_hcaer)=50.
417          count=count+1
418        endif
419
420      enddo ! of do iq=1,nq
421
422      ! 3. find condensable traceurs different from h2o and co2
423      do iq=1,nq
424        if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
425          count=count+1
426        endif
427        if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
428          count=count+1
429        endif
430
431      enddo ! of do iq=1,nq
432     
433      ! check that we identified all tracers:
434      if (count.ne.nq) then
435        write(*,*) "initracer: found only ",count," tracers"
436        write(*,*) "               expected ",nq
437        do iq=1,count
438          write(*,*)'      ',iq,' ',trim(noms(iq))
439        enddo
440!        stop
441      else
442        write(*,*) "initracer: found all expected tracers, namely:"
443        do iq=1,nq
444          write(*,*)'      ',iq,' ',trim(noms(iq))
445        enddo
446      endif
447
448      ! Get data of tracers. Need to rewind traceur.def first
449      if (is_master) then
450       rewind(407)
451       do
452        read(407,'(A)') tracline
453        if (index(tracline,"#") .ne. 1) then
454          exit
455        endif
456       enddo
457      endif
458      do iq=1,nqtot
459        if (is_master) read(407,'(A)') tracline
460        call bcast(tracline)
461        call get_tracdat(iq, tracline)
462      enddo
463
464      if (is_master) close(407)
465
466      ! Get specific data of condensable tracers
467      do iq=1,nq
468        if((is_condensable(iq)==1)) then
469                write(*,*) "There is a specie which is condensable, for generic condensation : ", noms(iq)
470                write(*,*) 'looking specie parameters for : ',noms(iq)(1:len(trim(noms(iq)))-4)
471                call specie_parameters_table(noms(iq)(1:len(trim(noms(iq)))-4))
472                constants_mass(iq)=m
473                constants_delta_vapH(iq)=delta_vapH
474                constants_Tref(iq)=Tref
475                constants_Pref(iq)=Pref
476                constants_epsi_generic(iq)=epsi_generic
477                constants_RLVTT_generic(iq)=RLVTT_generic
478                constants_metallicity_coeff(iq)=metallicity_coeff
479        else
480                write(*,*) "This tracer is not condensable, for generic condensation :  : ", noms(iq)
481                write(*,*) "We keep condensable constants at zero"
482        endif !(is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
483      enddo ! iq=1,nq
484
485      ! Calculate number of species in the chemistry
486      nesp = sum(is_chim)
487      write(*,*) 'Number of species in the chemistry nesp = ',nesp
488
489      ! Calculate number of generic tracers
490      ngt = sum(is_condensable)
491      write(*,*) 'Number of generic tracer is  ngt = ',ngt
492
493      ! Calculate number of radiative generic condensable species
494      n_rgcs = sum(is_rgcs)
495      write(*,*)'Number of Radiative Generic Condensable Species is n_rgcs = ',n_rgcs
496      if (n_rgcs> ngt/2) then
497        write(*,*) 'You have more Radiative Generic Condensable Species than Generic Condensable Species'
498        write(*,*)'This is not possible: check your Modern traceur.def'
499        call abort_physic("initracer, issue with # of RGCS and GCS")
500      endif
501
502!     Processing modern traceur options
503      if(moderntracdef) then
504        call ini_recombin
505      endif
506
507!------------------------------------------------------------
508!     Initialisation tracers ....
509!------------------------------------------------------------
510      ! rho_q(1:nq)=0
511
512      rho_dust=2500.  ! Mars dust density (kg.m-3)
513      rho_ice=920.    ! Water ice density (kg.m-3)
514      rho_co2=1620.   ! CO2 ice density (kg.m-3)
515
516!     Initialization for water vapor
517!     ------------------------------
518      if(water) then
519        radius(igcm_h2o_vap)=0.
520        Qext(igcm_h2o_vap)=0.
521        alpha_lift(igcm_h2o_vap) =0.
522        alpha_devil(igcm_h2o_vap)=0.
523        qextrhor(igcm_h2o_vap)= 0.
524
525         !! this is defined in surfdat_h.F90
526         IF (.not.ALLOCATED(dryness)) ALLOCATE(dryness(ngrid))
527
528         do ig=1,ngrid
529           dryness(ig) = 1.
530         enddo
531
532
533           radius(igcm_h2o_ice)=3.e-6
534           rho_q(igcm_h2o_ice)=rho_ice
535           Qext(igcm_h2o_ice)=0.
536!           alpha_lift(igcm_h2o_ice) =0.
537!           alpha_devil(igcm_h2o_ice)=0.
538           qextrhor(igcm_h2o_ice)= (3./4.)*Qext(igcm_h2o_ice) &
539            / (rho_ice*radius(igcm_h2o_ice))
540
541
542      end if  ! (water)
543
544
545!
546!     some extra (possibly redundant) sanity checks for tracers:
547!     ---------------------------------------------------------
548       if (water) then
549       ! verify that we indeed have h2o_vap and h2o_ice tracers
550         if (igcm_h2o_vap.eq.0) then
551           write(*,*) "initracer: error !!"
552           write(*,*) "  cannot use water option without ",&
553                     "an h2o_vap tracer !"
554           stop
555         endif
556         if (igcm_h2o_ice.eq.0) then
557           write(*,*) "initracer: error !!"
558           write(*,*) "  cannot use water option without ",&
559                     "an h2o_ice tracer !"
560           stop
561         endif
562       endif
563
564
565!     Output for records:
566!     ~~~~~~~~~~~~~~~~~~
567      write(*,*)
568      Write(*,*) '******** initracer : dust transport parameters :'
569      write(*,*) 'alpha_lift = ', alpha_lift
570      write(*,*) 'alpha_devil = ', alpha_devil
571      write(*,*) 'radius  = ', radius
572      write(*,*) 'Qext  = ', qext
573      write(*,*)
574
575      contains
576
577      subroutine get_tracdat(iq,tracline)
578        !-------------------ADDING NEW OPTIONS-------------------
579        ! Duplicate if sentence for adding new options
580        ! Do not forget to add the new variables in tracer_h.F90
581        ! Do not forget to allocate and initialize the new variables above
582        ! Please update list of options in "LMDZ.GENERIC/deftank/traceur.def.modern"
583        !-------------------------------------------------------
584          use tracer_h
585          implicit none
586          integer,           intent(in) :: iq       ! tracer index
587          character(len=500),intent(in) :: tracline ! traceur.def lines with parameters
588 
589          read(tracline,*) noms(iq)
590          ! JVO 20 : We should add a sanity check aborting when duplicates in names !
591          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
592          ! option mmol
593          if (index(tracline,'mmol='   ) /= 0) then
594              read(tracline(index(tracline,'mmol=')+len('mmol='):),*)&
595                  mmol(iq)
596              write(*,*) ' Parameter value (traceur.def) : mmol=', &
597                  mmol(iq)
598          else
599              write(*,*) ' Parameter value (default)     : mmol=', &
600                 mmol(iq)
601          end if
602          ! option aki
603          if (index(tracline,'aki='   ) /= 0) then
604              read(tracline(index(tracline,'aki=')+len('aki='):),*) &
605                  aki(iq)
606              write(*,*) ' Parameter value (traceur.def) : aki=', &
607                  aki(iq)
608          else
609              write(*,*) ' Parameter value (default)     : aki=', &
610                  aki(iq)
611          end if
612          ! option cpi
613          if (index(tracline,'cpi='   ) /= 0) then
614              read(tracline(index(tracline,'cpi=')+len('cpi='):),*) &
615                  cpi(iq)
616              write(*,*) ' Parameter value (traceur.def) : cpi=', &
617                  cpi(iq)
618          else
619              write(*,*) ' Parameter value (default)     : cpi=', &
620                  cpi(iq)
621          end if
622          ! option is_chim
623          if (index(tracline,'is_chim='   ) /= 0) then
624          read(tracline(index(tracline,'is_chim=')+len('is_chim='):),*) &
625                  is_chim(iq)
626              write(*,*) ' Parameter value (traceur.def) : is_chim=', &
627                  is_chim(iq)
628          else
629              write(*,*) ' Parameter value (default)     : is_chim=', &
630                  is_chim(iq)
631          end if
632          ! option is_rad
633          if (index(tracline,'is_rad=') /= 0) then
634              read(tracline(index(tracline,'is_rad=') &
635              +len('is_rad='):),*) is_rad(iq)
636              write(*,*) ' Parameter value (traceur.def) : is_rad=', &
637              is_rad(iq)
638          else
639              write(*,*) ' Parameter value (default)     : is_rad=',  &
640              is_rad(iq)
641          end if
642          ! option is_recomb
643          if (index(tracline,'is_recomb=') /= 0) then
644              read(tracline(index(tracline,'is_recomb=') &
645              +len('is_recomb='):),*) is_recomb(iq)
646              write(*,*) ' Parameter value (traceur.def) : is_recomb=', &
647              is_recomb(iq)
648          else
649              write(*,*) ' Parameter value (default)     : is_recomb=', &
650              is_recomb(iq)
651          end if
652          ! option is_recomb_qset
653          if (index(tracline,'is_recomb_qset=') /= 0) then
654              read(tracline(index(tracline,'is_recomb_qset=') &
655              +len('is_recomb_qset='):),*) is_recomb_qset(iq)
656              write(*,*) ' Parameter value (traceur.def) :'// &
657              ' is_recomb_qset=', &
658              is_recomb_qset(iq)
659          else
660              write(*,*) ' Parameter value (default)     :'// &
661              ' is_recomb_qset=', &
662              is_recomb_qset(iq)
663          endif
664          ! option is_recomb_qotf
665          if (index(tracline,'is_recomb_qotf=') /= 0) then
666              read(tracline(index(tracline,'is_recomb_qotf=') &
667              +len('is_recomb_qotf='):),*) is_recomb_qotf(iq)
668              write(*,*) ' Parameter value (traceur.def) :'// &
669              ' is_recomb_qotf=', &
670              is_recomb_qotf(iq)
671          else
672              write(*,*) ' Parameter value (default)     :'// &
673              ' is_recomb_qotf=',  &
674              is_recomb_qotf(iq)
675          end if
676          !option is_condensable (LT)
677          if (index(tracline,'is_condensable=') /=0) then
678            read(tracline(index(tracline,'is_condensable=') &
679              +len('is_condensable='):),*) is_condensable(iq)
680            write(*,*) ' Parameter value (traceur.def) :'// &
681              ' is_condensable=', is_condensable(iq)       
682          else
683              write(*,*) ' Parameter value (default)     :'// &
684              ' is_condensable=', is_condensable(iq)       
685          endif
686          !option radius
687          if (index(tracline,'radius=') .ne. 0) then
688            read(tracline(index(tracline,'radius=') &
689              +len('radius='):),*) radius(iq)
690            write(*,*)'Parameter value (traceur.def) :'// &
691              "radius=",radius(iq), "m"
692          else
693            write(*,*) ' Parameter value (default)     :'// &
694            ' radius=', radius(iq)," m"     
695          endif
696          !option rho
697          if (index(tracline,'rho=') .ne. 0) then
698            read(tracline(index(tracline,'rho=') &
699              +len('rho='):),*) rho_q(iq)
700            write(*,*)'Parameter value (traceur.def) :'// &
701              "rho=",rho_q(iq)
702          else
703            write(*,*) ' Parameter value (default)     :'// &
704              ' rho=', rho_q(iq)       
705          endif
706          !option is_rgcs
707          if (index(tracline,'is_rgcs') .ne. 0) then
708            read(tracline(index(tracline,'is_rgcs=') &
709              +len('is_rgcs='):),*) is_rgcs(iq)
710            write(*,*)'Parameter value (traceur.def) :'// &
711              'is_rgcs=',is_rgcs(iq)
712          else
713            write(*,*)'Parameter value (default) : '// &
714              'is_rgcs = ',is_rgcs(iq)
715          endif
716      end subroutine get_tracdat
717
718      end
719
Note: See TracBrowser for help on using the repository browser.