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