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

Last change on this file since 3992 was 3893, checked in by gmilcareck, 5 months ago

Remove all "call abort" and "stop" statement in LMDZ.GENERIC and replacing them by call abort_physic().

File size: 25.3 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", "Unable to find traceur.def file",1)
57        end if
58!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
59        READ(407,'(A)') tracline
60        IF (trim(tracline).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
61          READ(tracline,*) nqtot ! Try standard traceur.def
62        ELSE
63         moderntracdef = .true.
64         DO
65           READ(407,'(A)',iostat=ierr) tracline
66           IF (ierr==0) THEN
67             IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
68               READ(tracline,*) nqtot
69               ! Temporary not implemented solution
70               if (nqtot/=nq) then
71                  call abort_physic('initracer','Different number of '// &
72       'tracers in dynamics and physics not managed yet',1)
73               endif
74               EXIT
75             ENDIF
76           ELSE ! If pb, or if reached EOF without having found nqtot
77              call abort_physic('initracer','Unable to read numbers '// &
78        'of tracers in traceur.def',1)
79           ENDIF
80         ENDDO
81        ENDIF ! if modern or standard traceur.def
82       
83       endif ! of if (is_master)
84       
85       ! share the information with other procs/threads (if any)
86       CALL bcast(nqtot)
87       CALL bcast(moderntracdef)
88       
89!! -----------------------------------------------------------------------
90       !! For the moment number of tracers in dynamics and physics are equal
91       nqtot=nq
92       !! we allocate once for all arrays in common in tracer_h.F90
93       !! (supposedly those are not used before call to initracer)
94       IF (.NOT.ALLOCATED(noms))           ALLOCATE(noms(nq))
95       IF (.NOT.ALLOCATED(mmol))           ALLOCATE(mmol(nq))
96       IF (.NOT.ALLOCATED(aki))            ALLOCATE(aki(nqtot))
97       IF (.NOT.ALLOCATED(cpi))            ALLOCATE(cpi(nqtot))
98       IF (.NOT.ALLOCATED(radius))         ALLOCATE(radius(nq))
99       IF (.NOT.ALLOCATED(rho_q))          ALLOCATE(rho_q(nq))
100       IF (.NOT.ALLOCATED(qext))           ALLOCATE(qext(nq))
101       IF (.NOT.ALLOCATED(alpha_lift))     ALLOCATE(alpha_lift(nq))
102       IF (.NOT.ALLOCATED(alpha_devil))    ALLOCATE(alpha_devil(nq))
103       IF (.NOT.ALLOCATED(qextrhor))       ALLOCATE(qextrhor(nq))
104       IF (.NOT.ALLOCATED(igcm_dustbin))   ALLOCATE(igcm_dustbin(nq))
105       IF (.NOT.ALLOCATED(is_chim))        ALLOCATE(is_chim(nqtot))
106       IF (.NOT.ALLOCATED(is_rad))         ALLOCATE(is_rad(nqtot))
107       IF (.NOT.ALLOCATED(is_recomb))      ALLOCATE(is_recomb(nqtot))
108       IF (.NOT.ALLOCATED(is_recomb_qset)) THEN
109         ALLOCATE(is_recomb_qset(nqtot))
110       ENDIF
111       IF (.NOT.ALLOCATED(is_recomb_qotf)) THEN
112         ALLOCATE(is_recomb_qotf(nqtot))
113       ENDIF
114       IF (.NOT. allocated(is_condensable)) allocate(is_condensable(nq)) !LT
115       IF (.NOT. allocated(is_rgcs)) allocate(is_rgcs(nq)) !LT
116       IF (.NOT. allocated(constants_mass)) allocate(constants_mass(nq))
117       IF (.NOT. allocated(constants_delta_vapH)) allocate(constants_delta_vapH(nq))
118       IF (.NOT. allocated(constants_Tref)) allocate(constants_Tref(nq))
119       IF (.NOT. allocated(constants_Pref)) allocate(constants_Pref(nq))
120       IF (.NOT. allocated(constants_epsi_generic)) allocate(constants_epsi_generic(nq))
121       IF (.NOT. allocated(constants_RLVTT_generic)) allocate(constants_RLVTT_generic(nq))
122       IF (.NOT. allocated(constants_metallicity_coeff)) allocate(constants_metallicity_coeff(nq))
123       IF (.NOT. allocated(constants_RCPV_generic)) allocate(constants_RCPV_generic(nq))
124       IF (.NOT.ALLOCATED(half_life))         ALLOCATE(half_life(nqtot))
125       IF (.NOT.ALLOCATED(top_prod))         ALLOCATE(top_prod(nqtot))
126       IF (.NOT.ALLOCATED(bot_prod))         ALLOCATE(bot_prod(nqtot))
127
128       !! initialization
129       alpha_lift(:)     = 0.
130       alpha_devil(:)    = 0.
131       mmol(:)           = 0.
132       aki(:)            = 0.
133       cpi(:)            = 0.
134       is_chim(:)        = 0
135       is_rad(:)         = 0
136       is_recomb(:)      = 0
137       is_recomb_qset(:) = 0
138       is_recomb_qotf(:) = 0
139       half_life(:)      = 0.
140       top_prod(:)       = 0.
141       bot_prod(:)       = 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       constants_RCPV_generic(:)=0
160
161       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
162
163
164! Initialization: Read tracers names from traceur.def
165        do iq=1,nq
166          if (is_master) read(407,'(A)') tracline
167          call bcast(tracline)
168          blank = index(tracline,' ') ! Find position of 1st blank in tracline
169          noms(iq) = tracline(1:blank) !ensure that in Modern-trac case, noms is actually the name and not all properties
170        enddo
171
172! Identify tracers by their names: (and set corresponding values of mmol)
173      ! 0. initialize tracer indexes to zero:
174      ! NB: igcm_* indexes are commons in 'tracer.h'
175      do iq=1,nq
176        igcm_dustbin(iq)=0
177      enddo
178      igcm_dust_mass=0
179      igcm_dust_number=0
180      igcm_h2o_vap=0
181      igcm_h2o_ice=0
182      igcm_co2=0
183      igcm_co=0
184      igcm_o=0
185      igcm_o1d=0
186      igcm_o2=0
187      igcm_o3=0
188      igcm_h=0
189      igcm_h2=0
190      igcm_oh=0
191      igcm_ho2=0
192      igcm_h2o2=0
193      igcm_n2=0
194      igcm_n=0
195      igcm_n2d=0
196      igcm_no=0
197      igcm_no2=0
198      igcm_ar=0
199      igcm_ar_n2=0
200      igcm_co2_ice=0
201
202      igcm_ch4=0
203      igcm_ch3=0
204      igcm_ch=0
205      igcm_3ch2=0
206      igcm_1ch2=0
207      igcm_cho=0
208      igcm_ch2o=0
209      igcm_ch3o=0
210      igcm_c=0
211      igcm_c2=0
212      igcm_c2h=0
213      igcm_c2h2=0
214      igcm_c2h3=0
215      igcm_c2h4=0
216      igcm_c2h6=0
217      igcm_ch2co=0
218      igcm_ch3co=0
219      igcm_hcaer=0
220
221      ! 1. find dust tracers
222      count=0
223
224      ! 2. find chemistry and water tracers
225      do iq=1,nq
226        if (noms(iq).eq."co2") then
227          igcm_co2=iq
228          mmol(igcm_co2)=44.010
229          count=count+1
230!          write(*,*) 'co2: count=',count
231        endif
232        if (noms(iq).eq."co2_ice") then
233          igcm_co2_ice=iq
234          mmol(igcm_co2_ice)=44.010
235          count=count+1
236!          write(*,*) 'co2_ice: count=',count
237        endif
238        if (noms(iq).eq."h2o_vap") then
239          igcm_h2o_vap=iq
240          mmol(igcm_h2o_vap)=18.01528
241          count=count+1
242!          write(*,*) 'h2o_vap: count=',count
243        endif
244        if (noms(iq).eq."h2o_ice") then
245          igcm_h2o_ice=iq
246          mmol(igcm_h2o_ice)=18.01528
247          count=count+1
248!          write(*,*) 'h2o_ice: count=',count
249        endif
250        if (noms(iq).eq."co") then
251          igcm_co=iq
252          mmol(igcm_co)=28.010
253          count=count+1
254        endif
255        if (noms(iq).eq."o") then
256          igcm_o=iq
257          mmol(igcm_o)=15.9994
258          count=count+1
259        endif
260        if (noms(iq).eq."o1d") then
261          igcm_o1d=iq
262          mmol(igcm_o1d)=15.9994
263          count=count+1
264        endif
265        if (noms(iq).eq."o2") then
266          igcm_o2=iq
267          mmol(igcm_o2)=31.998
268          count=count+1
269        endif
270        if (noms(iq).eq."o3") then
271          igcm_o3=iq
272          mmol(igcm_o3)=47.9982
273          count=count+1
274        endif
275        if (noms(iq).eq."h") then
276          igcm_h=iq
277          mmol(igcm_h)=1.00784
278          count=count+1
279        endif
280        if (noms(iq).eq."h2") then
281          igcm_h2=iq
282          mmol(igcm_h2)=2.016
283          count=count+1
284        endif
285        if (noms(iq).eq."oh") then
286          igcm_oh=iq
287          mmol(igcm_oh)=17.008
288          count=count+1
289        endif
290        if (noms(iq).eq."ho2") then
291          igcm_ho2=iq
292          mmol(igcm_ho2)=33.0067
293          count=count+1
294        endif
295        if (noms(iq).eq."h2o2") then
296          igcm_h2o2=iq
297          mmol(igcm_h2o2)=34.0147
298          count=count+1
299        endif
300        if (noms(iq).eq."n2") then
301          igcm_n2=iq
302          mmol(igcm_n2)=28.0134
303          count=count+1
304        endif
305        if (noms(iq).eq."ch4") then
306          igcm_ch4=iq
307          mmol(igcm_ch4)=16.0425
308          count=count+1
309        endif
310        if (noms(iq).eq."ar") then
311          igcm_ar=iq
312          mmol(igcm_ar)=39.948
313          count=count+1
314        endif
315        if (noms(iq).eq."n") then
316          igcm_n=iq
317          mmol(igcm_n)=14.0067
318          count=count+1
319        endif
320        if (noms(iq).eq."no") then
321          igcm_no=iq
322          mmol(igcm_no)=30.0061
323          count=count+1
324        endif
325        if (noms(iq).eq."no2") then
326          igcm_no2=iq
327          mmol(igcm_no2)=46.0055
328          count=count+1
329        endif
330        if (noms(iq).eq."n2d") then
331          igcm_n2d=iq
332          mmol(igcm_n2d)=28.0134
333          count=count+1
334        endif
335        if (noms(iq).eq."ch3") then
336          igcm_ch3=iq
337          mmol(igcm_ch3)=15.034
338          count=count+1
339        endif
340        if (noms(iq).eq."ch") then
341          igcm_ch=iq
342          mmol(igcm_ch)=13.0186
343          count=count+1
344        endif
345        if (noms(iq).eq."3ch2") then
346          igcm_3ch2=iq
347          mmol(igcm_3ch2)=14.0266
348          count=count+1
349        endif
350        if (noms(iq).eq."1ch2") then
351          igcm_1ch2=iq
352          mmol(igcm_1ch2)=14.0266
353          count=count+1
354        endif
355        if (noms(iq).eq."cho") then
356          igcm_cho=iq
357          mmol(igcm_cho)=29.018
358          count=count+1
359        endif
360        if (noms(iq).eq."ch2o") then
361          igcm_ch2o=iq
362          mmol(igcm_ch2o)=30.031
363          count=count+1
364        endif
365        if (noms(iq).eq."ch3o") then
366          igcm_ch3o=iq
367          mmol(igcm_ch3o)=31.0339
368          count=count+1
369        endif
370        if (noms(iq).eq."c") then
371          igcm_c=iq
372          mmol(igcm_c)=12.0107
373          count=count+1
374        endif
375        if (noms(iq).eq."c2") then
376          igcm_c2=iq
377          mmol(igcm_c2)=24.0214
378          count=count+1
379        endif
380        if (noms(iq).eq."c2h") then
381          igcm_c2h=iq
382          mmol(igcm_c2h)=25.0293
383          count=count+1
384        endif
385        if (noms(iq).eq."c2h2") then
386          igcm_c2h2=iq
387          mmol(igcm_c2h2)=26.038
388          count=count+1
389        endif
390        if (noms(iq).eq."c2h3") then
391          igcm_c2h3=iq
392          mmol(igcm_c2h3)=27.0452
393          count=count+1
394        endif
395        if (noms(iq).eq."c2h4") then
396          igcm_c2h4=iq
397          mmol(igcm_c2h4)=28.053
398          count=count+1
399        endif
400        if (noms(iq).eq."c2h6") then
401          igcm_c2h6=iq
402          mmol(igcm_c2h6)=30.069
403          count=count+1
404        endif
405        if (noms(iq).eq."ch2co") then
406          igcm_ch2co=iq
407          mmol(igcm_ch2co)=42.0367
408          count=count+1
409        endif
410        if (noms(iq).eq."ch3co") then
411          igcm_ch3co=iq
412          mmol(igcm_ch3co)=43.0446
413          count=count+1
414        endif
415        if (noms(iq).eq."hcaer") then
416          igcm_hcaer=iq
417          mmol(igcm_hcaer)=50.
418          count=count+1
419        endif
420
421      enddo ! of do iq=1,nq
422
423      ! 3. find condensable traceurs different from h2o and co2
424      do iq=1,nq
425        if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
426          count=count+1
427        endif
428        if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
429          count=count+1
430        endif
431
432      enddo ! of do iq=1,nq
433     
434      ! check that we identified all tracers:
435      if (count.ne.nq) then
436        write(*,*) "initracer: found only ",count," tracers"
437        write(*,*) "               expected ",nq
438        do iq=1,count
439          write(*,*)'      ',iq,' ',trim(noms(iq))
440        enddo
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                constants_RCPV_generic(iq)=RCPV_generic
480        else
481                write(*,*) "This tracer is not condensable, for generic condensation :  : ", noms(iq)
482                write(*,*) "We keep condensable constants at zero"
483        endif !(is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
484      enddo ! iq=1,nq
485
486      ! Calculate number of species in the chemistry
487      nesp = sum(is_chim)
488      write(*,*) 'Number of species in the chemistry nesp = ',nesp
489
490      ! Calculate number of generic tracers
491      ngt = sum(is_condensable)
492      write(*,*) 'Number of generic tracer is  ngt = ',ngt
493
494      ! Calculate number of radiative generic condensable species
495      n_rgcs = sum(is_rgcs)
496      write(*,*)'Number of Radiative Generic Condensable Species is n_rgcs = ',n_rgcs
497      if (n_rgcs> ngt/2) then
498        write(*,*) 'You have more Radiative Generic Condensable Species than Generic Condensable Species'
499        write(*,*)'This is not possible: check your Modern traceur.def'
500        call abort_physic("initracer", "issue with # of RGCS and GCS",1)
501      endif
502
503!     Processing modern traceur options
504      if(moderntracdef) then
505        call ini_recombin
506      endif
507
508!------------------------------------------------------------
509!     Initialisation tracers ....
510!------------------------------------------------------------
511      ! rho_q(1:nq)=0
512
513      rho_dust=2500.  ! Mars dust density (kg.m-3)
514      rho_ice=920.    ! Water ice density (kg.m-3)
515      rho_co2=1620.   ! CO2 ice density (kg.m-3)
516
517!     Initialization for water vapor
518!     ------------------------------
519      if(water) then
520        radius(igcm_h2o_vap)=0.
521        Qext(igcm_h2o_vap)=0.
522        alpha_lift(igcm_h2o_vap) =0.
523        alpha_devil(igcm_h2o_vap)=0.
524        qextrhor(igcm_h2o_vap)= 0.
525
526         !! this is defined in surfdat_h.F90
527         IF (.not.ALLOCATED(dryness)) ALLOCATE(dryness(ngrid))
528
529         do ig=1,ngrid
530           dryness(ig) = 1.
531         enddo
532
533
534           radius(igcm_h2o_ice)=3.e-6
535           rho_q(igcm_h2o_ice)=rho_ice
536           Qext(igcm_h2o_ice)=0.
537!           alpha_lift(igcm_h2o_ice) =0.
538!           alpha_devil(igcm_h2o_ice)=0.
539           qextrhor(igcm_h2o_ice)= (3./4.)*Qext(igcm_h2o_ice) &
540            / (rho_ice*radius(igcm_h2o_ice))
541
542
543      end if  ! (water)
544
545
546!
547!     some extra (possibly redundant) sanity checks for tracers:
548!     ---------------------------------------------------------
549       if (water) then
550       ! verify that we indeed have h2o_vap and h2o_ice tracers
551         if (igcm_h2o_vap.eq.0) then
552           call abort_physic("initracer", "cannot use water option without an h2o_vap tracer", 1)
553         endif
554         if (igcm_h2o_ice.eq.0) then
555           call abort_physic("initracer", "cannot use water option without an h2o_ice tracer", 1)
556         endif
557       endif
558
559
560!     Output for records:
561!     ~~~~~~~~~~~~~~~~~~
562      write(*,*)
563      Write(*,*) '******** initracer : dust transport parameters :'
564      write(*,*) 'alpha_lift = ', alpha_lift
565      write(*,*) 'alpha_devil = ', alpha_devil
566      write(*,*) 'radius  = ', radius
567      write(*,*) 'Qext  = ', qext
568      write(*,*)
569
570      contains
571
572      subroutine get_tracdat(iq,tracline)
573        !-------------------ADDING NEW OPTIONS-------------------
574        ! Duplicate if sentence for adding new options
575        ! Do not forget to add the new variables in tracer_h.F90
576        ! Do not forget to allocate and initialize the new variables above
577        ! Please update list of options in "LMDZ.GENERIC/deftank/traceur.def.modern"
578        !-------------------------------------------------------
579          use tracer_h
580          implicit none
581          integer,           intent(in) :: iq       ! tracer index
582          character(len=500),intent(in) :: tracline ! traceur.def lines with parameters
583 
584          read(tracline,*) noms(iq)
585          ! JVO 20 : We should add a sanity check aborting when duplicates in names !
586          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
587          ! option mmol
588          if (index(tracline,'mmol='   ) /= 0) then
589              read(tracline(index(tracline,'mmol=')+len('mmol='):),*)&
590                  mmol(iq)
591              write(*,*) ' Parameter value (traceur.def) : mmol=', &
592                  mmol(iq)
593          else
594              write(*,*) ' Parameter value (default)     : mmol=', &
595                 mmol(iq)
596          end if
597          ! option aki
598          if (index(tracline,'aki='   ) /= 0) then
599              read(tracline(index(tracline,'aki=')+len('aki='):),*) &
600                  aki(iq)
601              write(*,*) ' Parameter value (traceur.def) : aki=', &
602                  aki(iq)
603          else
604              write(*,*) ' Parameter value (default)     : aki=', &
605                  aki(iq)
606          end if
607          ! option cpi
608          if (index(tracline,'cpi='   ) /= 0) then
609              read(tracline(index(tracline,'cpi=')+len('cpi='):),*) &
610                  cpi(iq)
611              write(*,*) ' Parameter value (traceur.def) : cpi=', &
612                  cpi(iq)
613          else
614              write(*,*) ' Parameter value (default)     : cpi=', &
615                  cpi(iq)
616          end if
617          ! option is_chim
618          if (index(tracline,'is_chim='   ) /= 0) then
619          read(tracline(index(tracline,'is_chim=')+len('is_chim='):),*) &
620                  is_chim(iq)
621              write(*,*) ' Parameter value (traceur.def) : is_chim=', &
622                  is_chim(iq)
623          else
624              write(*,*) ' Parameter value (default)     : is_chim=', &
625                  is_chim(iq)
626          end if
627          ! option is_rad
628          if (index(tracline,'is_rad=') /= 0) then
629              read(tracline(index(tracline,'is_rad=') &
630              +len('is_rad='):),*) is_rad(iq)
631              write(*,*) ' Parameter value (traceur.def) : is_rad=', &
632              is_rad(iq)
633          else
634              write(*,*) ' Parameter value (default)     : is_rad=',  &
635              is_rad(iq)
636          end if
637          ! option is_recomb
638          if (index(tracline,'is_recomb=') /= 0) then
639              read(tracline(index(tracline,'is_recomb=') &
640              +len('is_recomb='):),*) is_recomb(iq)
641              write(*,*) ' Parameter value (traceur.def) : is_recomb=', &
642              is_recomb(iq)
643          else
644              write(*,*) ' Parameter value (default)     : is_recomb=', &
645              is_recomb(iq)
646          end if
647          ! option is_recomb_qset
648          if (index(tracline,'is_recomb_qset=') /= 0) then
649              read(tracline(index(tracline,'is_recomb_qset=') &
650              +len('is_recomb_qset='):),*) is_recomb_qset(iq)
651              write(*,*) ' Parameter value (traceur.def) :'// &
652              ' is_recomb_qset=', &
653              is_recomb_qset(iq)
654          else
655              write(*,*) ' Parameter value (default)     :'// &
656              ' is_recomb_qset=', &
657              is_recomb_qset(iq)
658          endif
659          ! option is_recomb_qotf
660          if (index(tracline,'is_recomb_qotf=') /= 0) then
661              read(tracline(index(tracline,'is_recomb_qotf=') &
662              +len('is_recomb_qotf='):),*) is_recomb_qotf(iq)
663              write(*,*) ' Parameter value (traceur.def) :'// &
664              ' is_recomb_qotf=', &
665              is_recomb_qotf(iq)
666          else
667              write(*,*) ' Parameter value (default)     :'// &
668              ' is_recomb_qotf=',  &
669              is_recomb_qotf(iq)
670          end if
671          !option is_condensable (LT)
672          if (index(tracline,'is_condensable=') /=0) then
673            read(tracline(index(tracline,'is_condensable=') &
674              +len('is_condensable='):),*) is_condensable(iq)
675            write(*,*) ' Parameter value (traceur.def) :'// &
676              ' is_condensable=', is_condensable(iq)       
677          else
678              write(*,*) ' Parameter value (default)     :'// &
679              ' is_condensable=', is_condensable(iq)       
680          endif
681          !option radius
682          if (index(tracline,'radius=') .ne. 0) then
683            read(tracline(index(tracline,'radius=') &
684              +len('radius='):),*) radius(iq)
685            write(*,*)'Parameter value (traceur.def) :'// &
686              "radius=",radius(iq), "m"
687          else
688            write(*,*) ' Parameter value (default)     :'// &
689            ' radius=', radius(iq)," m"     
690          endif
691          !option rho
692          if (index(tracline,'rho=') .ne. 0) then
693            read(tracline(index(tracline,'rho=') &
694              +len('rho='):),*) rho_q(iq)
695            write(*,*)'Parameter value (traceur.def) :'// &
696              "rho=",rho_q(iq)
697          else
698            write(*,*) ' Parameter value (default)     :'// &
699              ' rho=', rho_q(iq)       
700          endif
701          !option is_rgcs
702          if (index(tracline,'is_rgcs') .ne. 0) then
703            read(tracline(index(tracline,'is_rgcs=') &
704              +len('is_rgcs='):),*) is_rgcs(iq)
705            write(*,*)'Parameter value (traceur.def) :'// &
706              'is_rgcs=',is_rgcs(iq)
707          else
708            write(*,*)'Parameter value (default) : '// &
709              'is_rgcs = ',is_rgcs(iq)
710          end if
711             ! option top_prod
712              if (index(tracline,'top_prod=') .ne. 0) then
713                  read(tracline(index(tracline,'top_prod=') &
714                  +len('top_prod='):),*) top_prod(iq)
715                  write(*,*) ' Parameter value (traceur.def) : top_prod=', &
716                  top_prod(iq)
717              else
718                  write(*,*) ' Parameter value (default)     : top_prod=',  &
719                  top_prod(iq)
720              end if
721              ! option bot_prod
722              if (index(tracline,'bot_prod=') .ne. 0) then
723                  read(tracline(index(tracline,'bot_prod=') &
724                  +len('bot_prod='):),*) bot_prod(iq)
725                  write(*,*) ' Parameter value (traceur.def) : bot_prod=', &
726                  bot_prod(iq)
727              else
728                  write(*,*) ' Parameter value (default)     : bot_prod=',  &
729                  bot_prod(iq)
730              end if
731      end subroutine get_tracdat
732
733      end
734
Note: See TracBrowser for help on using the repository browser.