3#include "../afivo/src/cpp_macros.h"
25 [character(len=20) ::
"ionization",
"attachment",
"recombination", &
26 "detachment",
"general"]
29 integer,
parameter :: rate_max_num_coeff = 9
33 integer,
allocatable :: ix_in(:)
34 integer,
allocatable :: ix_out(:)
35 integer,
allocatable :: multiplicity_out(:)
36 integer :: n_species_in
40 real(dp) :: rate_factor
43 real(dp) :: rate_data(rate_max_num_coeff) = -huge(1.0_dp)
44 integer :: lookup_table_index
45 real(dp),
allocatable :: x_data(:)
46 real(dp),
allocatable :: y_data(:)
47 character(len=50) :: description
52 integer,
allocatable :: ix_in(:)
53 integer,
allocatable :: ix_out(:)
54 integer,
allocatable :: multiplicity_out(:)
61 integer,
parameter :: rate_tabulated_field = 1
64 integer,
parameter :: rate_analytic_constant = 2
67 integer,
parameter :: rate_analytic_linear = 3
70 integer,
parameter :: rate_analytic_exp_v1 = 4
73 integer,
parameter :: rate_analytic_exp_v2 = 5
76 integer,
parameter :: rate_analytic_k1 = 6
79 integer,
parameter :: rate_analytic_k3 = 8
82 integer,
parameter :: rate_analytic_k4 = 9
85 integer,
parameter :: rate_analytic_k5 = 10
88 integer,
parameter :: rate_analytic_k6 = 11
91 integer,
parameter :: rate_analytic_k7 = 12
94 integer,
parameter :: rate_analytic_k8 = 13
97 integer,
parameter :: rate_analytic_k9 = 14
100 integer,
parameter :: rate_analytic_k10 = 15
103 integer,
parameter :: rate_analytic_k11 = 16
106 integer,
parameter :: rate_analytic_k12 = 17
109 integer,
parameter :: rate_analytic_k13 = 18
112 integer,
parameter :: rate_analytic_k14 = 19
115 integer,
parameter :: rate_analytic_k15 = 20
118 integer,
parameter :: max_num_species = 100
121 integer,
parameter :: max_num_reactions = 500
136 character(len=comp_len),
public,
protected ::
species_list(max_num_species)
145 type(reaction_t),
public,
protected ::
reactions(max_num_reactions)
148 type(tiny_react_t) :: tiny_react(max_num_reactions)
151 type(lt_t) :: chemtbl_fld
154 type(lt_t) :: chemtbl_ee
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
191 call cfg_get(cfg,
"input_data%file", reaction_file)
202 if (.not. read_success)
then
203 print *,
"m_chemistry: no reaction table found, using standard model"
217 reactions(1)%rate_type = rate_tabulated_field
223 reactions(1)%description =
"e + M -> e + e + M+"
230 reactions(2)%rate_type = rate_tabulated_field
238 error stop
"Varying gas density not yet supported"
258 tiny_react(n)%ix_out =
reactions(n)%ix_out
259 tiny_react(n)%multiplicity_out =
reactions(n)%multiplicity_out
286 call check_charge_conservation()
292 if (any(
reactions(n)%ix_in == i_elec) .and. &
293 .not. any(
reactions(n)%ix_out == i_elec) .and. &
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
306 else if (all(
reactions(n)%ix_in /= i_elec) .and. &
319 if (
reactions(n)%rate_type == rate_tabulated_field)
then
332 chemtbl_ee = lt_create(0.0_dp,
td_max_ev, &
338 if (
reactions(n)%rate_type == rate_tabulated_field)
then
364 print *,
"--- List of reactions ---"
366 write(*,
"(I4,' (',I0,') ',A15,A)") n,
reactions(n)%n_species_in, &
370 print *,
"-------------------------"
372 print *,
"--- List of gas species ---"
376 print *,
"-------------------------"
378 print *,
"--- List of plasma species ---"
382 print *,
"-------------------------"
385 call chemistry_modify_rates(cfg)
390 subroutine chemistry_modify_rates(cfg)
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(:)
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)
405 if (n /= n_modified) &
406 error stop
"size(modified_reaction_ix) /= size(modified_rate_factors)"
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)
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"
424 end subroutine chemistry_modify_rates
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
441 n_fields =
td_tbl%n_points
443 if (n_fields < 3) error stop
"Not enough data for linear extrapolation"
445 allocate(fields(n_fields))
446 allocate(energies(n_fields))
448 allocate(eta(n_fields), alpha(n_fields), src(n_fields), loss(n_fields))
454 call get_rates(fields, rates, n_fields, energies)
464 loss(:) = loss(:) + rates(:, n)
466 src(:) = src(:) + rates(:, n)
470 allocate(diff(n_fields))
473 allocate(mu(n_fields))
474 allocate(v(n_fields))
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)
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]"
492 write(my_unit, *) fields(n), &
503 subroutine check_charge_conservation()
504 integer :: n, q_in, q_out
510 if (q_in /= q_out)
then
512 error stop
"Charge is not conserved in the above reaction"
515 end subroutine check_charge_conservation
521 real(dp),
intent(out) :: field_td
523 real(dp),
intent(in) :: min_growth_rate
525 integer :: n, n_fields
526 real(dp),
allocatable :: energies(:)
527 real(dp),
allocatable :: fields(:), rates(:, :), src(:), loss(:)
529 n_fields =
td_tbl%n_points
530 allocate(fields(n_fields))
531 allocate(energies(n_fields))
533 allocate(src(n_fields), loss(n_fields))
538 call get_rates(fields, rates, n_fields, energies)
548 loss(:) = loss(:) + rates(:, n)
550 src(:) = src(:) + rates(:, n)
554 do n = n_fields, 1, -1
555 if (src(n) - loss(n) < min_growth_rate)
exit
559 if (n > 0) field_td = fields(n)
569 integer,
intent(in) :: n_cells
570 real(dp),
intent(in) :: fields(n_cells)
571 real(dp),
intent(out) :: rates(n_cells,
n_reactions)
572 real(dp),
intent(in),
optional :: energy_ev(n_cells)
573 integer :: n, n_coeff
574 real(dp) :: c0, c(rate_max_num_coeff)
575 real(dp) :: te(n_cells)
576 logical :: te_available
579 real(dp),
parameter :: electron_ev_to_k = 2 *
uc_elec_volt / &
582 te_available = .false.
591 c(1:n_coeff) =
reactions(n)%rate_data(1:n_coeff)
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, &
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
615 rates(:, n) = c0 * c(1) * (300 / te)**c(2)
616 case (rate_analytic_k3)
617 if (.not. te_available)
then
622 case (rate_analytic_k4)
624 case (rate_analytic_k5)
626 case (rate_analytic_k6)
628 case (rate_analytic_k7)
630 case (rate_analytic_k8)
632 case (rate_analytic_k9)
634 case (rate_analytic_k10)
636 case (rate_analytic_k11)
638 case (rate_analytic_k12)
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)
658 integer,
intent(in) :: n_cells
660 real(dp),
intent(in) :: dens(n_cells,
n_species)
662 real(dp),
intent(inout) :: rates(n_cells,
n_reactions)
664 real(dp),
intent(out) :: derivs(n_cells,
n_species)
667 derivs(:, :) = 0.0_dp
672 rates(:, n) = rates(:, n) * &
673 product(dens(:, tiny_react(n)%ix_in), dim=2)
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)
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)
691 subroutine read_ignored_species(filename, ignored_species)
692 character(len=*),
intent(in) :: filename
693 character(len=comp_len),
allocatable,
intent(inout) :: &
695 character(len=comp_len) :: tmp(max_num_species)
696 character(len=string_len) :: line
697 integer :: my_unit, n_ignored
701 open(newunit=my_unit, file=filename, action=
"read")
705 read(my_unit,
"(A)",
end=998) line
708 if (line ==
"ignored_species")
then
710 read(my_unit,
"(A)") line
711 if (line(1:5) /=
"-----") &
712 error stop
"ignored_species not followed by -----"
719 read(my_unit,
"(A)",
end=999) line
723 if (line(1:1) ==
"#") cycle
726 if (line(1:5) ==
"-----")
exit
728 n_ignored = n_ignored + 1
729 read(line, *) tmp(n_ignored)
733 ignored_species = tmp(1:n_ignored)
737999 error stop
"read_ignored_species: no closing dashes"
738 end subroutine read_ignored_species
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
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
758 character(len=comp_len),
allocatable :: ignored_species(:)
761 character(len=name_len) :: name
762 character(len=field_len),
allocatable :: members(:)
765 integer,
parameter :: max_groups = 10
766 integer :: i_group, group_size
767 type(group) :: groups(max_groups)
769 call read_ignored_species(filename, ignored_species)
771 open(newunit=my_unit, file=filename, action=
"read")
773 n_reactions_found = 0
776 read_success = .false.
780 read(my_unit,
"(A)",
end=998) line
783 if (line ==
"reaction_list")
then
785 read(my_unit,
"(A)") line
786 if (line(1:5) /=
"-----") &
787 error stop
"read_reactions: reaction_list not followed by -----"
794 read(my_unit,
"(A)",
end=999) line
798 if (line(1:1) ==
"#") cycle
801 if (line(1:5) ==
"-----")
exit
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
808 if (i_group > max_groups) error stop
"Too many groups"
810 if (i_group == 1)
then
811 group_size = n_found - 1
812 else if (n_found - 1 /= group_size)
then
814 error stop
"Groups for a reaction should have the same size"
817 groups(i_group)%name = line(i0(1):i1(1))
818 allocate(groups(i_group)%members(n_found - 1))
821 groups(i_group)%members(n-1) = adjustl(line(i0(n):i1(n)))
826 if (i_group > 0)
then
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
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))
850 if (count(reaction(lo:hi) == reaction(n)) > 1 .and. &
851 count(data_value(lo:hi) == data_value(n)) > 1)
then
853 print *, trim(reaction(k)),
",", data_value(k)
855 error stop
"Groups lead to duplicate reactions"
863 call get_fields_string(line,
",", n_fields_max, n_found, i0, i1)
865 if (n_found < 3 .or. n_found > 4)
then
867 error stop
"Invalid chemistry syntax"
870 if (n_reactions_found >= max_num_reactions) &
871 error stop
"Too many reactions, increase max_num_reactions"
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))
879 if (n_found > 3)
then
880 length_unit(n_reactions_found) = line(i0(4):i1(4))
882 length_unit(n_reactions_found) =
"m"
894 do n = 1, n_reactions_found
895 call parse_reaction(trim(reaction(n)), new_reaction, &
896 ignored_species, keep_reaction)
898 if (keep_reaction)
then
904 new_reaction%description = trim(reaction(n))
908 select case (how_to_get(n))
911 call read_reaction_table(filename, &
912 trim(data_value(n)), new_reaction)
913 new_reaction%n_coeff = 0
915 new_reaction%rate_type = rate_analytic_constant
916 new_reaction%n_coeff = 1
917 read(data_value(n), *) new_reaction%rate_data(1)
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)
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)
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"
997 error stop
"Unknown chemical reaction"
1001 select case (length_unit(n))
1006 new_reaction%rate_factor = new_reaction%rate_factor * &
1007 1e-6_dp**(new_reaction%n_species_in-1)
1009 print *,
"For reaction: ", trim(reaction(n))
1010 print *,
"Invalid length unit: ", trim(length_unit(n))
1021999 error stop
"read_reactions: no closing dashes for reaction list"
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
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
1036 subroutine parse_reaction(reaction_text, reaction, ignored_species, &
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
1055 call get_fields_string(reaction_text,
" ", max_components, n_found, i0, i1)
1058 keep_reaction = .true.
1062 reaction%n_species_in = 0
1065 component = reaction_text(i0(n):i1(n))
1067 if (component ==
"+") cycle
1069 if (component ==
"->")
then
1075 if (lge(component(1:1),
'1') .and. lle(component(1:1),
'9'))
then
1076 read(component(1:1), *) multiplicity
1077 component = component(2:)
1083 reaction%n_species_in = reaction%n_species_in + multiplicity
1095 if (component ==
"M")
then
1103 if (findloc(ignored_species, component, dim=1) > 0)
then
1104 is_gas_species = (
gas_index(component) > 0 .or. component ==
"M")
1106 if (left_side .and. .not. is_gas_species)
then
1109 keep_reaction = .false.
1120 error stop
"Too many species, increase max_num_species"
1127 do i = 1, multiplicity
1134 if (ix == ix_out(i))
exit
1137 if (i <= n_out)
then
1138 multiplicity_out(i) = multiplicity_out(i) + multiplicity
1143 multiplicity_out(n_out) = multiplicity
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
1154 print *,
"Error in the following reaction:"
1155 print *, trim(reaction_text)
1156 error stop
"No input species"
1158 end subroutine parse_reaction
1162 character(len=*),
intent(in) :: name
1170 subroutine get_fields_string(line, delims, n_max, n_found, ixs_start, ixs_end)
1172 character(len=*),
intent(in) :: line
1174 character(len=*),
intent(in) :: delims
1176 integer,
intent(in) :: n_max
1178 integer,
intent(inout) :: n_found
1180 integer,
intent(inout) :: ixs_start(n_max)
1182 integer,
intent(inout) :: ixs_end(n_max)
1184 integer :: ix, ix_prev
1189 do while (n_found < n_max)
1192 ix = verify(line(ix_prev+1:), delims)
1195 n_found = n_found + 1
1196 ixs_start(n_found) = ix_prev + ix
1199 ix = scan(line(ixs_start(n_found)+1:), delims) - 1
1202 ixs_end(n_found) = len(line)
1204 ixs_end(n_found) = ixs_start(n_found) + ix
1207 ix_prev = ixs_end(n_found)
1210 end subroutine get_fields_string
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
1225 nr = len_trim(replacement)
1228 i = index(buffer, text(:nt))
1230 buffer = buffer(:i-1) // replacement(:nr) // trim(buffer(i+nt:))
1233 len_outs = len_trim(buffer)
1234 allocate(
character(len=len_outs) :: outs)
1235 outs = buffer(1:len_outs)
1236 end function string_replace
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
1244 logical :: in_brackets = .false.
1249 do n = 1, len_trim(text)
1250 select case (text(n:n))
1252 in_brackets = .true.
1253 simple = trim(simple) //
"_"
1255 in_brackets = .false.
1257 simple = trim(simple) //
"_star"
1259 if (.not. in_brackets)
then
1262 simple = trim(simple) //
"_plus"
1264 if (.not. in_brackets)
then
1267 simple = trim(simple) //
"_min"
1269 simple = trim(simple) //
"_hat"
1271 simple = trim(simple) //
"p"
1273 simple = trim(simple) // text(n:n)
1278 if (simple ==
"e") charge = -1
1279 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.
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.
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