afivo-streamer 1.1
1D/2D/3D streamer simulations with AMR
Loading...
Searching...
No Matches
m_chemistry.f90
Go to the documentation of this file.
1!> Module for handling chemical reactions
3#include "../afivo/src/cpp_macros.h"
4 use m_types
5 use m_af_all
6 use m_lookup_table
8 use m_model
9 use m_random
10 use omp_lib
11
12 implicit none
13 private
14
15 !> Identifier for ionization reactions
16 integer, parameter, public :: ionization_reaction = 1
17 !> Identifier for attachment reactions
18 integer, parameter, public :: attachment_reaction = 2
19 !> Identifier for recombination reactions
20 integer, parameter, public :: recombination_reaction = 3
21 !> Identifier for detachment reactions
22 integer, parameter, public :: detachment_reaction = 4
23 !> Identifier for general reactions (not of any particular type)
24 integer, parameter, public :: general_reaction = 5
25
26 character(len=20), parameter, public :: reaction_names(*) = &
27 [character(len=20) :: "ionization", "attachment", "recombination", &
28 "detachment", "general"]
29
30 !> Maximum number of coefficients for a reaction rate function
31 integer, parameter :: rate_max_num_coeff = 9
32
33 !> Basic chemical reaction type
34 type reaction_t
35 integer, allocatable :: ix_in(:) !< Index of input species
36 integer, allocatable :: ix_out(:) !< Index of output species
37 integer, allocatable :: multiplicity_out(:) !< Multiplicity of output
38 integer :: n_species_in !< Number of input species
39 integer :: rate_type !< Type of reaction rate
40 !> Type of reaction (e.g. ionization)
41 integer :: reaction_type = general_reaction
42 real(dp) :: rate_factor !< Multiply rate by this factor
43 integer :: n_coeff !< Number of stored coefficients
44 !> Data for the reaction rate
45 real(dp) :: rate_data(rate_max_num_coeff) = -huge(1.0_dp)
46 integer :: lookup_table_index !< Index in lookup table
47 real(dp), allocatable :: x_data(:)
48 real(dp), allocatable :: y_data(:)
49 character(len=50) :: description
50 end type reaction_t
51
52 !> Compact reaction type, for efficiency
53 type tiny_react_t
54 integer, allocatable :: ix_in(:)
55 integer, allocatable :: ix_out(:)
56 integer, allocatable :: multiplicity_out(:)
57 end type tiny_react_t
58
59 !> Type for stochastic reaction
60 type stochastic_species_t
61 integer :: i_from !< Index of stochastic density (in tree)
62 integer :: i_to !< Index of regular density (in tree)
63 end type stochastic_species_t
64
65 !> Reaction with a field-dependent reaction rate
66 integer, parameter :: rate_tabulated_energy = 0
67
68 !> Reaction with a field-dependent reaction rate
69 integer, parameter :: rate_tabulated_field = 1
70
71 !> Reaction of the form c1
72 integer, parameter :: rate_analytic_constant = 2
73
74 !> Reaction of the form c1 * (Td - c2)
75 integer, parameter :: rate_analytic_linear = 3
76
77 !> Reaction of the form c1 * exp(-(c2/(c3 + Td))**2)
78 integer, parameter :: rate_analytic_exp_v1 = 4
79
80 !> Reaction of the form c1 * exp(-(Td/c2)**2)
81 integer, parameter :: rate_analytic_exp_v2 = 5
82
83 !> Reaction of the form c1 * (300 / Te)**c2
84 integer, parameter :: rate_analytic_k1 = 6
85
86 !> Reaction of the form (c1 * (kB_eV * Te + c2)**2 - c3) * c4
87 integer, parameter :: rate_analytic_k3 = 8
88
89 !> Reaction of the form c1 * (Tg / 300)**c2 * exp(-c3 / Tg)
90 integer, parameter :: rate_analytic_k4 = 9
91
92 !> Reaction of the form c1 * exp(-c2 / Tg)
93 integer, parameter :: rate_analytic_k5 = 10
94
95 !> Reaction of the form c1 * Tg**c2
96 integer, parameter :: rate_analytic_k6 = 11
97
98 !> Reaction of the form c1 * (Tg / c2)**c3
99 integer, parameter :: rate_analytic_k7 = 12
100
101 !> Reaction of the form c1 * (300 / Tg)**c2
102 integer, parameter :: rate_analytic_k8 = 13
103
104 !> Reaction of the form c1 * exp(-c2 * Tg)
105 integer, parameter :: rate_analytic_k9 = 14
106
107 !> Reaction of the form 10**(c1 + c2 * (Tg - 300))
108 integer, parameter :: rate_analytic_k10 = 15
109
110 !> Reaction of the form c1 * (300 / Tg)**c2 * exp(-c3 / Tg)
111 integer, parameter :: rate_analytic_k11 = 16
112
113 !> Reaction of the form c1 * Tg**c2 * exp(-c3 / Tg)
114 integer, parameter :: rate_analytic_k12 = 17
115
116 !> Reaction of the form c1 * exp(-(c2 / (c3 + Td))**c4)
117 integer, parameter :: rate_analytic_k13 = 18
118
119 !> Reaction of the form c1 * exp(-(Td / c2)**c3)
120 integer, parameter :: rate_analytic_k14 = 19
121
122 !> Reaction of the form c1 * exp(-(c2 /(kb * (Tg + Td/c3)))**c4)
123 integer, parameter :: rate_analytic_k15 = 20
124
125 !> Maximum number of species
126 integer, parameter :: max_num_species = 100
127
128 !> Maximum number of reactions
129 integer, parameter :: max_num_reactions = 500
130
131 !> Number of species present
132 integer, public, protected :: n_species = 0
133
134 !> Number of gas species present
135 integer, public, protected :: n_gas_species = 0
136
137 !> Number of plasma species present
138 integer, public, protected :: n_plasma_species = 0
139
140 !> Number of reactions present
141 integer, public, protected :: n_reactions = 0
142
143 !> Number of stochastic species
144 integer, public, protected :: n_stochastic_species = 0
145
146 !> List of the species
147 character(len=comp_len), public, protected :: species_list(max_num_species)
148
149 !> List of stochastic species
150 type(stochastic_species_t) :: stochastic_species_list(max_num_species)
151
152 !> Charge of the species
153 integer, public, protected :: species_charge(max_num_species) = 0
154
155 !> species_itree(n) holds the index of species n in the tree (cell-centered variables)
156 integer, public, protected :: species_itree(max_num_species)
157
158 !> List of reactions
159 type(reaction_t), public, protected :: reactions(max_num_reactions)
160
161 !> A copy of the list of reactions for performance reasons
162 type(tiny_react_t) :: tiny_react(max_num_reactions)
163
164 !> Lookup table with reaction rates versus field
165 type(lt_t) :: chemtbl_fld
166
167 !> Lookup table with reaction rates versus electron energy
168 type(lt_t) :: chemtbl_ee
169
170 !> List with indices of charged species
171 integer, allocatable, protected :: charged_species_itree(:)
172
173 !> List with charges of charged species
174 integer, allocatable, protected :: charged_species_charge(:)
175
176 public :: charged_species_itree
177 public :: charged_species_charge
178
179 public :: chemistry_initialize
183 public :: get_rates
184 public :: get_derivatives
185 public :: species_index
186
187 public :: read_reactions
188
189contains
190
191 !> Initialize module and load chemical reactions
192 subroutine chemistry_initialize(tree, cfg)
193 use m_config
195 use m_table_data
197 use m_gas
198 use m_dt
199 type(af_t), intent(inout) :: tree
200 type(cfg_t), intent(inout) :: cfg
201 integer :: n, i, j, i_elec, rtype
202 character(len=string_len) :: reaction_file
203 character(len=comp_len) :: tmp_name
204 logical :: read_success
205
206 call cfg_get(cfg, "input_data%file", reaction_file)
207
208 if (.not. gas_constant_density) then
209 ! Make sure the gas components are the first species
213 end if
214
215 call read_reactions(trim(reaction_file), read_success)
216
217 if (.not. read_success) then
218 print *, "m_chemistry: no reaction table found, using standard model"
219
220 species_list(1) = "e"
221 species_list(2) = "M+"
222 species_list(3) = "M-"
223 n_species = 3
224 n_reactions = 2
225
226 ! Ionization reaction
227 if (gas_constant_density) then
228 reactions(1)%ix_in = [1]
229 reactions(1)%ix_out = [1, 2]
230 reactions(1)%multiplicity_out = [2, 1]
231 reactions(1)%n_species_in = 2
232 reactions(1)%rate_type = rate_tabulated_field
233 reactions(1)%rate_factor = 1.0_dp
234 reactions(1)%x_data = td_tbl%x
235 reactions(1)%y_data = td_tbl%rows_cols(:, td_alpha) * &
236 td_tbl%rows_cols(:, td_mobility) * reactions(1)%x_data * &
238 reactions(1)%description = "e + M -> e + e + M+"
239
240 ! Attachment reaction
241 reactions(2)%ix_in = [1]
242 reactions(2)%ix_out = [3]
243 reactions(2)%multiplicity_out = [1]
244 reactions(2)%n_species_in = 2
245 reactions(2)%rate_type = rate_tabulated_field
246 reactions(2)%rate_factor = 1.0_dp
247 reactions(2)%x_data = td_tbl%x
248 reactions(2)%y_data = td_tbl%rows_cols(:, td_eta) * &
249 td_tbl%rows_cols(:, td_mobility) * reactions(2)%x_data * &
251 reactions(2)%description = "e + M -> M-"
252 else
253 error stop "Varying gas density not yet supported"
254 end if
255 end if
256
257 ! In case of energy equation, add an extra species
259 n_species = n_species + 1
260 species_list(n_species) = "e_energy"
261 end if
262
263 ! Convert names to simple ascii
264 do n = 1, n_species
265 tmp_name = species_list(n)
266 call to_simple_ascii(trim(tmp_name), species_list(n), &
268 end do
269
270 ! Also store in more memory-efficient structure
271 do n = 1, n_reactions
272 tiny_react(n)%ix_in = reactions(n)%ix_in
273 tiny_react(n)%ix_out = reactions(n)%ix_out
274 tiny_react(n)%multiplicity_out = reactions(n)%multiplicity_out
275 end do
276
277 ! Gas species are not stored in the tree for now
280
281 do n = n_gas_species+1, n_species
282 call af_add_cc_variable(tree, trim(species_list(n)), &
283 n_copies=af_advance_num_steps(time_integrator)+1, &
284 ix=species_itree(n))
285 end do
286
287 ! Store list with only charged species
288 n = count(species_charge(1:n_species) /= 0)
289 allocate(charged_species_itree(n))
290 allocate(charged_species_charge(n))
291
292 i = 0
293 do n = 1, n_species
294 if (species_charge(n) /= 0) then
295 i = i + 1
298 end if
299 end do
300
301 call check_charge_conservation()
302
303 ! Specify reaction types
304 i_elec = species_index("e")
305
306 do n = 1, n_reactions
307 if (any(reactions(n)%ix_in == i_elec) .and. &
308 .not. any(reactions(n)%ix_out == i_elec) .and. &
309 .not. any(species_charge(reactions(n)%ix_in) > 0)) then
310 ! In: an electron and no positive ions, out: no electrons
311 reactions(n)%reaction_type = attachment_reaction
312 else if (any(reactions(n)%ix_in == i_elec) .and. &
313 any(reactions(n)%ix_out == i_elec .and. &
314 reactions(n)%multiplicity_out == 2)) then
315 ! An electron in and out (with multiplicity 2)
316 reactions(n)%reaction_type = ionization_reaction
317 else if (any(species_charge(reactions(n)%ix_in) /= 0) .and. &
318 .not. any(species_charge(reactions(n)%ix_out) /= 0)) then
319 ! In: charged species, out: no charged species
320 reactions(n)%reaction_type = recombination_reaction
321 else if (all(reactions(n)%ix_in /= i_elec) .and. &
322 any(reactions(n)%ix_out == i_elec)) then
323 ! In: no electrons, out: an electron
324 reactions(n)%reaction_type = detachment_reaction
325 end if
326 end do
327
328 ! Create lookup tables for tabulated reaction data. First determine number
329 ! of reactions. In case of an energy equation, ionization and attachment
330 ! reactions are converted to use the electron energy instead of the field.
331 i = 0
332 j = count(reactions(1:n_reactions)%rate_type == rate_tabulated_energy)
333 do n = 1, n_reactions
334 if (reactions(n)%rate_type == rate_tabulated_field) then
335 rtype = reactions(n)%reaction_type
337 .or. rtype == attachment_reaction)) then
338 j = j + 1
339 else
340 i = i + 1
341 end if
342 end if
343 end do
344
345 chemtbl_fld = lt_create(td_tbl%x(1), td_tbl%x(td_tbl%n_points), &
347 chemtbl_ee = lt_create(0.0_dp, td_max_ev, &
349
350 i = 0
351 j = 0
352 do n = 1, n_reactions
353 if (reactions(n)%rate_type == rate_tabulated_field) then
354 rtype = reactions(n)%reaction_type
355
357 .or. rtype == attachment_reaction)) then
358 reactions(n)%rate_type = rate_tabulated_energy
359 j = j + 1
360 reactions(n)%lookup_table_index = j
361 ! Convert field to energy
362 call table_set_column(chemtbl_ee, j, &
363 lt_get_col(td_tbl, td_energy_ev, reactions(n)%x_data), &
364 reactions(n)%y_data)
365 else
366 i = i + 1
367 reactions(n)%lookup_table_index = i
368 call table_set_column(chemtbl_fld, i, reactions(n)%x_data, &
369 reactions(n)%y_data)
370 end if
371 else if (reactions(n)%rate_type == rate_tabulated_energy) then
372 j = j + 1
373 reactions(n)%lookup_table_index = j
374 call table_set_column(chemtbl_ee, j, reactions(n)%x_data, &
375 reactions(n)%y_data)
376 end if
377 end do
378
379 call store_stochastic_species()
380
381 print *, "--- List of reactions ---"
382 do n = 1, n_reactions
383 write(*, "(I4,' (',I0,') ',A15,A)") n, reactions(n)%n_species_in, &
384 reaction_names(reactions(n)%reaction_type), &
385 reactions(n)%description
386 end do
387 print *, "-------------------------"
388 print *, ""
389 print *, "--- List of gas species ---"
390 do n = 1, n_gas_species
391 write(*, "(I4,A)") n, ": " // species_list(n)
392 end do
393 print *, "-------------------------"
394 print *, ""
395 print *, "--- List of plasma species ---"
396 do n = n_gas_species+1, n_species
397 write(*, "(I4,A)") n, ": " // species_list(n)
398 end do
399 print *, "-------------------------"
400
401 ! Optionally modify rates
402 call chemistry_modify_rates(cfg)
403
404 end subroutine chemistry_initialize
405
406 !> Store stochastic species, for which species have a prefix "stoch_"
407 subroutine store_stochastic_species()
408 integer :: n, i_stoch, namelen, iv
409
410 i_stoch = 0
411 do n = n_gas_species+1, n_species
412 namelen = len_trim(species_list(n))
413
414 if (species_list(n)(1:6) == "stoch_") then
415 i_stoch = i_stoch + 1
416
417 stochastic_species_list(i_stoch)%i_from = species_itree(n)
418 iv = species_index(species_list(n)(7:))
419
420 if (iv == -1) then
421 print *, "Could not find [", species_list(n)(7:), "]"
422 print *, "which should be present for stochastic species (stoch_)"
423 error stop "Species for stochastic reaction not found"
424 else
425 stochastic_species_list(i_stoch)%i_to = species_itree(iv)
426 end if
427 end if
428 end do
429
430 n_stochastic_species = i_stoch
431
432 end subroutine store_stochastic_species
433
434 !> Convert stochastic species into regular species, using random numbers
436 type(af_t), intent(inout) :: tree !< Tree
437 type(rng_t), intent(inout) :: rng !< Random number generator
438 integer :: lvl, iblock, id, ijk, n
439 integer :: i_to, i_from, proc_id, n_procs
440 real(dp) :: dv, lambda, product_dr, n_produced
441#if NDIM == 2
442 real(dp), parameter :: two_pi = 2 * acos(-1.0_dp)
443#endif
444 type(prng_t) :: prng
445
446 ! Initialize parallel random numbers
447 n_procs = omp_get_max_threads()
448 call prng%init_parallel(n_procs, rng)
449
450 !$omp parallel private(lvl, iblock, id, proc_id, dV, i_to, &
451 !$omp &i_from, n, n_produced, IJK, lambda, product_dr)
452 proc_id = 1+omp_get_thread_num()
453
454 do lvl = 1, tree%highest_lvl
455 !$omp do
456 do iblock = 1, size(tree%lvls(lvl)%leaves)
457 id = tree%lvls(lvl)%leaves(iblock)
458 product_dr = product(tree%boxes(id)%dr)
459
460 do n = 1, n_stochastic_species
461 i_to = stochastic_species_list(n)%i_to
462 i_from = stochastic_species_list(n)%i_from
463
464 do kji_do(1, tree%n_cell)
465#if NDIM == 2
466 if (tree%coord_t == af_cyl) then
467 dv = two_pi * product_dr * af_cyl_radius_cc(tree%boxes(id), i)
468 else
469 dv = product_dr
470 end if
471#else
472 dv = product_dr
473#endif
474
475 lambda = dv * tree%boxes(id)%cc(ijk, i_from)
476 n_produced = prng%rngs(proc_id)%poisson(lambda)
477
478 tree%boxes(id)%cc(ijk, i_from) = 0.0_dp
479 tree%boxes(id)%cc(ijk, i_to) = tree%boxes(id)%cc(ijk, i_to) + &
480 n_produced/dv
481 end do; close_do
482 end do
483
484 end do
485 !$omp end do
486 end do
487 !$omp end parallel
488
489 call prng%update_seed(rng)
490
492
493 !> Modify reaction rates for sensitivity analysis
494 subroutine chemistry_modify_rates(cfg)
495 use m_config
496 type(cfg_t), intent(inout) :: cfg
497 integer :: n_modified, n, dummy_int(0), ix
498 real(dp) :: dummy_real(0)
499 integer, allocatable :: reaction_ix(:)
500 real(dp), allocatable :: rate_factors(:)
501
502 call cfg_add(cfg, "input_data%modified_reaction_ix", dummy_int, &
503 "Indices of reactions to be modified", .true.)
504 call cfg_add(cfg, "input_data%modified_rate_factors", dummy_real, &
505 "Reaction rate factors for modified reactions", .true.)
506 call cfg_get_size(cfg, "input_data%modified_reaction_ix", n_modified)
507 call cfg_get_size(cfg, "input_data%modified_rate_factors", n)
508
509 if (n /= n_modified) &
510 error stop "size(modified_reaction_ix) /= size(modified_rate_factors)"
511
512 if (n_modified > 0) then
513 allocate(reaction_ix(n_modified), rate_factors(n_modified))
514 call cfg_get(cfg, "input_data%modified_reaction_ix", reaction_ix)
515 call cfg_get(cfg, "input_data%modified_rate_factors", rate_factors)
516
517 if (minval(rate_factors) < 0) &
518 error stop "Negative value in modified_rate_factors"
519 if (minval(reaction_ix) < 1 .or. maxval(reaction_ix) > n_reactions) &
520 error stop "modified_reaction_ix outside valid range"
521
522 do n = 1, n_modified
523 ix = reaction_ix(n)
524 reactions(ix)%rate_factor = reactions(ix)%rate_factor * rate_factors(n)
525 end do
526 end if
527
528 end subroutine chemistry_modify_rates
529
530 !> Write a summary of the reactions (TODO) and the ionization and attachment
531 !> coefficients (if working at constant pressure)
532 subroutine chemistry_write_summary(fname)
533 use m_gas
535 character(len=*), intent(in) :: fname
536 real(dp), allocatable :: fields(:), energies(:)
537 real(dp), allocatable :: rates(:, :)
538 real(dp), allocatable :: eta(:), alpha(:), src(:), loss(:)
539 real(dp), allocatable :: v(:), mu(:), diff(:)
540 integer :: n, n_fields, i_elec
541 integer :: my_unit
542
543 if (gas_constant_density) then
544 i_elec = species_index("e")
545 n_fields = td_tbl%n_points
546
547 if (n_fields < 3) error stop "Not enough data for linear extrapolation"
548
549 allocate(fields(n_fields))
550 allocate(energies(n_fields))
551 allocate(rates(n_fields, n_reactions))
552 allocate(eta(n_fields), alpha(n_fields), src(n_fields), loss(n_fields))
553
554 fields = td_tbl%x
555
556 if (model_has_energy_equation) then
557 energies = lt_get_col(td_tbl, td_energy_ev, td_tbl%x)
558 call get_rates(fields, rates, n_fields, energies)
559 else
560 call get_rates(fields, rates, n_fields)
561 end if
562
563 loss(:) = 0.0_dp
564 src(:) = 0.0_dp
565
566 do n = 1, n_reactions
567 if (reactions(n)%reaction_type == attachment_reaction) then
568 loss(:) = loss(:) + rates(:, n)
569 else if (reactions(n)%reaction_type == ionization_reaction) then
570 src(:) = src(:) + rates(:, n)
571 end if
572 end do
573
574 allocate(diff(n_fields))
575 diff = lt_get_col(td_tbl, td_diffusion, fields)
576
577 allocate(mu(n_fields))
578 allocate(v(n_fields))
579 mu = lt_get_col(td_tbl, td_mobility, fields)
580 v = mu * fields * townsend_to_si
581
582 ! v(1) is zero, so extrapolate linearly
583 eta(2:) = loss(2:) / v(2:)
584 eta(1) = 2 * eta(2) - eta(3)
585 alpha(2:) = src(2:) / v(2:)
586 alpha(1) = 2 * alpha(2) - alpha(3)
587
588 ! Write to a file
589 open(newunit=my_unit, file=trim(fname), action="write")
590 write(my_unit, "(A)") &
591 "E/N[Td] E[V/m] Electron_mobility[m^2/(Vs)] &
592 &Electron_diffusion[m^2/s] Townsend_ioniz._coef._alpha[1/m] &
593 &Townsend_attach._coef._eta[1/m] Ionization_rate[1/s] &
594 &Attachment_rate[1/s]"
595 do n = 1, n_fields
596 write(my_unit, *) fields(n), &
597 fields(n) * townsend_to_si * gas_number_density, &
598 mu(n) / gas_number_density, &
599 diff(n) / gas_number_density, alpha(n), eta(n), &
600 src(n), loss(n)
601 end do
602 write(my_unit, *) ""
603 close(my_unit)
604 end if
605 end subroutine chemistry_write_summary
606
607 subroutine check_charge_conservation()
608 integer :: n, q_in, q_out
609
610 do n = 1, n_reactions
611 q_in = sum(species_charge(reactions(n)%ix_in))
612 q_out = sum(species_charge(reactions(n)%ix_out) * &
613 reactions(n)%multiplicity_out)
614 if (q_in /= q_out) then
615 print *, trim(reactions(n)%description)
616 error stop "Charge is not conserved in the above reaction"
617 end if
618 end do
619 end subroutine check_charge_conservation
620
621 !> Get the breakdown field in Townsend
622 subroutine chemistry_get_breakdown_field(field_td, min_growth_rate)
624 !> Breakdown field in Townsend
625 real(dp), intent(out) :: field_td
626 !> Minimal growth rate for breakdown
627 real(dp), intent(in) :: min_growth_rate
628
629 integer :: n, n_fields
630 real(dp), allocatable :: energies(:)
631 real(dp), allocatable :: fields(:), rates(:, :), src(:), loss(:)
632
633 n_fields = td_tbl%n_points
634 allocate(fields(n_fields))
635 allocate(energies(n_fields))
636 allocate(rates(n_fields, n_reactions))
637 allocate(src(n_fields), loss(n_fields))
638
639 fields = td_tbl%x
640 if (model_has_energy_equation) then
641 energies = lt_get_col(td_tbl, td_energy_ev, td_tbl%x)
642 call get_rates(fields, rates, n_fields, energies)
643 else
644 call get_rates(fields, rates, n_fields)
645 end if
646
647 loss(:) = 0.0_dp
648 src(:) = 0.0_dp
649
650 do n = 1, n_reactions
651 if (reactions(n)%reaction_type == attachment_reaction) then
652 loss(:) = loss(:) + rates(:, n)
653 else if (reactions(n)%reaction_type == ionization_reaction) then
654 src(:) = src(:) + rates(:, n)
655 end if
656 end do
657
658 do n = n_fields, 1, -1
659 if (src(n) - loss(n) < min_growth_rate) exit
660 end do
661
662 field_td = 0.0_dp
663 if (n > 0) field_td = fields(n)
664 end subroutine chemistry_get_breakdown_field
665
666 !> Compute reaction rates
667 !>
668 !> @todo These reactions do not take into account a variable gas_temperature
669 subroutine get_rates(fields, rates, n_cells, energy_eV)
671 use m_gas
673 integer, intent(in) :: n_cells !< Number of cells
674 real(dp), intent(in) :: fields(n_cells) !< The field (in Td) in the cells
675 real(dp), intent(out) :: rates(n_cells, n_reactions) !< The reaction rates
676 real(dp), intent(in), optional :: energy_ev(n_cells) !< Electron energy in eV
677 integer :: n, n_coeff
678 real(dp) :: c0, c(rate_max_num_coeff)
679 real(dp) :: te(n_cells) !> Electron Temperature in Kelvin
680 logical :: te_available
681
682 ! Conversion factor to go from eV to Kelvin
683 real(dp), parameter :: electron_ev_to_k = 2 * uc_elec_volt / &
685
686 te_available = .false.
687
688 do n = 1, n_reactions
689 ! A factor that the reaction rate is multiplied with, for example to
690 ! account for a constant gas number density and if the length unit is cm instead of m
691 c0 = reactions(n)%rate_factor
692
693 ! Coefficients for the reaction rate
694 n_coeff = reactions(n)%n_coeff
695 c(1:n_coeff) = reactions(n)%rate_data(1:n_coeff)
696
697 select case (reactions(n)%rate_type)
699 if (.not. present(energy_ev)) error stop "energy_eV required"
700 rates(:, n) = c0 * lt_get_col(chemtbl_ee, &
701 reactions(n)%lookup_table_index, energy_ev)
702 case (rate_tabulated_field)
703 rates(:, n) = c0 * lt_get_col(chemtbl_fld, &
704 reactions(n)%lookup_table_index, fields)
705 case (rate_analytic_constant)
706 rates(:, n) = c0 * c(1)
707 case (rate_analytic_linear)
708 rates(:, n) = c0 * c(1) * (fields - c(2))
709 case (rate_analytic_exp_v1)
710 rates(:, n) = c0 * c(1) * exp(-(c(2) / (c(3) + fields))**2)
711 case (rate_analytic_exp_v2)
712 rates(:, n) = c0 * c(1) * exp(-(fields/c(2))**2)
713 case (rate_analytic_k1)
714 if (.not. te_available) then
715 ! Note that we could use energy_eV if present, but this energy is
716 ! not guaranteed to be well-behaved
717 te = electron_ev_to_k * lt_get_col(td_tbl, td_energy_ev, fields)
718 end if
719 rates(:, n) = c0 * c(1) * (300 / te)**c(2)
720 case (rate_analytic_k3)
721 if (.not. te_available) then
722 te = electron_ev_to_k * lt_get_col(td_tbl, td_energy_ev, fields)
723 end if
724 ! We convert boltzmann_const from J / K to eV / K
725 rates(:, n) = c0 * (c(1) * ((uc_boltzmann_const / uc_elec_volt) * te + c(2))**2 - c(3)) * c(4)
726 case (rate_analytic_k4)
727 rates(:, n) = c0 * c(1) * (gas_temperature / 300)**c(2) * exp(-c(3) / gas_temperature)
728 case (rate_analytic_k5)
729 rates(:, n) = c0 * c(1) * exp(-c(2) / gas_temperature)
730 case (rate_analytic_k6)
731 rates(:, n) = c0 * c(1) * gas_temperature**c(2)
732 case (rate_analytic_k7)
733 rates(:, n) = c0 * c(1) * (gas_temperature / c(2))**c(3)
734 case (rate_analytic_k8)
735 rates(:, n) = c0 * c(1) * (300 / gas_temperature)**c(2)
736 case (rate_analytic_k9)
737 rates(:, n) = c0 * c(1) * exp(-c(2) * gas_temperature)
738 case (rate_analytic_k10)
739 rates(:, n) = c0 * 10**(c(1) + c(2) * (gas_temperature - 300))
740 case (rate_analytic_k11)
741 rates(:, n) = c0 * c(1) * (300 / gas_temperature)**c(2) * exp(-c(3) / gas_temperature)
742 case (rate_analytic_k12)
743 rates(:, n) = c0 * c(1) * gas_temperature**c(2) * exp(-c(3) / gas_temperature)
744 case (rate_analytic_k13)
745 rates(:, n) = c0 * c(1) * exp(-(c(2) / (c(3) + fields))**c(4))
746 case (rate_analytic_k14)
747 rates(:, n) = c0 * c(1) * exp(-(fields / c(2))**c(3))
748 case (rate_analytic_k15)
749 ! Note that this reaction depends on Ti, ionic temperature.
750 ! According to Galimberti(1979): Ti = T_gas + fields/c(3),
751 ! with c(3) = 0.18 Td/Kelvin, UC_boltzmann_const is in J/Kelvin, c(2)
752 ! is given in Joule in the input file
753 rates(:, n) = c0 * c(1) * exp(-(c(2) / (uc_boltzmann_const * &
754 (gas_temperature + fields/c(3))))**c(4))
755 end select
756 end do
757 end subroutine get_rates
758
759 !> Compute derivatives due to chemical reactions. Note that the 'rates'
760 !> argument is modified.
761 subroutine get_derivatives(dens, rates, derivs, n_cells)
762 integer, intent(in) :: n_cells
763 !> Species densities
764 real(dp), intent(in) :: dens(n_cells, n_species)
765 !> On input, reaction rate coefficients. On output, actual reaction rates.
766 real(dp), intent(inout) :: rates(n_cells, n_reactions)
767 !> Derivatives of the chemical species
768 real(dp), intent(out) :: derivs(n_cells, n_species)
769 integer :: n, i, ix
770
771 derivs(:, :) = 0.0_dp
772
773 ! Loop over reactions and add to derivative
774 do n = 1, n_reactions
775 ! Determine production rate of 'full' reaction
776 rates(:, n) = rates(:, n) * &
777 product(dens(:, tiny_react(n)%ix_in), dim=2)
778
779 ! Input species are removed
780 do i = 1, size(tiny_react(n)%ix_in)
781 ix = tiny_react(n)%ix_in(i)
782 derivs(:, ix) = derivs(:, ix) - rates(:, n)
783 end do
784
785 ! Output species are created
786 do i = 1, size(tiny_react(n)%ix_out)
787 ix = tiny_react(n)%ix_out(i)
788 derivs(:, ix) = derivs(:, ix) + rates(:, n) * &
789 tiny_react(n)%multiplicity_out(i)
790 end do
791 end do
792 end subroutine get_derivatives
793
794 !> Try to read a list species whose production will be ignored
795 subroutine read_ignored_species(filename, ignored_species)
796 character(len=*), intent(in) :: filename
797 character(len=comp_len), allocatable, intent(inout) :: &
798 ignored_species(:)
799 character(len=comp_len) :: tmp(max_num_species)
800 character(len=string_len) :: line
801 integer :: my_unit, n_ignored
802
803 n_ignored = 0
804
805 open(newunit=my_unit, file=filename, action="read")
806
807 ! Find list of reactions
808 do
809 read(my_unit, "(A)", end=998) line
810 line = adjustl(line)
811
812 if (line == "ignored_species") then
813 ! Read next line starting with at least 5 dashes
814 read(my_unit, "(A)") line
815 if (line(1:5) /= "-----") &
816 error stop "ignored_species not followed by -----"
817 exit
818 end if
819 end do
820
821 ! Read ignored species, one per line
822 do
823 read(my_unit, "(A)", end=999) line
824 line = adjustl(line)
825
826 ! Ignore comments
827 if (line(1:1) == "#") cycle
828
829 ! Exit when we read a line of dashes
830 if (line(1:5) == "-----") exit
831
832 n_ignored = n_ignored + 1
833 read(line, *) tmp(n_ignored)
834 end do
835
836998 close(my_unit)
837 ignored_species = tmp(1:n_ignored)
838 return
839
840 ! Error messages
841999 error stop "read_ignored_species: no closing dashes"
842 end subroutine read_ignored_species
843
844 !> Read reactions from a file
845 subroutine read_reactions(filename, read_success)
846 character(len=*), intent(in) :: filename
847 logical, intent(out) :: read_success
848 character(len=string_len) :: line
849 integer, parameter :: field_len = 50
850 character(len=field_len) :: data_value(max_num_reactions)
851 character(len=field_len) :: reaction(max_num_reactions)
852 character(len=field_len) :: how_to_get(max_num_reactions)
853 character(len=10) :: length_unit(max_num_reactions)
854 type(reaction_t) :: new_reaction
855 integer :: my_unit
856 integer :: n_reactions_found
857 integer, parameter :: n_fields_max = 40
858 integer :: i0(n_fields_max), i1(n_fields_max)
859 integer :: n, i, k, n_found, lo, hi
860 logical :: keep_reaction
861
862 character(len=comp_len), allocatable :: ignored_species(:)
863
864 type group
865 character(len=name_len) :: name
866 character(len=field_len), allocatable :: members(:)
867 end type group
868
869 integer, parameter :: max_groups = 10
870 integer :: i_group, group_size
871 type(group) :: groups(max_groups)
872
873 call read_ignored_species(filename, ignored_species)
874
875 open(newunit=my_unit, file=filename, action="read")
876
877 n_reactions_found = 0
878 i_group = 0
879 group_size = 0
880 read_success = .false.
881
882 ! Find list of reactions
883 do
884 read(my_unit, "(A)", end=998) line
885 line = adjustl(line)
886
887 if (line == "reaction_list") then
888 ! Read next line starting with at least 5 dashes
889 read(my_unit, "(A)") line
890 if (line(1:5) /= "-----") &
891 error stop "read_reactions: reaction_list not followed by -----"
892 exit
893 end if
894 end do
895
896 ! Start reading reactions
897 do
898 read(my_unit, "(A)", end=999) line
899 line = adjustl(line)
900
901 ! Ignore comments
902 if (line(1:1) == "#") cycle
903
904 ! Exit when we read a line of dashes
905 if (line(1:5) == "-----") exit
906
907 ! Group syntax for multiple species
908 if (line(1:1) == "@") then
909 call get_fields_string(line, "=,", n_fields_max, n_found, i0, i1)
910 i_group = i_group + 1
911
912 if (i_group > max_groups) error stop "Too many groups"
913
914 if (i_group == 1) then
915 group_size = n_found - 1
916 else if (n_found - 1 /= group_size) then
917 print *, trim(line)
918 error stop "Groups for a reaction should have the same size"
919 end if
920
921 groups(i_group)%name = line(i0(1):i1(1))
922 allocate(groups(i_group)%members(n_found - 1))
923
924 do n = 2, n_found
925 groups(i_group)%members(n-1) = adjustl(line(i0(n):i1(n)))
926 end do
927 cycle
928 end if
929
930 if (i_group > 0) then
931 ! Handle groups
932 lo = n_reactions_found
933 hi = n_reactions_found+group_size-1
934 reaction(lo+1:hi) = reaction(lo)
935 how_to_get(lo+1:hi) = how_to_get(lo)
936 data_value(lo+1:hi) = data_value(lo)
937 length_unit(lo+1:hi) = length_unit(lo)
938 n_reactions_found = hi
939
940 do k = 1, group_size
941 do i = 1, i_group
942 n = lo + k - 1
943 reaction(n) = string_replace(reaction(n), &
944 groups(i)%name, groups(i)%members(k))
945 how_to_get(n) = string_replace(how_to_get(n), &
946 groups(i)%name, groups(i)%members(k))
947 data_value(n) = string_replace(data_value(n), &
948 groups(i)%name, groups(i)%members(k))
949 end do
950 end do
951
952 ! Check if there are duplicate reactions
953 do n = lo, hi
954 if (count(reaction(lo:hi) == reaction(n)) > 1 .and. &
955 count(data_value(lo:hi) == data_value(n)) > 1) then
956 do k = lo, hi
957 print *, trim(reaction(k)), ",", data_value(k)
958 end do
959 error stop "Groups lead to duplicate reactions"
960 end if
961 end do
962
963 i_group = 0
964 group_size = 0
965 end if
966
967 call get_fields_string(line, ",", n_fields_max, n_found, i0, i1)
968
969 if (n_found < 3 .or. n_found > 4) then
970 print *, trim(line)
971 error stop "Invalid chemistry syntax"
972 end if
973
974 if (n_reactions_found >= max_num_reactions) &
975 error stop "Too many reactions, increase max_num_reactions"
976
977 n_reactions_found = n_reactions_found + 1
978 reaction(n_reactions_found) = line(i0(1):i1(1))
979 how_to_get(n_reactions_found) = line(i0(2):i1(2))
980 data_value(n_reactions_found) = line(i0(3):i1(3))
981
982 ! Fourth entry can hold a custom length unit, default is meter
983 if (n_found > 3) then
984 length_unit(n_reactions_found) = line(i0(4):i1(4))
985 else
986 length_unit(n_reactions_found) = "m"
987 end if
988 end do
989
990998 continue
991
992 ! Close the file (so that we can re-open it for reading data)
993 close(my_unit)
994
995 n_reactions = 0
996
997 ! Now parse the reactions
998 do n = 1, n_reactions_found
999 call parse_reaction(trim(reaction(n)), new_reaction, &
1000 ignored_species, keep_reaction)
1001
1002 if (keep_reaction) then
1004 else
1005 cycle
1006 end if
1007
1008 new_reaction%description = trim(reaction(n))
1009
1010 ! IMPORTANT: If you change the reactions below, do not forget to update
1011 ! documentation/chemistry.md accordingly!
1012 select case (how_to_get(n))
1013 case ("field_table")
1014 ! Reaction data should be present in the same file
1015 call read_reaction_table(filename, &
1016 trim(data_value(n)), new_reaction)
1017 new_reaction%n_coeff = 0
1018 case ("c1")
1019 new_reaction%rate_type = rate_analytic_constant
1020 new_reaction%n_coeff = 1
1021 read(data_value(n), *) new_reaction%rate_data(1)
1022 case ("c1*(Td-c2)")
1023 new_reaction%rate_type = rate_analytic_linear
1024 new_reaction%n_coeff = 2
1025 read(data_value(n), *) new_reaction%rate_data(1:2)
1026 case ("c1*exp(-(c2/(c3+Td))**2)")
1027 new_reaction%rate_type = rate_analytic_exp_v1
1028 new_reaction%n_coeff = 3
1029 read(data_value(n), *) new_reaction%rate_data(1:3)
1030 case ("c1*exp(-(Td/c2)**2)")
1031 new_reaction%rate_type = rate_analytic_exp_v2
1032 new_reaction%n_coeff = 2
1033 read(data_value(n), *) new_reaction%rate_data(1:2)
1034 case ("c1*(300/Te)**c2")
1035 new_reaction%rate_type = rate_analytic_k1
1036 new_reaction%n_coeff = 2
1037 read(data_value(n), *) new_reaction%rate_data(1:2)
1038 case ("(c1*(kB_eV*Te+c2)**2-c3)*c4")
1039 new_reaction%rate_type = rate_analytic_k3
1040 new_reaction%n_coeff = 4
1041 read(data_value(n), *) new_reaction%rate_data(1:4)
1042 case ("c1*(Tg/300)**c2*exp(-c3/Tg)")
1043 new_reaction%rate_type = rate_analytic_k4
1044 new_reaction%n_coeff = 3
1045 read(data_value(n), *) new_reaction%rate_data(1:3)
1046 case ("c1*exp(-c2/Tg)")
1047 new_reaction%rate_type = rate_analytic_k5
1048 new_reaction%n_coeff = 2
1049 read(data_value(n), *) new_reaction%rate_data(1:2)
1050 case ("c1*Tg**c2")
1051 new_reaction%rate_type = rate_analytic_k6
1052 new_reaction%n_coeff = 2
1053 read(data_value(n), *) new_reaction%rate_data(1:2)
1054 case ("c1*(Tg/c2)**c3")
1055 new_reaction%rate_type = rate_analytic_k7
1056 new_reaction%n_coeff = 3
1057 read(data_value(n), *) new_reaction%rate_data(1:3)
1058 case ("c1*(300/Tg)**c2")
1059 new_reaction%rate_type = rate_analytic_k8
1060 new_reaction%n_coeff = 2
1061 read(data_value(n), *) new_reaction%rate_data(1:2)
1062 case ("c1*exp(-c2*Tg)")
1063 new_reaction%rate_type = rate_analytic_k9
1064 new_reaction%n_coeff = 2
1065 read(data_value(n), *) new_reaction%rate_data(1:2)
1066 case ("10**(c1+c2*(Tg-300))")
1067 new_reaction%rate_type = rate_analytic_k10
1068 new_reaction%n_coeff = 2
1069 read(data_value(n), *) new_reaction%rate_data(1:2)
1070 case ("c1*(300/Tg)**c2*exp(-c3/Tg)")
1071 new_reaction%rate_type = rate_analytic_k11
1072 new_reaction%n_coeff = 3
1073 read(data_value(n), *) new_reaction%rate_data(1:3)
1074 case ("c1*Tg**c2*exp(-c3/Tg)")
1075 new_reaction%rate_type = rate_analytic_k12
1076 new_reaction%n_coeff = 3
1077 read(data_value(n), *) new_reaction%rate_data(1:3)
1078 case ("c1*exp(-(c2/(c3+Td))**c4)")
1079 new_reaction%rate_type = rate_analytic_k13
1080 new_reaction%n_coeff = 4
1081 read(data_value(n), *) new_reaction%rate_data(1:4)
1082 case ("c1*exp(-(Td/c2)**c3)")
1083 new_reaction%rate_type = rate_analytic_k14
1084 new_reaction%n_coeff = 3
1085 read(data_value(n), *) new_reaction%rate_data(1:3)
1086 case ("c1*exp(-(c2/(kb*(Tg+Td/c3)))**c4)")
1087 new_reaction%rate_type = rate_analytic_k15
1088 new_reaction%n_coeff = 4
1089 read(data_value(n), *) new_reaction%rate_data(1:4)
1090 case default
1091 print *, "Unknown rate type: ", trim(how_to_get(n))
1092 print *, "For reaction: ", trim(reaction(n))
1093 print *, "In file: ", trim(filename)
1094 print *, "See documentation/chemistry.md"
1095 if (how_to_get(n) /= "field_table" .and. &
1096 index(how_to_get(n), "c1") == 0) then
1097 print *, "You probably use the old reaction format"
1098 print *, "Try to use tools/chemistry_update_reactions.sh"
1099 print *, "See also the chemistry documentation"
1100 end if
1101 error stop "Unknown chemical reaction"
1102 end select
1103
1104 ! Correct for length unit in the rate function (e.g. [k] = cm3/s)
1105 select case (length_unit(n))
1106 case ("m")
1107 continue
1108 case ("cm")
1109 ! 1 cm^3 = 1e-6 m^3
1110 new_reaction%rate_factor = new_reaction%rate_factor * &
1111 1e-6_dp**(new_reaction%n_species_in-1)
1112 case default
1113 print *, "For reaction: ", trim(reaction(n))
1114 print *, "Invalid length unit: ", trim(length_unit(n))
1115 error stop
1116 end select
1117
1118 reactions(n_reactions) = new_reaction
1119 end do
1120
1121 if (n_reactions > 0) read_success = .true.
1122 return
1123
1124 ! Error messages
1125999 error stop "read_reactions: no closing dashes for reaction list"
1126 end subroutine read_reactions
1127
1128 !> Read a reaction table
1129 subroutine read_reaction_table(filename, dataname, rdata)
1131 character(len=*), intent(in) :: filename
1132 character(len=*), intent(in) :: dataname
1133 type(reaction_t), intent(inout) :: rdata
1134
1135 rdata%rate_type = rate_tabulated_field
1136 call table_from_file(filename, dataname, rdata%x_data, rdata%y_data)
1137 end subroutine read_reaction_table
1138
1139 !> Parse a reaction and store it
1140 subroutine parse_reaction(reaction_text, reaction, ignored_species, &
1141 keep_reaction)
1142 use m_gas
1143 character(len=*), intent(in) :: reaction_text
1144 type(reaction_t), intent(out) :: reaction
1145 character(len=comp_len), intent(in) :: ignored_species(:)
1146 logical, intent(out) :: keep_reaction
1147 integer, parameter :: max_components = 100
1148 character(len=comp_len) :: component
1149 integer :: i, ix, n, n_found, multiplicity
1150 integer :: n_in, n_out
1151 integer :: i0(max_components)
1152 integer :: i1(max_components)
1153 integer :: ix_in(max_components)
1154 integer :: ix_out(max_components)
1155 integer :: multiplicity_out(max_components)
1156 logical :: left_side, is_gas_species
1157 real(dp) :: rfactor
1158
1159 call get_fields_string(reaction_text, " ", max_components, n_found, i0, i1)
1160
1161 left_side = .true.
1162 keep_reaction = .true.
1163 n_in = 0
1164 n_out = 0
1165 rfactor = 1.0_dp
1166 reaction%n_species_in = 0
1167
1168 do n = 1, n_found
1169 component = reaction_text(i0(n):i1(n))
1170
1171 if (component == "+") cycle
1172
1173 if (component == "->") then
1174 left_side = .false.
1175 cycle
1176 end if
1177
1178 ! Assume we have a multiplicity less than 10
1179 if (lge(component(1:1), '1') .and. lle(component(1:1), '9')) then
1180 read(component(1:1), *) multiplicity
1181 component = component(2:)
1182 else
1183 multiplicity = 1
1184 end if
1185
1186 if (left_side) then
1187 reaction%n_species_in = reaction%n_species_in + multiplicity
1188 end if
1189
1190 ! If the gas density is constant, remove gas species from the reaction
1191 if (gas_constant_density) then
1192 ix = gas_index(component)
1193 if (ix /= -1) then
1194 ! Multiply reaction by constant density
1195 if (left_side) rfactor = rfactor * gas_densities(ix)
1196 cycle
1197 end if
1198
1199 if (component == "M") then
1200 ! Assume this stands for 'any molecule'
1201 if (left_side) rfactor = rfactor * gas_number_density
1202 cycle
1203 end if
1204 end if
1205
1206 ! Handle ignored species
1207 if (findloc(ignored_species, component, dim=1) > 0) then
1208 is_gas_species = (gas_index(component) > 0 .or. component == "M")
1209
1210 if (left_side .and. .not. is_gas_species) then
1211 ! Ignore the whole reaction, since this species will not be
1212 ! produced (and will thus have zero density)
1213 keep_reaction = .false.
1214 return
1215 else
1216 ! Ignore the production of this species, but keep the reaction
1217 cycle
1218 end if
1219 end if
1220
1221 ix = species_index(component)
1222 if (ix == -1) then
1223 if (n_species >= max_num_species) &
1224 error stop "Too many species, increase max_num_species"
1225 n_species = n_species + 1
1226 ix = n_species
1227 species_list(ix) = trim(component)
1228 end if
1229
1230 if (left_side) then
1231 do i = 1, multiplicity
1232 n_in = n_in + 1
1233 ix_in(n_in) = ix
1234 end do
1235 else
1236 ! Check if species is already present in right-hand side
1237 do i = 1, n_out
1238 if (ix == ix_out(i)) exit
1239 end do
1240
1241 if (i <= n_out) then
1242 multiplicity_out(i) = multiplicity_out(i) + multiplicity
1243 else
1244 ! If not already present, add the species
1245 n_out = n_out + 1
1246 ix_out(n_out) = ix
1247 multiplicity_out(n_out) = multiplicity
1248 end if
1249 end if
1250 end do
1251
1252 reaction%ix_in = ix_in(1:n_in)
1253 reaction%ix_out = ix_out(1:n_out)
1254 reaction%multiplicity_out = multiplicity_out(1:n_out)
1255 reaction%rate_factor = rfactor
1256
1257 if (n_in == 0) then
1258 print *, "Error in the following reaction:"
1259 print *, trim(reaction_text)
1260 error stop "No input species"
1261 end if
1262 end subroutine parse_reaction
1263
1264 !> Find index of a species, return -1 if not found
1265 elemental integer function species_index(name)
1266 character(len=*), intent(in) :: name
1267 do species_index = 1, n_species
1268 if (species_list(species_index) == name) exit
1269 end do
1271 end function species_index
1272
1273 !> Routine to find the indices of entries in a string
1274 subroutine get_fields_string(line, delims, n_max, n_found, ixs_start, ixs_end)
1275 !> The line from which we want to read
1276 character(len=*), intent(in) :: line
1277 !> A string with delimiters. For example delims = " ,'"""//char(9)
1278 character(len=*), intent(in) :: delims
1279 !> Maximum number of entries to read in
1280 integer, intent(in) :: n_max
1281 !> Number of entries found
1282 integer, intent(inout) :: n_found
1283 !> On return, ix_start(i) holds the starting point of entry i
1284 integer, intent(inout) :: ixs_start(n_max)
1285 !> On return, ix_end(i) holds the end point of entry i
1286 integer, intent(inout) :: ixs_end(n_max)
1287
1288 integer :: ix, ix_prev
1289
1290 ix_prev = 0
1291 n_found = 0
1292
1293 do while (n_found < n_max)
1294
1295 ! Find the starting point of the next entry (a non-delimiter value)
1296 ix = verify(line(ix_prev+1:), delims)
1297 if (ix == 0) exit
1298
1299 n_found = n_found + 1
1300 ixs_start(n_found) = ix_prev + ix ! This is the absolute position in 'line'
1301
1302 ! Get the end point of the current entry (next delimiter index minus one)
1303 ix = scan(line(ixs_start(n_found)+1:), delims) - 1
1304
1305 if (ix == -1) then ! If there is no last delimiter,
1306 ixs_end(n_found) = len(line) ! the end of the line is the endpoint
1307 else
1308 ixs_end(n_found) = ixs_start(n_found) + ix
1309 end if
1310
1311 ix_prev = ixs_end(n_found) ! We continue to search from here
1312 end do
1313
1314 end subroutine get_fields_string
1315
1316 !> Replace text in string, inspired by
1317 !> http://fortranwiki.org/fortran/show/String_Functions
1318 pure function string_replace(string, text, replacement) result(outs)
1319 character(len=*), intent(in) :: string
1320 character(len=*), intent(in) :: text
1321 character(len=*), intent(in) :: replacement
1322 character(len=:), allocatable :: outs
1323 integer, parameter :: buffer_space = 256
1324 character(len=(len(string)+buffer_space)) :: buffer
1325 integer :: i, nt, nr, len_outs
1326
1327 buffer = string
1328 nt = len_trim(text)
1329 nr = len_trim(replacement)
1330
1331 do
1332 i = index(buffer, text(:nt))
1333 if (i == 0) exit
1334 buffer = buffer(:i-1) // replacement(:nr) // trim(buffer(i+nt:))
1335 end do
1336
1337 len_outs = len_trim(buffer)
1338 allocate(character(len=len_outs) :: outs)
1339 outs = buffer(1:len_outs)
1340 end function string_replace
1341
1342 !> An inefficient routine to replace *^+- characters in a string
1343 subroutine to_simple_ascii(text, simple, charge)
1344 character(len=*), intent(in) :: text
1345 character(len=*), intent(inout) :: simple
1346 integer, intent(out) :: charge
1347 integer :: n
1348 logical :: in_brackets = .false.
1349
1350 charge = 0
1351 simple = ""
1352
1353 do n = 1, len_trim(text)
1354 select case (text(n:n))
1355 case ('(')
1356 in_brackets = .true.
1357 simple = trim(simple) // "_"
1358 case (')')
1359 in_brackets = .false.
1360 case ('*')
1361 simple = trim(simple) // "_star"
1362 case ('+')
1363 if (.not. in_brackets) then
1364 charge = charge + 1
1365 end if
1366 simple = trim(simple) // "_plus"
1367 case ('-')
1368 if (.not. in_brackets) then
1369 charge = charge - 1
1370 end if
1371 simple = trim(simple) // "_min"
1372 case ('^')
1373 simple = trim(simple) // "_hat"
1374 case ("'")
1375 simple = trim(simple) // "p"
1376 case default
1377 simple = trim(simple) // text(n:n)
1378 end select
1379 end do
1380
1381 ! Handle some species separately
1382 if (simple == "e" .or. simple == "stoch_e") charge = -1
1383 end subroutine to_simple_ascii
1384
1385end module m_chemistry
Module for handling chemical reactions.
character(len=20), dimension(*), parameter, public reaction_names
integer, parameter rate_tabulated_energy
Reaction with a field-dependent reaction rate.
subroutine, public read_reactions(filename, read_success)
Read reactions from a file.
subroutine, public chemistry_write_summary(fname)
Write a summary of the reactions (TODO) and the ionization and attachment coefficients (if working at...
integer, dimension(max_num_species), public, protected species_charge
Charge of the species.
integer, parameter, public attachment_reaction
Identifier for attachment reactions.
integer, parameter, public ionization_reaction
Identifier for ionization reactions.
subroutine, public get_rates(fields, rates, n_cells, energy_ev)
Compute reaction rates.
integer, parameter, public detachment_reaction
Identifier for detachment reactions.
integer, parameter, public general_reaction
Identifier for general reactions (not of any particular type)
subroutine, public chemistry_initialize(tree, cfg)
Initialize module and load chemical reactions.
integer, dimension(:), allocatable, public, protected charged_species_itree
List with indices of charged species.
integer, parameter, public recombination_reaction
Identifier for recombination reactions.
integer, public, protected n_stochastic_species
Number of stochastic species.
subroutine, public get_derivatives(dens, rates, derivs, n_cells)
Compute derivatives due to chemical reactions. Note that the 'rates' argument is modified.
integer, dimension(max_num_species), public, protected species_itree
species_itree(n) holds the index of species n in the tree (cell-centered variables)
subroutine, public chemistry_convert_stochastic_species(tree, rng)
Convert stochastic species into regular species, using random numbers.
subroutine, public chemistry_get_breakdown_field(field_td, min_growth_rate)
Get the breakdown field in Townsend.
character(len=comp_len), dimension(max_num_species), public, protected species_list
List of the species.
integer, public, protected n_species
Number of species present.
integer, public, protected n_plasma_species
Number of plasma species present.
elemental integer function, public species_index(name)
Find index of a species, return -1 if not found.
integer, dimension(:), allocatable, public, protected charged_species_charge
List with charges of charged species.
type(reaction_t), dimension(max_num_reactions), public, protected reactions
List of reactions.
integer, public, protected n_gas_species
Number of gas species present.
integer, public, protected n_reactions
Number of reactions present.
Module to set the time step.
Definition m_dt.f90:2
integer, public, protected time_integrator
Which time integrator is used.
Definition m_dt.f90:49
Module that stores parameters related to the gas.
Definition m_gas.f90:2
character(len=comp_len), dimension(:), allocatable, public, protected gas_components
Definition m_gas.f90:30
elemental integer function, public gas_index(name)
Find index of a gas component, return -1 if not found.
Definition m_gas.f90:285
real(dp), public, protected gas_temperature
Definition m_gas.f90:21
real(dp), parameter, public townsend_to_si
Definition m_gas.f90:42
real(dp), dimension(:), allocatable, public, protected gas_densities
Definition m_gas.f90:27
logical, public, protected gas_constant_density
Whether the gas has a constant density.
Definition m_gas.f90:12
real(dp), public, protected gas_number_density
Definition m_gas.f90:33
Module to set the type of model.
Definition m_model.f90:2
logical, public, protected model_has_energy_equation
Whether the model has an energy equation.
Definition m_model.f90:20
Module with settings and routines for tabulated data.
subroutine, public table_set_column(tbl, i_col, x, y, max_err)
Interpolate data and store in lookup table.
integer, public, protected table_xspacing
X-spacing for lookup table.
integer, public, protected table_size
How large lookup tables should be.
Module that provides routines for reading in arbritrary transport data.
real(dp), public, protected td_max_ev
Maximal energy (eV) in input data (automatically updated)
type(lt_t), public, protected td_tbl
integer, parameter, public td_diffusion
Electron diffusion constant.
integer, parameter, public td_eta
Attachment coefficient.
integer, public, protected td_energy_ev
Electron energy in eV (used with chemistry)
integer, parameter, public td_alpha
Ionization coefficient.
integer, parameter, public td_mobility
Electron mobility.
Module with basic types.
Definition m_types.f90:2
Module that contains physical and numerical constants.
real(dp), parameter uc_boltzmann_const
real(dp), parameter uc_elec_volt