3#include "../afivo/src/cpp_macros.h"
21 character(len=20) :: photoi_method =
'helmholtz'
24 character(len=20) :: photoi_source_type =
'Zheleznyak'
27 character(len=20) :: photoi_species =
'O2_plus'
33 character(len=string_len) :: photoi_excited_species =
'UNDEFINED'
36 integer :: i_excited_species = -1
40 real(dp) :: photoi_eta = 0.05_dp
43 real(dp) :: photoi_quenching_pressure = 40e-3
46 real(dp) :: photoi_photoemission_time = 0.0_dp
55 integer,
public,
protected ::
i_photo = -1
66 type(af_t),
intent(inout) :: tree
67 type(cfg_t),
intent(inout) :: cfg
71 "Whether photoionization is enabled")
73 "Update photoionization every N time steps")
74 call cfg_add_get(cfg,
"photoi%method", photoi_method, &
75 "Which photoionization method to use (helmholtz, montecarlo)")
76 call cfg_add_get(cfg,
"photoi%eta", photoi_eta, &
77 "Photoionization efficiency factor, typically around 0.05-0.1")
78 call cfg_add_get(cfg,
"photoi%quenching_pressure", photoi_quenching_pressure, &
79 "Photoionization quenching pressure (bar)")
80 call cfg_add_get(cfg,
"photoi%source_type", photoi_source_type, &
81 "How to compute the photoi. source (Zheleznyak, from_species)")
82 call cfg_add_get(cfg,
"photoi%excited_species", photoi_excited_species, &
83 "Which excited species to use when photoi%source_type = from_species")
84 call cfg_add_get(cfg,
"photoi%species", photoi_species, &
85 "Which species is ionized by photoionization")
88 call cfg_add_get(cfg,
"photoi%photoemission_time", photoi_photoemission_time, &
89 "Photoemission time delay in case photoi_source_type is 'from_species'")
91 "Whether photoemission is enabled")
93 "Update photoemission every N time step")
96 if (photoi_eta <= 0.0_dp) error stop
"photoi%eta <= 0.0"
97 if (photoi_eta > 1.0_dp) error stop
"photoi%eta > 1.0"
98 if (photoi_quenching_pressure <= 0.0_dp) &
99 error stop
"photoi%quenching_pressure <= 0.0"
103 print *,
"photoi%species not found (", trim(photoi_species),
")"
104 error stop
"photoi%species not present"
107 call af_add_cc_variable(tree,
"photo", ix=
i_photo)
110 select case (photoi_source_type)
113 error stop
"Cannot use photoi%excited_species with Zheleznyak"
114 case (
'from_species')
116 error stop
"Please define photoi%excited_species"
118 i_excited_species = &
119 af_find_cc_variable(tree, trim(photoi_excited_species))
121 error stop
"Unknown value for photoi%source_type"
125 select case (photoi_method)
133 print *,
"Unknown photoi_method: ", trim(photoi_method)
144 type(af_t),
intent(inout) :: tree
145 real(dp),
intent(in),
optional :: dt
146 real(dp) :: quench_fac, decay_fraction, eff_decay_time, decay_rate
150 quench_fac = photoi_quenching_pressure / &
155 select case (photoi_source_type)
157 call af_loop_box_arg(tree, photoionization_rate_from_alpha, &
158 [photoi_eta * quench_fac], .true.)
159 case (
'from_species')
161 eff_decay_time = photoi_photoemission_time
162 decay_fraction = 1 - exp(-dt / eff_decay_time)
163 if (dt > 1e-6_dp * eff_decay_time)
then
164 decay_rate = decay_fraction / dt
166 decay_rate = 1.0 / eff_decay_time
169 call af_loop_box_arg(tree, photoionization_rate_from_species, &
170 [quench_fac * decay_rate, 1-decay_fraction], .true.)
173 select case (photoi_method)
233 subroutine photoionization_rate_from_alpha(box, coeff)
238 type(box_t),
intent(inout) :: box
239 real(dp),
intent(in) :: coeff(:)
241 real(dp) :: fld, td, alpha, mobility, tmp, gas_dens
242 type(lt_loc_t) :: loc
256 loc = lt_get_loc(
td_tbl, td)
261 tmp = fld * mobility * alpha * box%cc(ijk,
i_electron) * coeff(1)
263 box%cc(ijk,
i_rhs) = tmp
265 end subroutine photoionization_rate_from_alpha
268 subroutine photoionization_rate_from_species(box, coeff)
273 type(box_t),
intent(inout) :: box
274 real(dp),
intent(in) :: coeff(:)
276 real(dp) :: source_factor, decay_factor
279 source_factor = coeff(1)
280 decay_factor = coeff(2)
282 box%cc(dtimes(1:nc),
i_rhs) = source_factor * &
283 box%cc(dtimes(1:nc), i_excited_species)
284 box%cc(dtimes(1:nc), i_excited_species) = decay_factor * &
285 box%cc(dtimes(1:nc), i_excited_species)
286 end subroutine photoionization_rate_from_species
Module for handling chemical reactions.
elemental integer function, public species_index(name)
Find index of a species, return -1 if not found.
Module that stores parameters related to the gas.
integer, public, protected i_gas_dens
real(dp), parameter, public si_to_townsend
real(dp), public, protected gas_pressure
logical, public, protected gas_constant_density
Whether the gas has a constant density.
real(dp), public, protected gas_number_density
Module for photoionization with the Helmholtz approximation.
subroutine, public photoi_helmh_compute(tree, i_photo)
subroutine, public photoi_helmh_bc(box, nb, iv, coords, bc_val, bc_type)
subroutine, public photoi_helmh_initialize(tree, cfg, is_used, eta)
Initialize options for Helmholtz photoionization.
Module that provides routines for Monte Carlo photoionization.
subroutine, public phmc_set_src(tree, rng, i_src, i_photo, use_cyl, dt)
Set the source term due to photoionization. At most phmc_num_photons discrete photons are produced.
subroutine, public phmc_initialize(cfg, is_used)
Initialize photoionization parameters.
logical, public, protected phmc_physical_photons
Whether physical photons are used.
Top-module for photoionization, which can make use of different methods.
integer, public, protected photoi_species_index
Index of species ionized by photoionization (in list of chemical species)
subroutine, public photoi_set_src(tree, dt)
Sets the photoionization.
integer, public, protected photoe_per_steps
Update photoemission every N time step.
logical, public, protected photoe_enabled
Whether photoemission is enabled.
logical, public, protected photoi_enabled
Whether photoionization is enabled.
integer, public, protected i_photo
Optional variable (when using photoionization)
subroutine, public photoi_initialize(tree, cfg)
Initialize photoionization parameters.
integer, public, protected photoi_per_steps
Update photoionization every N time step.
This module contains several pre-defined variables like:
logical, public, protected st_cylindrical
Whether cylindrical coordinates are used.
integer, public, protected i_electron
Index of electron density.
integer, public, protected i_rhs
Index of source term Poisson.
integer, public, protected i_electric_fld
Index of electric field norm.
type(rng_t), public st_rng
Random number generator.
Module that provides routines for reading in arbritrary transport data.
type(lt_t), public, protected td_tbl
integer, parameter, public td_alpha
Ionization coefficient.
integer, parameter, public td_mobility
Electron mobility.
character(len= *), parameter undefined_str
Undefined string.
Module that contains physical and numerical constants.