afivo-streamer 1.1
1D/2D/3D streamer simulations with AMR
Loading...
Searching...
No Matches
m_photoi.f90
Go to the documentation of this file.
1!> Top-module for photoionization, which can make use of different methods
2module m_photoi
3#include "../afivo/src/cpp_macros.h"
6 use m_af_all
7 use m_streamer
8 use m_types
10
11 implicit none
12 private
13
14 !> Whether photoionization is enabled
15 logical, protected, public :: photoi_enabled = .false.
16
17 !> Whether photoemission is enabled
18 logical, protected, public :: photoe_enabled = .false.
19
20 !> Which photoionization method to use (helmholtz, montecarlo)
21 character(len=20) :: photoi_method = 'helmholtz'
22
23 !> Which photoionization source to use (Zheleznyak, from_species)
24 character(len=20) :: photoi_source_type = 'Zheleznyak'
25
26 !> Which species is ionized by photoionization
27 character(len=20) :: photoi_species = 'O2_plus'
28
29 !> Index of species ionized by photoionization (in list of chemical species)
30 integer, public, protected :: photoi_species_index = -1
31
32 !> Name of the excited species in case photoi_source_type is 'from_species'
33 character(len=string_len) :: photoi_excited_species = 'UNDEFINED'
34
35 !> Index of the excited species
36 integer :: i_excited_species = -1
37
38 !> Photoionization efficiency factor, typically around 0.05-0.1, not for
39 !> Helmholtz-Luque should be 1.0
40 real(dp) :: photoi_eta = 0.05_dp
41
42 !> Quenching pressure
43 real(dp) :: photoi_quenching_pressure = 40e-3 ! bar
44
45 !> Decay time in case photoi_source_type is 'from_species'
46 real(dp) :: photoi_photoemission_time = 0.0_dp
47
48 !> Update photoionization every N time step
49 integer, protected, public :: photoi_per_steps = 5
50
51 !> Update photoemission every N time step
52 integer, protected, public :: photoe_per_steps = 10
53
54 !> Optional variable (when using photoionization)
55 integer, public, protected :: i_photo = -1 ! Photoionization rate
56
57 public :: photoi_initialize
58 public :: photoi_set_src
59 ! public :: photoe_set_src
60
61contains
62
63 !> Initialize photoionization parameters
64 subroutine photoi_initialize(tree, cfg)
65 use m_config
66 type(af_t), intent(inout) :: tree
67 type(cfg_t), intent(inout) :: cfg !< The configuration for the simulation
68
69 !> [photoi_general_parameters]
70 call cfg_add_get(cfg, "photoi%enabled", photoi_enabled, &
71 "Whether photoionization is enabled")
72 call cfg_add_get(cfg, "photoi%per_steps", photoi_per_steps, &
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")
86 !> [photoi_general_parameters]
87
88 call cfg_add_get(cfg, "photoi%photoemission_time", photoi_photoemission_time, &
89 "Photoemission time delay in case photoi_source_type is 'from_species'")
90 call cfg_add_get(cfg, "photoe%enabled", photoe_enabled, &
91 "Whether photoemission is enabled")
92 call cfg_add_get(cfg, "photoe%per_steps", photoe_per_steps, &
93 "Update photoemission every N time step")
94
95 if (photoi_enabled) then
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"
100
101 photoi_species_index = species_index(photoi_species)
102 if (photoi_species_index == -1) then
103 print *, "photoi%species not found (", trim(photoi_species), ")"
104 error stop "photoi%species not present"
105 end if
106
107 call af_add_cc_variable(tree, "photo", ix=i_photo)
108 call af_set_cc_methods(tree, i_photo, photoi_helmh_bc, af_gc_interp)
109
110 select case (photoi_source_type)
111 case ('Zheleznyak')
112 if (photoi_excited_species /= undefined_str) &
113 error stop "Cannot use photoi%excited_species with Zheleznyak"
114 case ('from_species')
115 if (photoi_excited_species == undefined_str) &
116 error stop "Please define photoi%excited_species"
117
118 i_excited_species = &
119 af_find_cc_variable(tree, trim(photoi_excited_species))
120 case default
121 error stop "Unknown value for photoi%source_type"
122 end select
123 end if
124
125 select case (photoi_method)
126 case ("helmholtz")
127 call photoi_helmh_initialize(tree, cfg, photoi_enabled, photoi_eta)
128 call phmc_initialize(cfg, .false.)
129 case ("montecarlo")
130 call photoi_helmh_initialize(tree, cfg, .false., photoi_eta)
132 case default
133 print *, "Unknown photoi_method: ", trim(photoi_method)
134 error stop
135 end select
136
137 end subroutine photoi_initialize
138
139 !> Sets the photoionization
140 subroutine photoi_set_src(tree, dt)
142 use m_gas
143
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
147
148 ! Compute quench factor, because some excited species will be quenched by
149 ! collisions, preventing the emission of a UV photon
150 quench_fac = photoi_quenching_pressure / &
151 (gas_pressure + photoi_quenching_pressure)
152
153 ! Set photon production rate per cell, which is proportional to the
154 ! ionization rate.
155 select case (photoi_source_type)
156 case ('Zheleznyak')
157 call af_loop_box_arg(tree, photoionization_rate_from_alpha, &
158 [photoi_eta * quench_fac], .true.)
159 case ('from_species')
160 ! effective decay time when we have quenching
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
165 else
166 decay_rate = 1.0 / eff_decay_time
167 end if
168
169 call af_loop_box_arg(tree, photoionization_rate_from_species, &
170 [quench_fac * decay_rate, 1-decay_fraction], .true.)
171 end select
172
173 select case (photoi_method)
174 case ("helmholtz")
175 ! Use Helmholtz approximation
176 call photoi_helmh_compute(tree, i_photo)
177 case ("montecarlo")
178 if (phmc_physical_photons) then
179 call phmc_set_src(tree, st_rng, i_rhs, &
181 else
182 call phmc_set_src(tree, st_rng, i_rhs, &
184 end if
185 end select
186
187 end subroutine photoi_set_src
188
189 ! !> Sets the photoemission
190 ! subroutine photoe_set_src(tree, dt)
191 ! use m_units_constants
192 ! use m_gas
193
194 ! type(af_t), intent(inout) :: tree
195 ! real(dp), intent(in), optional :: dt
196 ! real(dp), parameter :: p_quench = 40.0e-3_dp
197 ! real(dp) :: quench_fac, decay_fraction, decay_rate
198
199 ! ! Compute quench factor, because some excited species will be quenched by
200 ! ! collisions, preventing the emission of a UV photon
201 ! quench_fac = p_quench / (gas_pressure + p_quench)
202
203 ! ! Set photon production rate per cell, which is proportional to the
204 ! ! ionization rate.
205 ! select case (photoi_source_type)
206 ! case ('Zheleznyak')
207 ! call af_loop_box_arg(tree, photoionization_rate_from_alpha, &
208 ! [photoi_eta * quench_fac], .true.)
209 ! case ('from_species')
210 ! decay_fraction = 1 - exp(-dt / eff_decay_time)
211
212 ! if (dt > 1e-6_dp * photoi_decay_time) then
213 ! decay_rate = decay_fraction / dt
214 ! else
215 ! decay_rate = 1 / photoi_decay_time
216 ! end if
217
218 ! call af_loop_box_arg(tree, photoionization_rate_from_species, &
219 ! [photoi_eta * quench_fac * decay_rate, 1-decay_fraction], .true.)
220 ! end select
221
222 ! if (phmc_physical_photons) then
223 ! call phe_mc_set_src(tree, ST_rng, i_rhs, &
224 ! i_photo, ST_cylindrical, dt)
225 ! else
226 ! call phe_mc_set_src(tree, ST_rng, i_rhs, &
227 ! i_photo, ST_cylindrical)
228 ! end if
229
230 ! end subroutine photoe_set_src
231
232 !> Sets the photoionization_rate
233 subroutine photoionization_rate_from_alpha(box, coeff)
234 use m_streamer
235 use m_lookup_table
237 use m_gas
238 type(box_t), intent(inout) :: box
239 real(dp), intent(in) :: coeff(:)
240 integer :: ijk, nc
241 real(dp) :: fld, td, alpha, mobility, tmp, gas_dens
242 type(lt_loc_t) :: loc
243
244 nc = box%n_cell
245
246 do kji_do(1,nc)
247 fld = box%cc(ijk, i_electric_fld)
248
249 if (gas_constant_density) then
250 gas_dens = gas_number_density
251 else
252 gas_dens = box%cc(ijk, i_gas_dens)
253 end if
254
255 td = fld * si_to_townsend / gas_dens
256 loc = lt_get_loc(td_tbl, td)
257 ! No normalization required: (alpha/N) * (mu/N)
258 alpha = lt_get_col_at_loc(td_tbl, td_alpha, loc)
259 mobility = lt_get_col_at_loc(td_tbl, td_mobility, loc)
260
261 tmp = fld * mobility * alpha * box%cc(ijk, i_electron) * coeff(1)
262 if (tmp < 0) tmp = 0
263 box%cc(ijk, i_rhs) = tmp
264 end do; close_do
265 end subroutine photoionization_rate_from_alpha
266
267 !> Sets the photoionization_rate
268 subroutine photoionization_rate_from_species(box, coeff)
269 use m_streamer
270 use m_lookup_table
272 use m_gas
273 type(box_t), intent(inout) :: box
274 real(dp), intent(in) :: coeff(:)
275 integer :: nc
276 real(dp) :: source_factor, decay_factor
277
278 nc = box%n_cell
279 source_factor = coeff(1)
280 decay_factor = coeff(2)
281
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
287
288end module m_photoi
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.
Definition m_gas.f90:2
integer, public, protected i_gas_dens
Definition m_gas.f90:45
real(dp), parameter, public si_to_townsend
Definition m_gas.f90:39
real(dp), public, protected gas_pressure
Definition m_gas.f90:18
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 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.
Definition m_photoi.f90:2
integer, public, protected photoi_species_index
Index of species ionized by photoionization (in list of chemical species)
Definition m_photoi.f90:30
subroutine, public photoi_set_src(tree, dt)
Sets the photoionization.
Definition m_photoi.f90:141
integer, public, protected photoe_per_steps
Update photoemission every N time step.
Definition m_photoi.f90:52
logical, public, protected photoe_enabled
Whether photoemission is enabled.
Definition m_photoi.f90:18
logical, public, protected photoi_enabled
Whether photoionization is enabled.
Definition m_photoi.f90:15
integer, public, protected i_photo
Optional variable (when using photoionization)
Definition m_photoi.f90:55
subroutine, public photoi_initialize(tree, cfg)
Initialize photoionization parameters.
Definition m_photoi.f90:65
integer, public, protected photoi_per_steps
Update photoionization every N time step.
Definition m_photoi.f90:49
This module contains several pre-defined variables like:
Definition m_streamer.f90:6
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.
Module with basic types.
Definition m_types.f90:2
character(len= *), parameter undefined_str
Undefined string.
Definition m_types.f90:10
Module that contains physical and numerical constants.