source: trunk/LMDZ.PLUTO/libf/phypluto/initracer.F90 @ 3572

Last change on this file since 3572 was 3572, checked in by debatzbr, 3 weeks ago

Remove generic_aerosols and generic_condensation, along with their related variables (useless). RENAME THE VARIABLE AEROHAZE TO OPTICHAZE.

File size: 23.3 KB
Line 
1      SUBROUTINE initracer(ngrid,nq)
2
3      use surfdat_h, ONLY: dryness
4      USE tracer_h
5      USE callkeys_mod, only: optichaze,nb_monomer,haze,fractal,fasthaze,rad_haze,callmufi
6      USE recombin_corrk_mod, ONLY: ini_recombin
7      USE mod_phys_lmdz_para, only: is_master, bcast
8      use aerosol_mod, only: iaero_haze,i_haze
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!            B. de Batz de Trenquelléon (2024): specific microphysical tracers   
27!=======================================================================
28
29      integer,intent(in) :: ngrid,nq
30
31      character(len=500) :: tracline   ! to read traceur.def lines
32      integer :: blank      !to store the index of 1st blank when reading tracers names
33      integer iq,ig,count,ierr
34      real r0_lift , reff_lift, rho_haze
35      integer nqhaze(nq)               ! to store haze tracers
36      integer i, ia, block, j
37      character(len=20) :: txt ! to store some text
38      character(LEN=20) :: tracername ! to temporarily store text
39      character(LEN=20) :: str
40
41!-----------------------------------------------------------------------
42!  radius(nq)      ! aerosol particle radius (m)
43!  rho_q(nq)       ! tracer densities (kg.m-3)
44!  qext(nq)        ! Single Scat. Extinction coeff at 0.67 um
45!  alpha_lift(nq)  ! saltation vertical flux/horiz flux ratio (m-1)
46!  alpha_devil(nq) ! lifting coeeficient by dust devil
47!  rho_dust          ! Mars dust density
48!  rho_ice           ! Water ice density
49!  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
50!  varian            ! Characteristic variance of log-normal distribution
51!-----------------------------------------------------------------------
52
53      if (is_master) then ! only the master proc/thread needs do this
54
55        moderntracdef=.false. ! For modern traceur.def (default false, old type)
56
57        open(407, form = 'formatted', status = 'old', &
58             file = 'traceur.def', iostat=ierr)
59        if (ierr /=0) then
60          ! call abort_physic('initracer',&
61          ! 'Problem in opening traceur.def',1)
62          print*,'Problem in opening traceur.def'
63          stop
64        end if
65!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
66        READ(407,'(A)') tracline
67        IF (trim(tracline).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
68          READ(tracline,*) nqtot ! Try standard traceur.def
69        ELSE
70         moderntracdef = .true.
71         DO
72           READ(407,'(A)',iostat=ierr) tracline
73           IF (ierr==0) THEN
74             IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
75               READ(tracline,*) nqtot
76               ! Temporary not implemented solution
77               if (nqtot/=nq) then
78      !            call abort_physic('initracer','Different number of '// &
79      ! 'tracers in dynamics and physics not managed yet',1)
80                print*,'!= nbr oftracers in dynamics and physics not managed yet'
81                stop
82               endif
83               EXIT
84             ENDIF
85           ELSE ! If pb, or if reached EOF without having found nqtot
86      !        call abort_physic('initracer','Unable to read numbers '// &
87      !  'of tracers in traceur.def',1)
88            print*,"unable to read numbers of tracer in tracer.def"
89            stop
90           ENDIF
91         ENDDO
92        ENDIF ! if modern or standard traceur.def
93
94       endif ! of if (is_master)
95
96       ! share the information with other procs/threads (if any)
97       CALL bcast(nqtot)
98       CALL bcast(moderntracdef)
99
100!! -----------------------------------------------------------------------
101       !! For the moment number of tracers in dynamics and physics are equal
102       nqtot=nq
103       !! we allocate once for all arrays in common in tracer_h.F90
104       !! (supposedly those are not used before call to initracer)
105       IF (.NOT.ALLOCATED(noms))           ALLOCATE(noms(nq))
106       IF (.NOT.ALLOCATED(mmol))           ALLOCATE(mmol(nq))
107       IF (.NOT.ALLOCATED(aki))            ALLOCATE(aki(nqtot))
108       IF (.NOT.ALLOCATED(cpi))            ALLOCATE(cpi(nqtot))
109       IF (.NOT.ALLOCATED(radius))         ALLOCATE(radius(nq))
110       IF (.NOT.ALLOCATED(rho_q))          ALLOCATE(rho_q(nq))
111       IF (.NOT.ALLOCATED(qext))           ALLOCATE(qext(nq))
112       IF (.NOT.ALLOCATED(alpha_lift))     ALLOCATE(alpha_lift(nq))
113       IF (.NOT.ALLOCATED(alpha_devil))    ALLOCATE(alpha_devil(nq))
114       IF (.NOT.ALLOCATED(qextrhor))       ALLOCATE(qextrhor(nq))
115      !  IF (.NOT.ALLOCATED(igcm_dustbin))   ALLOCATE(igcm_dustbin(nq))
116       IF (.NOT.ALLOCATED(is_chim))        ALLOCATE(is_chim(nqtot))
117       IF (.NOT.ALLOCATED(is_rad))         ALLOCATE(is_rad(nqtot))
118       IF (.NOT.ALLOCATED(is_recomb))      ALLOCATE(is_recomb(nqtot))
119       IF (.NOT.ALLOCATED(is_recomb_qset)) THEN
120         ALLOCATE(is_recomb_qset(nqtot))
121       ENDIF
122       IF (.NOT.ALLOCATED(is_recomb_qotf)) THEN
123         ALLOCATE(is_recomb_qotf(nqtot))
124       ENDIF
125       IF (.NOT. allocated(is_condensable)) allocate(is_condensable(nq)) !LT
126       IF (.NOT. allocated(is_rgcs)) allocate(is_rgcs(nq)) !LT
127       IF (.NOT. allocated(constants_mass)) allocate(constants_mass(nq))
128       IF (.NOT. allocated(constants_delta_gasH)) allocate(constants_delta_gasH(nq))
129       IF (.NOT. allocated(constants_Tref)) allocate(constants_Tref(nq))
130       IF (.NOT. allocated(constants_Pref)) allocate(constants_Pref(nq))
131       IF (.NOT. allocated(constants_RLVTT_generic)) allocate(constants_RLVTT_generic(nq))
132       IF (.NOT. allocated(constants_RCPV_generic)) allocate(constants_RCPV_generic(nq))
133
134       !! initialization
135       alpha_lift(:)     = 0.
136       alpha_devil(:)    = 0.
137       mmol(:)           = 0.
138       aki(:)            = 0.
139       cpi(:)            = 0.
140       is_chim(:)        = 0
141       is_rad(:)         = 0
142       is_recomb(:)      = 0
143       is_recomb_qset(:) = 0
144       is_recomb_qotf(:) = 0
145
146       ! Added by JVO 2017 : these arrays are handled later
147       ! -> initialization is the least we can do, please !!!
148       radius(:)=0.
149       qext(:)=0.
150
151       ! For condensable tracers, by Lucas Teinturier and Noé Clément (2022)
152
153       is_condensable(:)= 0
154       is_rgcs(:) = 0
155       constants_mass(:)=0
156       constants_delta_gasH(:)=0
157       constants_Tref(:)=0
158       constants_Pref(:)=0
159       constants_RLVTT_generic(:)=0
160       constants_RCPV_generic(:)=0
161
162       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
163
164
165! Initialization: Read tracers names from traceur.def
166        do iq=1,nq
167          if (is_master) read(407,'(A)') tracline
168          call bcast(tracline)
169          blank = index(tracline,' ') ! Find position of 1st blank in tracline
170          noms(iq) = tracline(1:blank) !ensure that in Modern-trac case, noms is actually the name and not all properties
171        enddo
172
173! Identify tracers by their names: (and set corresponding values of mmol)
174      ! 0. initialize tracer indexes to zero:
175      ! NB: igcm_* indexes are commons in 'tracer.h'
176      igcm_n2=0
177      igcm_ch4_gas=0
178      igcm_ch4_ice=0
179      igcm_prec_haze=0
180      igcm_co_gas=0
181      igcm_co_ice=0
182
183      nqhaze(:)=0
184      i=0
185      DO iq=1,nq
186         txt=noms(iq)
187         IF (txt(1:4).eq."haze") THEN
188            i=i+1
189            nqhaze(i)=iq
190         ENDIF
191      ENDDO
192      if ((haze.or.fasthaze).and.i==0) then
193         print*, 'Haze active but no haze tracer in traceur.def'
194         stop
195      endif
196
197      igcm_haze=0
198      igcm_haze10=0
199      igcm_haze30=0
200      igcm_haze50=0
201      igcm_haze100=0
202
203!     Eddy diffusion tracers
204      igcm_eddy1e6=0
205      igcm_eddy1e7=0
206      igcm_eddy5e7=0
207      igcm_eddy1e8=0
208      igcm_eddy5e8=0
209      write(*,*) 'initracer: noms() ', noms
210      rho_n2=1000        ! n2 ice
211      rho_ch4_ice=520.       ! rho ch4 ice
212      rho_co_ice=520.       ! rho ch4 ice
213      rho_haze=800.     ! haze
214
215      ! 1. find dust tracers
216      count=0
217
218      ! 2. find chemistry and water tracers
219      do iq=1,nq
220        if (noms(iq).eq."n2") then
221          igcm_n2=iq
222          mmol(igcm_n2)=28.
223          count=count+1
224          write(*,*) 'Tracer ',count,' = n2'
225        endif
226        if (noms(iq).eq."ch4_gas") then
227          igcm_ch4_gas=iq
228          mmol(igcm_ch4_gas)=16.
229          count=count+1
230          write(*,*) 'Tracer ',count,' = ch4 gas'
231        endif
232        if (noms(iq).eq."ch4_ice") then
233          igcm_ch4_ice=iq
234          mmol(igcm_ch4_ice)=16.
235          radius(iq)=3.e-6
236          rho_q(iq)=rho_ch4_ice
237          count=count+1
238          write(*,*) 'Tracer ',count,' = ch4 ice'
239        endif
240        if (noms(iq).eq."co_gas") then
241          igcm_co_gas=iq
242          mmol(igcm_co_gas)=28.
243          count=count+1
244          write(*,*) 'Tracer ',count,' = co gas'
245        endif
246        if (noms(iq).eq."co_ice") then
247          igcm_co_ice=iq
248          mmol(igcm_co_ice)=28.
249          radius(iq)=3.e-6
250          rho_q(iq)=rho_co_ice
251          count=count+1
252          write(*,*) 'Tracer ',count,' = co ice'
253        endif
254        if (noms(iq).eq."prec_haze") then
255          igcm_prec_haze=iq
256          count=count+1
257          write(*,*) 'Tracer ',count,' = prec haze'
258        endif
259        if (noms(iq).eq."haze") then
260          igcm_haze=iq
261          count=count+1
262          radius(iq)=rad_haze
263          rho_q(iq)=rho_haze
264          write(*,*) 'Tracer ',count,' = haze'
265        endif
266        if (noms(iq).eq."haze10") then
267          igcm_haze10=iq
268          count=count+1
269          radius(iq)=10.e-9
270          rho_q(iq)=rho_haze
271          write(*,*) 'Tracer ',count,' = haze10'
272        endif
273        if (noms(iq).eq."haze30") then
274          igcm_haze30=iq
275          count=count+1
276          radius(iq)=30.e-9
277          rho_q(iq)=rho_haze
278          write(*,*) 'Tracer ',count,' = haze30'
279        endif
280        if (noms(iq).eq."haze50") then
281          igcm_haze50=iq
282          count=count+1
283          radius(iq)=50.e-9
284          rho_q(iq)=rho_haze
285          write(*,*) 'Tracer ',count,' = haze50'
286        endif
287        if (noms(iq).eq."haze100") then
288          igcm_haze100=iq
289          count=count+1
290          radius(iq)=100.e-9
291          rho_q(iq)=rho_haze
292          write(*,*) 'Tracer ',count,' = haze100'
293        endif
294!       Eddy diffusion tracers
295        if (noms(iq).eq."eddy1e6") then
296          igcm_eddy1e6=iq
297          count=count+1
298          write(*,*) 'Tracer ',count,' = eddy1e6'
299        endif
300        if (noms(iq).eq."eddy1e7") then
301          igcm_eddy1e7=iq
302          count=count+1
303          write(*,*) 'Tracer ',count,' = eddy1e7'
304        endif
305        if (noms(iq).eq."eddy5e7") then
306          igcm_eddy5e7=iq
307          count=count+1
308          write(*,*) 'Tracer ',count,' = eddy5e7'
309        endif
310        if (noms(iq).eq."eddy1e8") then
311          igcm_eddy1e8=iq
312          count=count+1
313          write(*,*) 'Tracer ',count,' = eddy1e8'
314        endif
315        if (noms(iq).eq."eddy5e8") then
316          igcm_eddy5e8=iq
317          count=count+1
318          write(*,*) 'Tracer ',count,' = eddy5e8'
319        endif
320      enddo ! of do iq=1,nq
321
322      ! ! 3. find condensable traceurs different from h2o and n2
323      ! do iq=1,nq
324      !   if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"n2") .eq. 0)) then
325      !     count=count+1
326      !   endif
327      !   if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"n2") .eq. 0)) then
328      !     count=count+1
329      !   endif
330
331      ! enddo ! of do iq=1,nq
332
333      ! check that we identified all tracers:
334      if (count.ne.nq) then
335        write(*,*) "initracer: found only ",count," tracers"
336        write(*,*) "               expected ",nq
337        do iq=1,count
338          write(*,*)'      ',iq,' ',trim(noms(iq))
339        enddo
340!        stop
341      else
342        write(*,*) "initracer: found all expected tracers, namely:"
343        do iq=1,nq
344          write(*,*)'      ',iq,' ',trim(noms(iq))
345        enddo
346      endif
347     
348      ! Compute number of microphysics tracers:
349      ! By convention they all have the prefix "mu_" (case sensitive !)
350      nmicro = 0
351      IF (callmufi) THEN
352         DO iq=1,nq
353            str = noms(iq)
354            IF (str(1:3) == "mu_") THEN
355               nmicro = nmicro+1
356               count = count+1
357            ENDIF
358         ENDDO
359
360         ! Checking the expected number of tracers:
361         ! Microphysics moment model: nmicro = 4
362         IF (nmicro < 4) THEN
363            WRITE(*,*) "initracer:error:"," Inconsistent number of microphysical tracers"
364            WRITE(*,*) "expected at least 4 tracers,", nmicro, " given"
365            CALL abort
366         ELSE IF (nmicro > 4) THEN
367            WRITE(*,*) "!!! WARNING !!! initracer: I was expecting only four tracers, you gave me more."
368            CALL abort
369         ENDIF
370
371         ! microphysics indexes share the same values than original tracname.
372         IF (.NOT.ALLOCATED(micro_indx)) ALLOCATE(micro_indx(nmicro))
373         j = 1
374         DO i=1,nq
375            str = noms(i)
376            IF (str(1:3) == "mu_") THEN
377               micro_indx(j) = i
378               j=j+1
379            ENDIF
380         ENDDO
381     
382      ELSE
383         IF (.NOT.ALLOCATED(micro_indx)) ALLOCATE(micro_indx(nmicro))
384     
385      ENDIF ! end of callmufi
386
387      ! Get data of tracers. Need to rewind traceur.def first
388      if (is_master) then
389       rewind(407)
390       do
391        read(407,'(A)') tracline
392        if (index(tracline,"#") .ne. 1) then
393          exit
394        endif
395       enddo
396      endif
397      do iq=1,nqtot
398        if (is_master) read(407,'(A)') tracline
399        call bcast(tracline)
400        call get_tracdat(iq, tracline)
401      enddo
402
403      if (is_master) close(407)
404
405      ! Get specific data of condensable tracers
406      do iq=1,nq
407        if((is_condensable(iq)==1)) then
408                write(*,*) "There is a specie which is condensable, for generic condensation : ", noms(iq)
409                write(*,*) 'looking specie parameters for : ',noms(iq)(1:len(trim(noms(iq)))-4)
410                ! call specie_parameters_table(noms(iq)(1:len(trim(noms(iq)))-4))
411                ! constants_mass(iq)=m
412                ! constants_delta_gasH(iq)=delta_gasH
413                ! constants_Tref(iq)=Tref
414                ! constants_Pref(iq)=Pref
415                ! constants_RLVTT_generic(iq)=RLVTT_generic
416                ! constants_RCPV_generic(iq)=RCPV_generic
417        else
418                write(*,*) "This tracer is not condensable, for generic condensation :  : ", noms(iq)
419                write(*,*) "We keep condensable constants at zero"
420        endif !(is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
421      enddo ! iq=1,nq
422
423      ! Calculate number of species in the chemistry
424      nesp = sum(is_chim)
425      write(*,*) 'Number of species in the chemistry nesp = ',nesp
426
427      ! Calculate number of generic tracers
428      ngt = sum(is_condensable)
429      write(*,*) 'Number of generic tracer is  ngt = ',ngt
430
431      ! Calculate number of radiative generic condensable species
432      n_rgcs = sum(is_rgcs)
433      write(*,*)'Number of Radiative Generic Condensable Species is n_rgcs = ',n_rgcs
434      if (n_rgcs> ngt/2) then
435        write(*,*) 'You have more Radiative Generic Condensable Species than Generic Condensable Species'
436        write(*,*)'This is not possible: check your Modern traceur.def'
437        call abort_physic("initracer, issue with # of RGCS and GCS")
438      endif
439
440      ! Calculate number of microphysical tracer
441      write(*,*) 'Number of microphysical tracer nmicro = ',nmicro
442      IF (callmufi) THEN
443         call dumptracers(micro_indx)
444      ENDIF
445
446!     Processing modern traceur options
447      if(moderntracdef) then
448        call ini_recombin
449      endif
450
451!------------------------------------------------------------
452!     Initialisation tracers ....
453!------------------------------------------------------------
454      ! rho_q(1:nq)=0
455
456      rho_n2=1000.   ! N2 ice density (kg.m-3)
457
458      lw_co=274000.
459      lw_ch4=586700.
460      lw_n2=2.5e5
461      write(*,*) "lw_n2 = ", lw_n2
462
463      if (haze) then
464        ! the sedimentation radius remains radius(igcm_haze)
465        if (fractal) then
466           nmono=nb_monomer
467        else
468           nmono=1
469        endif ! end fractal
470
471        ia=0
472        if (optichaze) then
473          ia=ia+1
474          iaero_haze=ia
475          write(*,*) '--- number of haze aerosol = ', iaero_haze
476
477          block=0  ! Only one type of haze is active : the first one set in traceur.def
478          do iq=1,nq
479            tracername=noms(iq)
480            write(*,*) "--> tracername ",iq,'/',nq,' = ',tracername
481            if (tracername(1:4).eq."haze".and.block.eq.0) then
482              i_haze=iq
483              block=1
484              write(*,*) "i_haze=",i_haze
485              write(*,*) "Careful: if you set many haze traceurs in &
486                traceur.def,only ",tracername," will be radiatively active &
487                (first one in traceur.def)"
488            endif
489          enddo
490        endif ! end optichaze
491      endif ! end callmufi or haze
492
493!     Output for records:
494!     ~~~~~~~~~~~~~~~~~~
495      write(*,*)
496      Write(*,*) '******** initracer : dust transport parameters :'
497      write(*,*) 'alpha_lift = ', alpha_lift
498      write(*,*) 'alpha_devil = ', alpha_devil
499      write(*,*) 'radius  = ', radius
500      write(*,*) 'Qext  = ', qext
501      write(*,*)
502
503      contains
504
505      subroutine get_tracdat(iq,tracline)
506        !-------------------ADDING NEW OPTIONS-------------------
507        ! Duplicate if sentence for adding new options
508        ! Do not forget to add the new variables in tracer_h.F90
509        ! Do not forget to allocate and initialize the new variables above
510        ! Please update list of options in "LMDZ.GENERIC/deftank/traceur.def.modern"
511        !-------------------------------------------------------
512          use tracer_h
513          implicit none
514          integer,           intent(in) :: iq       ! tracer index
515          character(len=500),intent(in) :: tracline ! traceur.def lines with parameters
516
517          read(tracline,*) noms(iq)
518          ! JVO 20 : We should add a sanity check aborting when duplicates in names !
519          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
520          ! option mmol
521          if (index(tracline,'mmol='   ) /= 0) then
522              read(tracline(index(tracline,'mmol=')+len('mmol='):),*)&
523                  mmol(iq)
524              write(*,*) ' Parameter value (traceur.def) : mmol=', &
525                  mmol(iq)
526          else
527              write(*,*) ' Parameter value (default)     : mmol=', &
528                 mmol(iq)
529          end if
530          ! option aki
531          if (index(tracline,'aki='   ) /= 0) then
532              read(tracline(index(tracline,'aki=')+len('aki='):),*) &
533                  aki(iq)
534              write(*,*) ' Parameter value (traceur.def) : aki=', &
535                  aki(iq)
536          else
537              write(*,*) ' Parameter value (default)     : aki=', &
538                  aki(iq)
539          end if
540          ! option cpi
541          if (index(tracline,'cpi='   ) /= 0) then
542              read(tracline(index(tracline,'cpi=')+len('cpi='):),*) &
543                  cpi(iq)
544              write(*,*) ' Parameter value (traceur.def) : cpi=', &
545                  cpi(iq)
546          else
547              write(*,*) ' Parameter value (default)     : cpi=', &
548                  cpi(iq)
549          end if
550          ! option is_chim
551          if (index(tracline,'is_chim='   ) /= 0) then
552          read(tracline(index(tracline,'is_chim=')+len('is_chim='):),*) &
553                  is_chim(iq)
554              write(*,*) ' Parameter value (traceur.def) : is_chim=', &
555                  is_chim(iq)
556          else
557              write(*,*) ' Parameter value (default)     : is_chim=', &
558                  is_chim(iq)
559          end if
560          ! option is_rad
561          if (index(tracline,'is_rad=') /= 0) then
562              read(tracline(index(tracline,'is_rad=') &
563              +len('is_rad='):),*) is_rad(iq)
564              write(*,*) ' Parameter value (traceur.def) : is_rad=', &
565              is_rad(iq)
566          else
567              write(*,*) ' Parameter value (default)     : is_rad=',  &
568              is_rad(iq)
569          end if
570          ! option is_recomb
571          if (index(tracline,'is_recomb=') /= 0) then
572              read(tracline(index(tracline,'is_recomb=') &
573              +len('is_recomb='):),*) is_recomb(iq)
574              write(*,*) ' Parameter value (traceur.def) : is_recomb=', &
575              is_recomb(iq)
576          else
577              write(*,*) ' Parameter value (default)     : is_recomb=', &
578              is_recomb(iq)
579          end if
580          ! option is_recomb_qset
581          if (index(tracline,'is_recomb_qset=') /= 0) then
582              read(tracline(index(tracline,'is_recomb_qset=') &
583              +len('is_recomb_qset='):),*) is_recomb_qset(iq)
584              write(*,*) ' Parameter value (traceur.def) :'// &
585              ' is_recomb_qset=', &
586              is_recomb_qset(iq)
587          else
588              write(*,*) ' Parameter value (default)     :'// &
589              ' is_recomb_qset=', &
590              is_recomb_qset(iq)
591          endif
592          ! option is_recomb_qotf
593          if (index(tracline,'is_recomb_qotf=') /= 0) then
594              read(tracline(index(tracline,'is_recomb_qotf=') &
595              +len('is_recomb_qotf='):),*) is_recomb_qotf(iq)
596              write(*,*) ' Parameter value (traceur.def) :'// &
597              ' is_recomb_qotf=', &
598              is_recomb_qotf(iq)
599          else
600              write(*,*) ' Parameter value (default)     :'// &
601              ' is_recomb_qotf=',  &
602              is_recomb_qotf(iq)
603          end if
604          !option is_condensable (LT)
605          if (index(tracline,'is_condensable=') /=0) then
606            read(tracline(index(tracline,'is_condensable=') &
607              +len('is_condensable='):),*) is_condensable(iq)
608            write(*,*) ' Parameter value (traceur.def) :'// &
609              ' is_condensable=', is_condensable(iq)
610          else
611              write(*,*) ' Parameter value (default)     :'// &
612              ' is_condensable=', is_condensable(iq)
613          endif
614          !option radius
615          if (index(tracline,'radius=') .ne. 0) then
616            read(tracline(index(tracline,'radius=') &
617              +len('radius='):),*) radius(iq)
618            write(*,*)'Parameter value (traceur.def) :'// &
619              "radius=",radius(iq), "m"
620          else
621            write(*,*) ' Parameter value (default)     :'// &
622            ' radius=', radius(iq)," m"
623          endif
624          !option rho
625          if (index(tracline,'rho=') .ne. 0) then
626            read(tracline(index(tracline,'rho=') &
627              +len('rho='):),*) rho_q(iq)
628            write(*,*)'Parameter value (traceur.def) :'// &
629              "rho=",rho_q(iq)
630          else
631            write(*,*) ' Parameter value (default)     :'// &
632              ' rho=', rho_q(iq)
633          endif
634          !option is_rgcs
635          if (index(tracline,'is_rgcs') .ne. 0) then
636            read(tracline(index(tracline,'is_rgcs=') &
637              +len('is_rgcs='):),*) is_rgcs(iq)
638            write(*,*)'Parameter value (traceur.def) :'// &
639              'is_rgcs=',is_rgcs(iq)
640          else
641            write(*,*)'Parameter value (default) : '// &
642              'is_rgcs = ',is_rgcs(iq)
643          endif
644      end subroutine get_tracdat
645
646      end
647
Note: See TracBrowser for help on using the repository browser.