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

Last change on this file since 3781 was 3724, checked in by mmaurice, 8 months ago

Generic PCM

Add radioactive tracers in order to compare tracers mixing between 3D and 1D.

MM

File size: 25.7 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       IF (.NOT. allocated(constants_RCPV_generic)) allocate(constants_RCPV_generic(nq))
131       IF (.NOT.ALLOCATED(half_life))         ALLOCATE(half_life(nqtot))
132       IF (.NOT.ALLOCATED(top_prod))         ALLOCATE(top_prod(nqtot))
133       IF (.NOT.ALLOCATED(bot_prod))         ALLOCATE(bot_prod(nqtot))
134
135       !! initialization
136       alpha_lift(:)     = 0.
137       alpha_devil(:)    = 0.
138       mmol(:)           = 0.
139       aki(:)            = 0.
140       cpi(:)            = 0.
141       is_chim(:)        = 0
142       is_rad(:)         = 0
143       is_recomb(:)      = 0
144       is_recomb_qset(:) = 0
145       is_recomb_qotf(:) = 0
146       half_life(:)      = 0.
147       top_prod(:)       = 0.
148       bot_prod(:)       = 0.
149       
150       ! Added by JVO 2017 : these arrays are handled later
151       ! -> initialization is the least we can do, please !!!
152       radius(:)=0.
153       qext(:)=0.
154
155       ! For condensable tracers, by Lucas Teinturier and Noé Clément (2022)
156
157       is_condensable(:)= 0
158       is_rgcs(:) = 0
159       constants_mass(:)=0
160       constants_delta_vapH(:)=0
161       constants_Tref(:)=0
162       constants_Pref(:)=0
163       constants_epsi_generic(:)=0
164       constants_RLVTT_generic(:)=0
165       constants_metallicity_coeff(:)=0
166       constants_RCPV_generic(:)=0
167
168       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
169
170
171! Initialization: Read tracers names from traceur.def
172        do iq=1,nq
173          if (is_master) read(407,'(A)') tracline
174          call bcast(tracline)
175          blank = index(tracline,' ') ! Find position of 1st blank in tracline
176          noms(iq) = tracline(1:blank) !ensure that in Modern-trac case, noms is actually the name and not all properties
177        enddo
178
179! Identify tracers by their names: (and set corresponding values of mmol)
180      ! 0. initialize tracer indexes to zero:
181      ! NB: igcm_* indexes are commons in 'tracer.h'
182      do iq=1,nq
183        igcm_dustbin(iq)=0
184      enddo
185      igcm_dust_mass=0
186      igcm_dust_number=0
187      igcm_h2o_vap=0
188      igcm_h2o_ice=0
189      igcm_co2=0
190      igcm_co=0
191      igcm_o=0
192      igcm_o1d=0
193      igcm_o2=0
194      igcm_o3=0
195      igcm_h=0
196      igcm_h2=0
197      igcm_oh=0
198      igcm_ho2=0
199      igcm_h2o2=0
200      igcm_n2=0
201      igcm_n=0
202      igcm_n2d=0
203      igcm_no=0
204      igcm_no2=0
205      igcm_ar=0
206      igcm_ar_n2=0
207      igcm_co2_ice=0
208
209      igcm_ch4=0
210      igcm_ch3=0
211      igcm_ch=0
212      igcm_3ch2=0
213      igcm_1ch2=0
214      igcm_cho=0
215      igcm_ch2o=0
216      igcm_ch3o=0
217      igcm_c=0
218      igcm_c2=0
219      igcm_c2h=0
220      igcm_c2h2=0
221      igcm_c2h3=0
222      igcm_c2h4=0
223      igcm_c2h6=0
224      igcm_ch2co=0
225      igcm_ch3co=0
226      igcm_hcaer=0
227
228      ! 1. find dust tracers
229      count=0
230
231      ! 2. find chemistry and water tracers
232      do iq=1,nq
233        if (noms(iq).eq."co2") then
234          igcm_co2=iq
235          mmol(igcm_co2)=44.010
236          count=count+1
237!          write(*,*) 'co2: count=',count
238        endif
239        if (noms(iq).eq."co2_ice") then
240          igcm_co2_ice=iq
241          mmol(igcm_co2_ice)=44.010
242          count=count+1
243!          write(*,*) 'co2_ice: count=',count
244        endif
245        if (noms(iq).eq."h2o_vap") then
246          igcm_h2o_vap=iq
247          mmol(igcm_h2o_vap)=18.01528
248          count=count+1
249!          write(*,*) 'h2o_vap: count=',count
250        endif
251        if (noms(iq).eq."h2o_ice") then
252          igcm_h2o_ice=iq
253          mmol(igcm_h2o_ice)=18.01528
254          count=count+1
255!          write(*,*) 'h2o_ice: count=',count
256        endif
257        if (noms(iq).eq."co") then
258          igcm_co=iq
259          mmol(igcm_co)=28.010
260          count=count+1
261        endif
262        if (noms(iq).eq."o") then
263          igcm_o=iq
264          mmol(igcm_o)=15.9994
265          count=count+1
266        endif
267        if (noms(iq).eq."o1d") then
268          igcm_o1d=iq
269          mmol(igcm_o1d)=15.9994
270          count=count+1
271        endif
272        if (noms(iq).eq."o2") then
273          igcm_o2=iq
274          mmol(igcm_o2)=31.998
275          count=count+1
276        endif
277        if (noms(iq).eq."o3") then
278          igcm_o3=iq
279          mmol(igcm_o3)=47.9982
280          count=count+1
281        endif
282        if (noms(iq).eq."h") then
283          igcm_h=iq
284          mmol(igcm_h)=1.00784
285          count=count+1
286        endif
287        if (noms(iq).eq."h2") then
288          igcm_h2=iq
289          mmol(igcm_h2)=2.016
290          count=count+1
291        endif
292        if (noms(iq).eq."oh") then
293          igcm_oh=iq
294          mmol(igcm_oh)=17.008
295          count=count+1
296        endif
297        if (noms(iq).eq."ho2") then
298          igcm_ho2=iq
299          mmol(igcm_ho2)=33.0067
300          count=count+1
301        endif
302        if (noms(iq).eq."h2o2") then
303          igcm_h2o2=iq
304          mmol(igcm_h2o2)=34.0147
305          count=count+1
306        endif
307        if (noms(iq).eq."n2") then
308          igcm_n2=iq
309          mmol(igcm_n2)=28.0134
310          count=count+1
311        endif
312        if (noms(iq).eq."ch4") then
313          igcm_ch4=iq
314          mmol(igcm_ch4)=16.0425
315          count=count+1
316        endif
317        if (noms(iq).eq."ar") then
318          igcm_ar=iq
319          mmol(igcm_ar)=39.948
320          count=count+1
321        endif
322        if (noms(iq).eq."n") then
323          igcm_n=iq
324          mmol(igcm_n)=14.0067
325          count=count+1
326        endif
327        if (noms(iq).eq."no") then
328          igcm_no=iq
329          mmol(igcm_no)=30.0061
330          count=count+1
331        endif
332        if (noms(iq).eq."no2") then
333          igcm_no2=iq
334          mmol(igcm_no2)=46.0055
335          count=count+1
336        endif
337        if (noms(iq).eq."n2d") then
338          igcm_n2d=iq
339          mmol(igcm_n2d)=28.0134
340          count=count+1
341        endif
342        if (noms(iq).eq."ch3") then
343          igcm_ch3=iq
344          mmol(igcm_ch3)=15.034
345          count=count+1
346        endif
347        if (noms(iq).eq."ch") then
348          igcm_ch=iq
349          mmol(igcm_ch)=13.0186
350          count=count+1
351        endif
352        if (noms(iq).eq."3ch2") then
353          igcm_3ch2=iq
354          mmol(igcm_3ch2)=14.0266
355          count=count+1
356        endif
357        if (noms(iq).eq."1ch2") then
358          igcm_1ch2=iq
359          mmol(igcm_1ch2)=14.0266
360          count=count+1
361        endif
362        if (noms(iq).eq."cho") then
363          igcm_cho=iq
364          mmol(igcm_cho)=29.018
365          count=count+1
366        endif
367        if (noms(iq).eq."ch2o") then
368          igcm_ch2o=iq
369          mmol(igcm_ch2o)=30.031
370          count=count+1
371        endif
372        if (noms(iq).eq."ch3o") then
373          igcm_ch3o=iq
374          mmol(igcm_ch3o)=31.0339
375          count=count+1
376        endif
377        if (noms(iq).eq."c") then
378          igcm_c=iq
379          mmol(igcm_c)=12.0107
380          count=count+1
381        endif
382        if (noms(iq).eq."c2") then
383          igcm_c2=iq
384          mmol(igcm_c2)=24.0214
385          count=count+1
386        endif
387        if (noms(iq).eq."c2h") then
388          igcm_c2h=iq
389          mmol(igcm_c2h)=25.0293
390          count=count+1
391        endif
392        if (noms(iq).eq."c2h2") then
393          igcm_c2h2=iq
394          mmol(igcm_c2h2)=26.038
395          count=count+1
396        endif
397        if (noms(iq).eq."c2h3") then
398          igcm_c2h3=iq
399          mmol(igcm_c2h3)=27.0452
400          count=count+1
401        endif
402        if (noms(iq).eq."c2h4") then
403          igcm_c2h4=iq
404          mmol(igcm_c2h4)=28.053
405          count=count+1
406        endif
407        if (noms(iq).eq."c2h6") then
408          igcm_c2h6=iq
409          mmol(igcm_c2h6)=30.069
410          count=count+1
411        endif
412        if (noms(iq).eq."ch2co") then
413          igcm_ch2co=iq
414          mmol(igcm_ch2co)=42.0367
415          count=count+1
416        endif
417        if (noms(iq).eq."ch3co") then
418          igcm_ch3co=iq
419          mmol(igcm_ch3co)=43.0446
420          count=count+1
421        endif
422        if (noms(iq).eq."hcaer") then
423          igcm_hcaer=iq
424          mmol(igcm_hcaer)=50.
425          count=count+1
426        endif
427
428      enddo ! of do iq=1,nq
429
430      ! 3. find condensable traceurs different from h2o and co2
431      do iq=1,nq
432        if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
433          count=count+1
434        endif
435        if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
436          count=count+1
437        endif
438
439      enddo ! of do iq=1,nq
440     
441      ! check that we identified all tracers:
442      if (count.ne.nq) then
443        write(*,*) "initracer: found only ",count," tracers"
444        write(*,*) "               expected ",nq
445        do iq=1,count
446          write(*,*)'      ',iq,' ',trim(noms(iq))
447        enddo
448!        stop
449      else
450        write(*,*) "initracer: found all expected tracers, namely:"
451        do iq=1,nq
452          write(*,*)'      ',iq,' ',trim(noms(iq))
453        enddo
454      endif
455
456      ! Get data of tracers. Need to rewind traceur.def first
457      if (is_master) then
458       rewind(407)
459       do
460        read(407,'(A)') tracline
461        if (index(tracline,"#") .ne. 1) then
462          exit
463        endif
464       enddo
465      endif
466      do iq=1,nqtot
467        if (is_master) read(407,'(A)') tracline
468        call bcast(tracline)
469        call get_tracdat(iq, tracline)
470      enddo
471
472      if (is_master) close(407)
473
474      ! Get specific data of condensable tracers
475      do iq=1,nq
476        if((is_condensable(iq)==1)) then
477                write(*,*) "There is a specie which is condensable, for generic condensation : ", noms(iq)
478                write(*,*) 'looking specie parameters for : ',noms(iq)(1:len(trim(noms(iq)))-4)
479                call specie_parameters_table(noms(iq)(1:len(trim(noms(iq)))-4))
480                constants_mass(iq)=m
481                constants_delta_vapH(iq)=delta_vapH
482                constants_Tref(iq)=Tref
483                constants_Pref(iq)=Pref
484                constants_epsi_generic(iq)=epsi_generic
485                constants_RLVTT_generic(iq)=RLVTT_generic
486                constants_metallicity_coeff(iq)=metallicity_coeff
487                constants_RCPV_generic(iq)=RCPV_generic
488        else
489                write(*,*) "This tracer is not condensable, for generic condensation :  : ", noms(iq)
490                write(*,*) "We keep condensable constants at zero"
491        endif !(is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
492      enddo ! iq=1,nq
493
494      ! Calculate number of species in the chemistry
495      nesp = sum(is_chim)
496      write(*,*) 'Number of species in the chemistry nesp = ',nesp
497
498      ! Calculate number of generic tracers
499      ngt = sum(is_condensable)
500      write(*,*) 'Number of generic tracer is  ngt = ',ngt
501
502      ! Calculate number of radiative generic condensable species
503      n_rgcs = sum(is_rgcs)
504      write(*,*)'Number of Radiative Generic Condensable Species is n_rgcs = ',n_rgcs
505      if (n_rgcs> ngt/2) then
506        write(*,*) 'You have more Radiative Generic Condensable Species than Generic Condensable Species'
507        write(*,*)'This is not possible: check your Modern traceur.def'
508        call abort_physic("initracer, issue with # of RGCS and GCS")
509      endif
510
511!     Processing modern traceur options
512      if(moderntracdef) then
513        call ini_recombin
514      endif
515
516!------------------------------------------------------------
517!     Initialisation tracers ....
518!------------------------------------------------------------
519      ! rho_q(1:nq)=0
520
521      rho_dust=2500.  ! Mars dust density (kg.m-3)
522      rho_ice=920.    ! Water ice density (kg.m-3)
523      rho_co2=1620.   ! CO2 ice density (kg.m-3)
524
525!     Initialization for water vapor
526!     ------------------------------
527      if(water) then
528        radius(igcm_h2o_vap)=0.
529        Qext(igcm_h2o_vap)=0.
530        alpha_lift(igcm_h2o_vap) =0.
531        alpha_devil(igcm_h2o_vap)=0.
532        qextrhor(igcm_h2o_vap)= 0.
533
534         !! this is defined in surfdat_h.F90
535         IF (.not.ALLOCATED(dryness)) ALLOCATE(dryness(ngrid))
536
537         do ig=1,ngrid
538           dryness(ig) = 1.
539         enddo
540
541
542           radius(igcm_h2o_ice)=3.e-6
543           rho_q(igcm_h2o_ice)=rho_ice
544           Qext(igcm_h2o_ice)=0.
545!           alpha_lift(igcm_h2o_ice) =0.
546!           alpha_devil(igcm_h2o_ice)=0.
547           qextrhor(igcm_h2o_ice)= (3./4.)*Qext(igcm_h2o_ice) &
548            / (rho_ice*radius(igcm_h2o_ice))
549
550
551      end if  ! (water)
552
553
554!
555!     some extra (possibly redundant) sanity checks for tracers:
556!     ---------------------------------------------------------
557       if (water) then
558       ! verify that we indeed have h2o_vap and h2o_ice tracers
559         if (igcm_h2o_vap.eq.0) then
560           write(*,*) "initracer: error !!"
561           write(*,*) "  cannot use water option without ",&
562                     "an h2o_vap tracer !"
563           stop
564         endif
565         if (igcm_h2o_ice.eq.0) then
566           write(*,*) "initracer: error !!"
567           write(*,*) "  cannot use water option without ",&
568                     "an h2o_ice tracer !"
569           stop
570         endif
571       endif
572
573
574!     Output for records:
575!     ~~~~~~~~~~~~~~~~~~
576      write(*,*)
577      Write(*,*) '******** initracer : dust transport parameters :'
578      write(*,*) 'alpha_lift = ', alpha_lift
579      write(*,*) 'alpha_devil = ', alpha_devil
580      write(*,*) 'radius  = ', radius
581      write(*,*) 'Qext  = ', qext
582      write(*,*)
583
584      contains
585
586      subroutine get_tracdat(iq,tracline)
587        !-------------------ADDING NEW OPTIONS-------------------
588        ! Duplicate if sentence for adding new options
589        ! Do not forget to add the new variables in tracer_h.F90
590        ! Do not forget to allocate and initialize the new variables above
591        ! Please update list of options in "LMDZ.GENERIC/deftank/traceur.def.modern"
592        !-------------------------------------------------------
593          use tracer_h
594          implicit none
595          integer,           intent(in) :: iq       ! tracer index
596          character(len=500),intent(in) :: tracline ! traceur.def lines with parameters
597 
598          read(tracline,*) noms(iq)
599          ! JVO 20 : We should add a sanity check aborting when duplicates in names !
600          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
601          ! option mmol
602          if (index(tracline,'mmol='   ) /= 0) then
603              read(tracline(index(tracline,'mmol=')+len('mmol='):),*)&
604                  mmol(iq)
605              write(*,*) ' Parameter value (traceur.def) : mmol=', &
606                  mmol(iq)
607          else
608              write(*,*) ' Parameter value (default)     : mmol=', &
609                 mmol(iq)
610          end if
611          ! option aki
612          if (index(tracline,'aki='   ) /= 0) then
613              read(tracline(index(tracline,'aki=')+len('aki='):),*) &
614                  aki(iq)
615              write(*,*) ' Parameter value (traceur.def) : aki=', &
616                  aki(iq)
617          else
618              write(*,*) ' Parameter value (default)     : aki=', &
619                  aki(iq)
620          end if
621          ! option cpi
622          if (index(tracline,'cpi='   ) /= 0) then
623              read(tracline(index(tracline,'cpi=')+len('cpi='):),*) &
624                  cpi(iq)
625              write(*,*) ' Parameter value (traceur.def) : cpi=', &
626                  cpi(iq)
627          else
628              write(*,*) ' Parameter value (default)     : cpi=', &
629                  cpi(iq)
630          end if
631          ! option is_chim
632          if (index(tracline,'is_chim='   ) /= 0) then
633          read(tracline(index(tracline,'is_chim=')+len('is_chim='):),*) &
634                  is_chim(iq)
635              write(*,*) ' Parameter value (traceur.def) : is_chim=', &
636                  is_chim(iq)
637          else
638              write(*,*) ' Parameter value (default)     : is_chim=', &
639                  is_chim(iq)
640          end if
641          ! option is_rad
642          if (index(tracline,'is_rad=') /= 0) then
643              read(tracline(index(tracline,'is_rad=') &
644              +len('is_rad='):),*) is_rad(iq)
645              write(*,*) ' Parameter value (traceur.def) : is_rad=', &
646              is_rad(iq)
647          else
648              write(*,*) ' Parameter value (default)     : is_rad=',  &
649              is_rad(iq)
650          end if
651          ! option is_recomb
652          if (index(tracline,'is_recomb=') /= 0) then
653              read(tracline(index(tracline,'is_recomb=') &
654              +len('is_recomb='):),*) is_recomb(iq)
655              write(*,*) ' Parameter value (traceur.def) : is_recomb=', &
656              is_recomb(iq)
657          else
658              write(*,*) ' Parameter value (default)     : is_recomb=', &
659              is_recomb(iq)
660          end if
661          ! option is_recomb_qset
662          if (index(tracline,'is_recomb_qset=') /= 0) then
663              read(tracline(index(tracline,'is_recomb_qset=') &
664              +len('is_recomb_qset='):),*) is_recomb_qset(iq)
665              write(*,*) ' Parameter value (traceur.def) :'// &
666              ' is_recomb_qset=', &
667              is_recomb_qset(iq)
668          else
669              write(*,*) ' Parameter value (default)     :'// &
670              ' is_recomb_qset=', &
671              is_recomb_qset(iq)
672          endif
673          ! option is_recomb_qotf
674          if (index(tracline,'is_recomb_qotf=') /= 0) then
675              read(tracline(index(tracline,'is_recomb_qotf=') &
676              +len('is_recomb_qotf='):),*) is_recomb_qotf(iq)
677              write(*,*) ' Parameter value (traceur.def) :'// &
678              ' is_recomb_qotf=', &
679              is_recomb_qotf(iq)
680          else
681              write(*,*) ' Parameter value (default)     :'// &
682              ' is_recomb_qotf=',  &
683              is_recomb_qotf(iq)
684          end if
685          !option is_condensable (LT)
686          if (index(tracline,'is_condensable=') /=0) then
687            read(tracline(index(tracline,'is_condensable=') &
688              +len('is_condensable='):),*) is_condensable(iq)
689            write(*,*) ' Parameter value (traceur.def) :'// &
690              ' is_condensable=', is_condensable(iq)       
691          else
692              write(*,*) ' Parameter value (default)     :'// &
693              ' is_condensable=', is_condensable(iq)       
694          endif
695          !option radius
696          if (index(tracline,'radius=') .ne. 0) then
697            read(tracline(index(tracline,'radius=') &
698              +len('radius='):),*) radius(iq)
699            write(*,*)'Parameter value (traceur.def) :'// &
700              "radius=",radius(iq), "m"
701          else
702            write(*,*) ' Parameter value (default)     :'// &
703            ' radius=', radius(iq)," m"     
704          endif
705          !option rho
706          if (index(tracline,'rho=') .ne. 0) then
707            read(tracline(index(tracline,'rho=') &
708              +len('rho='):),*) rho_q(iq)
709            write(*,*)'Parameter value (traceur.def) :'// &
710              "rho=",rho_q(iq)
711          else
712            write(*,*) ' Parameter value (default)     :'// &
713              ' rho=', rho_q(iq)       
714          endif
715          !option is_rgcs
716          if (index(tracline,'is_rgcs') .ne. 0) then
717            read(tracline(index(tracline,'is_rgcs=') &
718              +len('is_rgcs='):),*) is_rgcs(iq)
719            write(*,*)'Parameter value (traceur.def) :'// &
720              'is_rgcs=',is_rgcs(iq)
721          else
722            write(*,*)'Parameter value (default) : '// &
723              'is_rgcs = ',is_rgcs(iq)
724          end if
725             ! option top_prod
726              if (index(tracline,'top_prod=') .ne. 0) then
727                  read(tracline(index(tracline,'top_prod=') &
728                  +len('top_prod='):),*) top_prod(iq)
729                  write(*,*) ' Parameter value (traceur.def) : top_prod=', &
730                  top_prod(iq)
731              else
732                  write(*,*) ' Parameter value (default)     : top_prod=',  &
733                  top_prod(iq)
734              end if
735              ! option bot_prod
736              if (index(tracline,'bot_prod=') .ne. 0) then
737                  read(tracline(index(tracline,'bot_prod=') &
738                  +len('bot_prod='):),*) bot_prod(iq)
739                  write(*,*) ' Parameter value (traceur.def) : bot_prod=', &
740                  bot_prod(iq)
741              else
742                  write(*,*) ' Parameter value (default)     : bot_prod=',  &
743                  bot_prod(iq)
744              end if
745      end subroutine get_tracdat
746
747      end
748
Note: See TracBrowser for help on using the repository browser.