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

Last change on this file since 3567 was 3101, checked in by emillour, 15 months ago

correction integer to real + added RCPV_generic

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