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