3#include "../afivo/src/cpp_macros.h"
27 [character(len=20) ::
"ionization",
"attachment",
"recombination", &
28 "detachment",
"general"]
31 integer,
parameter :: rate_max_num_coeff = 9
35 integer,
allocatable :: ix_in(:)
36 integer,
allocatable :: ix_out(:)
37 integer,
allocatable :: multiplicity_out(:)
38 integer :: n_species_in
42 real(dp) :: rate_factor
45 real(dp) :: rate_data(rate_max_num_coeff) = -huge(1.0_dp)
46 integer :: lookup_table_index
47 real(dp),
allocatable :: x_data(:)
48 real(dp),
allocatable :: y_data(:)
49 character(len=50) :: description
54 integer,
allocatable :: ix_in(:)
55 integer,
allocatable :: ix_out(:)
56 integer,
allocatable :: multiplicity_out(:)
60 type stochastic_species_t
63 end type stochastic_species_t
69 integer,
parameter :: rate_tabulated_field = 1
72 integer,
parameter :: rate_analytic_constant = 2
75 integer,
parameter :: rate_analytic_linear = 3
78 integer,
parameter :: rate_analytic_exp_v1 = 4
81 integer,
parameter :: rate_analytic_exp_v2 = 5
84 integer,
parameter :: rate_analytic_k1 = 6
87 integer,
parameter :: rate_analytic_k3 = 8
90 integer,
parameter :: rate_analytic_k4 = 9
93 integer,
parameter :: rate_analytic_k5 = 10
96 integer,
parameter :: rate_analytic_k6 = 11
99 integer,
parameter :: rate_analytic_k7 = 12
102 integer,
parameter :: rate_analytic_k8 = 13
105 integer,
parameter :: rate_analytic_k9 = 14
108 integer,
parameter :: rate_analytic_k10 = 15
111 integer,
parameter :: rate_analytic_k11 = 16
114 integer,
parameter :: rate_analytic_k12 = 17
117 integer,
parameter :: rate_analytic_k13 = 18
120 integer,
parameter :: rate_analytic_k14 = 19
123 integer,
parameter :: rate_analytic_k15 = 20
126 integer,
parameter :: max_num_species = 100
129 integer,
parameter :: max_num_reactions = 500
147 character(len=comp_len),
public,
protected ::
species_list(max_num_species)
150 type(stochastic_species_t) :: stochastic_species_list(max_num_species)
159 type(reaction_t),
public,
protected ::
reactions(max_num_reactions)
162 type(tiny_react_t) :: tiny_react(max_num_reactions)
165 type(lt_t) :: chemtbl_fld
168 type(lt_t) :: chemtbl_ee
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
206 call cfg_get(cfg,
"input_data%file", reaction_file)
217 if (.not. read_success)
then
218 print *,
"m_chemistry: no reaction table found, using standard model"
232 reactions(1)%rate_type = rate_tabulated_field
238 reactions(1)%description =
"e + M -> e + e + M+"
245 reactions(2)%rate_type = rate_tabulated_field
253 error stop
"Varying gas density not yet supported"
273 tiny_react(n)%ix_out =
reactions(n)%ix_out
274 tiny_react(n)%multiplicity_out =
reactions(n)%multiplicity_out
301 call check_charge_conservation()
307 if (any(
reactions(n)%ix_in == i_elec) .and. &
308 .not. any(
reactions(n)%ix_out == i_elec) .and. &
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
321 else if (all(
reactions(n)%ix_in /= i_elec) .and. &
334 if (
reactions(n)%rate_type == rate_tabulated_field)
then
347 chemtbl_ee = lt_create(0.0_dp,
td_max_ev, &
353 if (
reactions(n)%rate_type == rate_tabulated_field)
then
379 call store_stochastic_species()
381 print *,
"--- List of reactions ---"
383 write(*,
"(I4,' (',I0,') ',A15,A)") n,
reactions(n)%n_species_in, &
387 print *,
"-------------------------"
389 print *,
"--- List of gas species ---"
393 print *,
"-------------------------"
395 print *,
"--- List of plasma species ---"
399 print *,
"-------------------------"
402 call chemistry_modify_rates(cfg)
407 subroutine store_stochastic_species()
408 integer :: n, i_stoch, namelen, iv
415 i_stoch = i_stoch + 1
422 print *,
"which should be present for stochastic species (stoch_)"
423 error stop
"Species for stochastic reaction not found"
432 end subroutine store_stochastic_species
436 type(af_t),
intent(inout) :: tree
437 type(rng_t),
intent(inout) :: rng
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
442 real(dp),
parameter :: two_pi = 2 * acos(-1.0_dp)
447 n_procs = omp_get_max_threads()
448 call prng%init_parallel(n_procs, rng)
452 proc_id = 1+omp_get_thread_num()
454 do lvl = 1, tree%highest_lvl
456 do iblock = 1,
size(tree%lvls(lvl)%leaves)
457 id = tree%lvls(lvl)%leaves(iblock)
458 product_dr = product(tree%boxes(id)%dr)
461 i_to = stochastic_species_list(n)%i_to
462 i_from = stochastic_species_list(n)%i_from
464 do kji_do(1, tree%n_cell)
466 if (tree%coord_t == af_cyl)
then
467 dv = two_pi * product_dr * af_cyl_radius_cc(tree%boxes(id), i)
475 lambda = dv * tree%boxes(id)%cc(ijk, i_from)
476 n_produced = prng%rngs(proc_id)%poisson(lambda)
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) + &
489 call prng%update_seed(rng)
494 subroutine chemistry_modify_rates(cfg)
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(:)
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)
509 if (n /= n_modified) &
510 error stop
"size(modified_reaction_ix) /= size(modified_rate_factors)"
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)
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"
528 end subroutine chemistry_modify_rates
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
545 n_fields =
td_tbl%n_points
547 if (n_fields < 3) error stop
"Not enough data for linear extrapolation"
549 allocate(fields(n_fields))
550 allocate(energies(n_fields))
552 allocate(eta(n_fields), alpha(n_fields), src(n_fields), loss(n_fields))
556 if (model_has_energy_equation)
then
558 call get_rates(fields, rates, n_fields, energies)
568 loss(:) = loss(:) + rates(:, n)
570 src(:) = src(:) + rates(:, n)
574 allocate(diff(n_fields))
577 allocate(mu(n_fields))
578 allocate(v(n_fields))
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)
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]"
596 write(my_unit, *) fields(n), &
607 subroutine check_charge_conservation()
608 integer :: n, q_in, q_out
614 if (q_in /= q_out)
then
616 error stop
"Charge is not conserved in the above reaction"
619 end subroutine check_charge_conservation
625 real(dp),
intent(out) :: field_td
627 real(dp),
intent(in) :: min_growth_rate
629 integer :: n, n_fields
630 real(dp),
allocatable :: energies(:)
631 real(dp),
allocatable :: fields(:), rates(:, :), src(:), loss(:)
633 n_fields =
td_tbl%n_points
634 allocate(fields(n_fields))
635 allocate(energies(n_fields))
637 allocate(src(n_fields), loss(n_fields))
640 if (model_has_energy_equation)
then
642 call get_rates(fields, rates, n_fields, energies)
652 loss(:) = loss(:) + rates(:, n)
654 src(:) = src(:) + rates(:, n)
658 do n = n_fields, 1, -1
659 if (src(n) - loss(n) < min_growth_rate)
exit
663 if (n > 0) field_td = fields(n)
673 integer,
intent(in) :: n_cells
674 real(dp),
intent(in) :: fields(n_cells)
675 real(dp),
intent(out) :: rates(n_cells,
n_reactions)
676 real(dp),
intent(in),
optional :: energy_ev(n_cells)
677 integer :: n, n_coeff
678 real(dp) :: c0, c(rate_max_num_coeff)
679 real(dp) :: te(n_cells)
680 logical :: te_available
683 real(dp),
parameter :: electron_ev_to_k = 2 *
uc_elec_volt / &
686 te_available = .false.
695 c(1:n_coeff) =
reactions(n)%rate_data(1:n_coeff)
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, &
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
718 te_available = .true.
720 rates(:, n) = c0 * c(1) * (300 / te)**c(2)
721 case (rate_analytic_k3)
722 if (.not. te_available)
then
724 te_available = .true.
728 case (rate_analytic_k4)
730 case (rate_analytic_k5)
732 case (rate_analytic_k6)
734 case (rate_analytic_k7)
736 case (rate_analytic_k8)
738 case (rate_analytic_k9)
740 case (rate_analytic_k10)
742 case (rate_analytic_k11)
744 case (rate_analytic_k12)
746 case (rate_analytic_k13)
747 rates(:, n) = c0 * c(1) * exp(-(c(2) / (c(3) + fields))**c(4))
748 case (rate_analytic_k14)
749 rates(:, n) = c0 * c(1) * exp(-(fields / c(2))**c(3))
750 case (rate_analytic_k15)
764 integer,
intent(in) :: n_cells
766 real(dp),
intent(in) :: dens(n_cells,
n_species)
768 real(dp),
intent(inout) :: rates(n_cells,
n_reactions)
770 real(dp),
intent(out) :: derivs(n_cells,
n_species)
773 derivs(:, :) = 0.0_dp
778 rates(:, n) = rates(:, n) * &
779 product(dens(:, tiny_react(n)%ix_in), dim=2)
782 do i = 1,
size(tiny_react(n)%ix_in)
783 ix = tiny_react(n)%ix_in(i)
784 derivs(:, ix) = derivs(:, ix) - rates(:, n)
788 do i = 1,
size(tiny_react(n)%ix_out)
789 ix = tiny_react(n)%ix_out(i)
790 derivs(:, ix) = derivs(:, ix) + rates(:, n) * &
791 tiny_react(n)%multiplicity_out(i)
797 subroutine read_ignored_species(filename, ignored_species)
798 character(len=*),
intent(in) :: filename
799 character(len=comp_len),
allocatable,
intent(inout) :: &
801 character(len=comp_len) :: tmp(max_num_species)
802 character(len=string_len) :: line
803 integer :: my_unit, n_ignored
807 open(newunit=my_unit, file=filename, action=
"read")
811 read(my_unit,
"(A)",
end=998) line
814 if (line ==
"ignored_species")
then
816 read(my_unit,
"(A)") line
817 if (line(1:5) /=
"-----") &
818 error stop
"ignored_species not followed by -----"
825 read(my_unit,
"(A)",
end=999) line
829 if (line(1:1) ==
"#") cycle
832 if (line(1:5) ==
"-----")
exit
834 n_ignored = n_ignored + 1
835 read(line, *) tmp(n_ignored)
839 ignored_species = tmp(1:n_ignored)
843999 error stop
"read_ignored_species: no closing dashes"
844 end subroutine read_ignored_species
848 character(len=*),
intent(in) :: filename
849 logical,
intent(out) :: read_success
850 character(len=string_len) :: line
851 integer,
parameter :: field_len = 50
852 character(len=field_len) :: data_value(max_num_reactions)
853 character(len=field_len) :: reaction(max_num_reactions)
854 character(len=field_len) :: how_to_get(max_num_reactions)
855 character(len=10) :: length_unit(max_num_reactions)
856 type(reaction_t) :: new_reaction
858 integer :: n_reactions_found
859 integer,
parameter :: n_fields_max = 40
860 integer :: i0(n_fields_max), i1(n_fields_max)
861 integer :: n, i, k, n_found, lo, hi
862 logical :: keep_reaction
864 character(len=comp_len),
allocatable :: ignored_species(:)
867 character(len=name_len) :: name
868 character(len=field_len),
allocatable :: members(:)
871 integer,
parameter :: max_groups = 10
872 integer :: i_group, group_size
873 type(group) :: groups(max_groups)
875 call read_ignored_species(filename, ignored_species)
877 open(newunit=my_unit, file=filename, action=
"read")
879 n_reactions_found = 0
882 read_success = .false.
886 read(my_unit,
"(A)",
end=998) line
889 if (line ==
"reaction_list")
then
891 read(my_unit,
"(A)") line
892 if (line(1:5) /=
"-----") &
893 error stop
"read_reactions: reaction_list not followed by -----"
900 read(my_unit,
"(A)",
end=999) line
904 if (line(1:1) ==
"#") cycle
907 if (line(1:5) ==
"-----")
exit
910 if (line(1:1) ==
"@")
then
911 call get_fields_string(line,
"=,", n_fields_max, n_found, i0, i1)
912 i_group = i_group + 1
914 if (i_group > max_groups) error stop
"Too many groups"
916 if (i_group == 1)
then
917 group_size = n_found - 1
918 else if (n_found - 1 /= group_size)
then
920 error stop
"Groups for a reaction should have the same size"
923 groups(i_group)%name = line(i0(1):i1(1))
924 allocate(groups(i_group)%members(n_found - 1))
927 groups(i_group)%members(n-1) = adjustl(line(i0(n):i1(n)))
932 if (i_group > 0)
then
934 lo = n_reactions_found
935 hi = n_reactions_found+group_size-1
936 reaction(lo+1:hi) = reaction(lo)
937 how_to_get(lo+1:hi) = how_to_get(lo)
938 data_value(lo+1:hi) = data_value(lo)
939 length_unit(lo+1:hi) = length_unit(lo)
940 n_reactions_found = hi
945 reaction(n) = string_replace(reaction(n), &
946 groups(i)%name, groups(i)%members(k))
947 how_to_get(n) = string_replace(how_to_get(n), &
948 groups(i)%name, groups(i)%members(k))
949 data_value(n) = string_replace(data_value(n), &
950 groups(i)%name, groups(i)%members(k))
956 if (count(reaction(lo:hi) == reaction(n)) > 1 .and. &
957 count(data_value(lo:hi) == data_value(n)) > 1)
then
959 print *, trim(reaction(k)),
",", data_value(k)
961 error stop
"Groups lead to duplicate reactions"
969 call get_fields_string(line,
",", n_fields_max, n_found, i0, i1)
971 if (n_found < 3 .or. n_found > 4)
then
973 error stop
"Invalid chemistry syntax"
976 if (n_reactions_found >= max_num_reactions) &
977 error stop
"Too many reactions, increase max_num_reactions"
979 n_reactions_found = n_reactions_found + 1
980 reaction(n_reactions_found) = line(i0(1):i1(1))
981 how_to_get(n_reactions_found) = line(i0(2):i1(2))
982 data_value(n_reactions_found) = line(i0(3):i1(3))
985 if (n_found > 3)
then
986 length_unit(n_reactions_found) = line(i0(4):i1(4))
988 length_unit(n_reactions_found) =
"m"
1000 do n = 1, n_reactions_found
1001 call parse_reaction(trim(reaction(n)), new_reaction, &
1002 ignored_species, keep_reaction)
1004 if (keep_reaction)
then
1010 new_reaction%description = trim(reaction(n))
1014 select case (how_to_get(n))
1015 case (
"field_table")
1017 call read_reaction_table(filename, &
1018 trim(data_value(n)), new_reaction)
1019 new_reaction%n_coeff = 0
1021 new_reaction%rate_type = rate_analytic_constant
1022 new_reaction%n_coeff = 1
1023 read(data_value(n), *) new_reaction%rate_data(1)
1025 new_reaction%rate_type = rate_analytic_linear
1026 new_reaction%n_coeff = 2
1027 read(data_value(n), *) new_reaction%rate_data(1:2)
1028 case (
"c1*exp(-(c2/(c3+Td))**2)")
1029 new_reaction%rate_type = rate_analytic_exp_v1
1030 new_reaction%n_coeff = 3
1031 read(data_value(n), *) new_reaction%rate_data(1:3)
1032 case (
"c1*exp(-(Td/c2)**2)")
1033 new_reaction%rate_type = rate_analytic_exp_v2
1034 new_reaction%n_coeff = 2
1035 read(data_value(n), *) new_reaction%rate_data(1:2)
1036 case (
"c1*(300/Te)**c2")
1037 new_reaction%rate_type = rate_analytic_k1
1038 new_reaction%n_coeff = 2
1039 read(data_value(n), *) new_reaction%rate_data(1:2)
1040 case (
"(c1*(kB_eV*Te+c2)**2-c3)*c4")
1041 new_reaction%rate_type = rate_analytic_k3
1042 new_reaction%n_coeff = 4
1043 read(data_value(n), *) new_reaction%rate_data(1:4)
1044 case (
"c1*(Tg/300)**c2*exp(-c3/Tg)")
1045 new_reaction%rate_type = rate_analytic_k4
1046 new_reaction%n_coeff = 3
1047 read(data_value(n), *) new_reaction%rate_data(1:3)
1048 case (
"c1*exp(-c2/Tg)")
1049 new_reaction%rate_type = rate_analytic_k5
1050 new_reaction%n_coeff = 2
1051 read(data_value(n), *) new_reaction%rate_data(1:2)
1053 new_reaction%rate_type = rate_analytic_k6
1054 new_reaction%n_coeff = 2
1055 read(data_value(n), *) new_reaction%rate_data(1:2)
1056 case (
"c1*(Tg/c2)**c3")
1057 new_reaction%rate_type = rate_analytic_k7
1058 new_reaction%n_coeff = 3
1059 read(data_value(n), *) new_reaction%rate_data(1:3)
1060 case (
"c1*(300/Tg)**c2")
1061 new_reaction%rate_type = rate_analytic_k8
1062 new_reaction%n_coeff = 2
1063 read(data_value(n), *) new_reaction%rate_data(1:2)
1064 case (
"c1*exp(-c2*Tg)")
1065 new_reaction%rate_type = rate_analytic_k9
1066 new_reaction%n_coeff = 2
1067 read(data_value(n), *) new_reaction%rate_data(1:2)
1068 case (
"10**(c1+c2*(Tg-300))")
1069 new_reaction%rate_type = rate_analytic_k10
1070 new_reaction%n_coeff = 2
1071 read(data_value(n), *) new_reaction%rate_data(1:2)
1072 case (
"c1*(300/Tg)**c2*exp(-c3/Tg)")
1073 new_reaction%rate_type = rate_analytic_k11
1074 new_reaction%n_coeff = 3
1075 read(data_value(n), *) new_reaction%rate_data(1:3)
1076 case (
"c1*Tg**c2*exp(-c3/Tg)")
1077 new_reaction%rate_type = rate_analytic_k12
1078 new_reaction%n_coeff = 3
1079 read(data_value(n), *) new_reaction%rate_data(1:3)
1080 case (
"c1*exp(-(c2/(c3+Td))**c4)")
1081 new_reaction%rate_type = rate_analytic_k13
1082 new_reaction%n_coeff = 4
1083 read(data_value(n), *) new_reaction%rate_data(1:4)
1084 case (
"c1*exp(-(Td/c2)**c3)")
1085 new_reaction%rate_type = rate_analytic_k14
1086 new_reaction%n_coeff = 3
1087 read(data_value(n), *) new_reaction%rate_data(1:3)
1088 case (
"c1*exp(-(c2/(kb*(Tg+Td/c3)))**c4)")
1089 new_reaction%rate_type = rate_analytic_k15
1090 new_reaction%n_coeff = 4
1091 read(data_value(n), *) new_reaction%rate_data(1:4)
1093 print *,
"Unknown rate type: ", trim(how_to_get(n))
1094 print *,
"For reaction: ", trim(reaction(n))
1095 print *,
"In file: ", trim(filename)
1096 print *,
"See documentation/chemistry.md"
1097 if (how_to_get(n) /=
"field_table" .and. &
1098 index(how_to_get(n),
"c1") == 0)
then
1099 print *,
"You probably use the old reaction format"
1100 print *,
"Try to use tools/chemistry_update_reactions.sh"
1101 print *,
"See also the chemistry documentation"
1103 error stop
"Unknown chemical reaction"
1107 select case (length_unit(n))
1112 new_reaction%rate_factor = new_reaction%rate_factor * &
1113 1e-6_dp**(new_reaction%n_species_in-1)
1115 print *,
"For reaction: ", trim(reaction(n))
1116 print *,
"Invalid length unit: ", trim(length_unit(n))
1127999 error stop
"read_reactions: no closing dashes for reaction list"
1131 subroutine read_reaction_table(filename, dataname, rdata)
1133 character(len=*),
intent(in) :: filename
1134 character(len=*),
intent(in) :: dataname
1135 type(reaction_t),
intent(inout) :: rdata
1137 rdata%rate_type = rate_tabulated_field
1138 call table_from_file(filename, dataname, rdata%x_data, rdata%y_data)
1139 end subroutine read_reaction_table
1142 subroutine parse_reaction(reaction_text, reaction, ignored_species, &
1145 character(len=*),
intent(in) :: reaction_text
1146 type(reaction_t),
intent(out) :: reaction
1147 character(len=comp_len),
intent(in) :: ignored_species(:)
1148 logical,
intent(out) :: keep_reaction
1149 integer,
parameter :: max_components = 100
1150 character(len=comp_len) :: component
1151 integer :: i, ix, n, n_found, multiplicity
1152 integer :: n_in, n_out
1153 integer :: i0(max_components)
1154 integer :: i1(max_components)
1155 integer :: ix_in(max_components)
1156 integer :: ix_out(max_components)
1157 integer :: multiplicity_out(max_components)
1158 logical :: left_side, is_gas_species
1161 call get_fields_string(reaction_text,
" ", max_components, n_found, i0, i1)
1164 keep_reaction = .true.
1168 reaction%n_species_in = 0
1171 component = reaction_text(i0(n):i1(n))
1173 if (component ==
"+") cycle
1175 if (component ==
"->")
then
1181 if (lge(component(1:1),
'1') .and. lle(component(1:1),
'9'))
then
1182 read(component(1:1), *) multiplicity
1183 component = component(2:)
1189 reaction%n_species_in = reaction%n_species_in + multiplicity
1201 if (component ==
"M")
then
1209 if (findloc(ignored_species, component, dim=1) > 0)
then
1210 is_gas_species = (
gas_index(component) > 0 .or. component ==
"M")
1212 if (left_side .and. .not. is_gas_species)
then
1215 keep_reaction = .false.
1226 error stop
"Too many species, increase max_num_species"
1233 do i = 1, multiplicity
1240 if (ix == ix_out(i))
exit
1243 if (i <= n_out)
then
1244 multiplicity_out(i) = multiplicity_out(i) + multiplicity
1249 multiplicity_out(n_out) = multiplicity
1254 reaction%ix_in = ix_in(1:n_in)
1255 reaction%ix_out = ix_out(1:n_out)
1256 reaction%multiplicity_out = multiplicity_out(1:n_out)
1257 reaction%rate_factor = rfactor
1260 print *,
"Error in the following reaction:"
1261 print *, trim(reaction_text)
1262 error stop
"No input species"
1264 end subroutine parse_reaction
1268 character(len=*),
intent(in) :: name
1276 subroutine get_fields_string(line, delims, n_max, n_found, ixs_start, ixs_end)
1278 character(len=*),
intent(in) :: line
1280 character(len=*),
intent(in) :: delims
1282 integer,
intent(in) :: n_max
1284 integer,
intent(inout) :: n_found
1286 integer,
intent(inout) :: ixs_start(n_max)
1288 integer,
intent(inout) :: ixs_end(n_max)
1290 integer :: ix, ix_prev
1295 do while (n_found < n_max)
1298 ix = verify(line(ix_prev+1:), delims)
1301 n_found = n_found + 1
1302 ixs_start(n_found) = ix_prev + ix
1305 ix = scan(line(ixs_start(n_found)+1:), delims) - 1
1308 ixs_end(n_found) = len(line)
1310 ixs_end(n_found) = ixs_start(n_found) + ix
1313 ix_prev = ixs_end(n_found)
1316 end subroutine get_fields_string
1320 pure function string_replace(string, text, replacement)
result(outs)
1321 character(len=*),
intent(in) :: string
1322 character(len=*),
intent(in) :: text
1323 character(len=*),
intent(in) :: replacement
1324 character(len=:),
allocatable :: outs
1325 integer,
parameter :: buffer_space = 256
1326 character(len=(len(string)+buffer_space)) :: buffer
1327 integer :: i, nt, nr, len_outs
1331 nr = len_trim(replacement)
1334 i = index(buffer, text(:nt))
1336 buffer = buffer(:i-1) // replacement(:nr) // trim(buffer(i+nt:))
1339 len_outs = len_trim(buffer)
1340 allocate(
character(len=len_outs) :: outs)
1341 outs = buffer(1:len_outs)
1342 end function string_replace
1345 subroutine to_simple_ascii(text, simple, charge)
1346 character(len=*),
intent(in) :: text
1347 character(len=*),
intent(inout) :: simple
1348 integer,
intent(out) :: charge
1350 logical :: in_brackets = .false.
1355 do n = 1, len_trim(text)
1356 select case (text(n:n))
1358 in_brackets = .true.
1359 simple = trim(simple) //
"_"
1361 in_brackets = .false.
1363 simple = trim(simple) //
"_star"
1365 if (.not. in_brackets)
then
1368 simple = trim(simple) //
"_plus"
1370 if (.not. in_brackets)
then
1373 simple = trim(simple) //
"_min"
1375 simple = trim(simple) //
"_hat"
1377 simple = trim(simple) //
"p"
1379 simple = trim(simple) // text(n:n)
1384 if (simple ==
"e" .or. simple ==
"stoch_e") charge = -1
1385 end subroutine to_simple_ascii
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.
integer, public, protected time_integrator
Which time integrator is used.
Module that stores parameters related to the gas.
character(len=comp_len), dimension(:), allocatable, public, protected gas_components
elemental integer function, public gas_index(name)
Find index of a gas component, return -1 if not found.
real(dp), public, protected gas_temperature
real(dp), parameter, public townsend_to_si
real(dp), dimension(:), allocatable, public, protected gas_densities
logical, public, protected gas_constant_density
Whether the gas has a constant density.
real(dp), public, protected gas_number_density
Module to set the type of model.
logical, public, protected model_has_energy_equation
Whether the model has an energy equation.
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 that contains physical and numerical constants.
real(dp), parameter uc_boltzmann_const
real(dp), parameter uc_elec_volt