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

Last change on this file since 2725 was 2725, checked in by aslmd, 3 years ago

switch to epsi_generic & RLVTT_generic

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