afivo-streamer 1.1
1D/2D/3D streamer simulations with AMR
Loading...
Searching...
No Matches
m_photoi_mc.f90
Go to the documentation of this file.
1!> Module that provides routines for Monte Carlo photoionization
3#include "../afivo/src/cpp_macros.h"
4 use m_streamer
5 use m_lookup_table
6 use m_af_all
8
9 implicit none
10 private
11
12 !> Type to quickly look up absorption lengths from a table
13 type, public :: phmc_tbl_t
14 type(lt_t) :: tbl !< The lookup table
15 !> Value of full integral over the absorption function (can be less than
16 !> one if non-ionizing absorption is present)
17 real(dp) :: full_integral
18 end type phmc_tbl_t
19
20 !> Whether physical photons are used
21 logical, public, protected :: phmc_physical_photons = .true.
22
23 !> Whether a constant or adaptive grid spacing is used
24 logical, protected :: phmc_const_dx = .true.
25
26 !> Minimum grid spacing for photoionization
27 real(dp), protected :: phmc_min_dx = 1e-9_dp
28
29 !> At which grid spacing photons are absorbed compared to their mean distance
30 real(dp), protected :: phmc_absorp_fac = 0.25_dp
31
32 !> Number of photons to use
33 integer, protected :: phmc_num_photons = 5000*1000
34
35 !> Minimal photon weight
36 real(dp), protected :: phmc_min_weight = 1.0_dp
37
38 !> Non-ionizing absorption coefficient (1/m), 0 if none
39 real(dp), protected :: phmc_k3 = 0.0_dp
40
41 !> Table for photoionization
42 type(phmc_tbl_t), public, protected :: phmc_tbl
43
44 interface
45 real(dp) function absfunc_t(dist, p_o2, p_h2o)
46 import
47 real(dp), intent(in) :: dist !< Distance
48 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
49 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
50 end function absfunc_t
51 end interface
52
53 ! Public methods
54 public :: phmc_initialize
56 public :: phmc_get_table_air
57 public :: phmc_do_absorption
58 public :: phmc_set_src
59 ! public :: phe_mc_set_src
60
61contains
62
63 !> Initialize photoionization parameters
64 subroutine phmc_initialize(cfg, is_used, source_type)
65 use m_config
66 use m_streamer
67 use m_gas
68 type(cfg_t), intent(inout) :: cfg !< The configuration for the simulation
69 logical, intent(in) :: is_used !< Whether Monte Carlo photoionization is used
70 character(len=*), intent(in) :: source_type
71 integer :: ix
72 real(dp) :: p_o2, p_h2o
73
74 !< [photoi_mc_parameters]
75 call cfg_add_get(cfg, "photoi_mc%physical_photons", phmc_physical_photons, &
76 "Whether physical photons are used")
77 call cfg_add_get(cfg, "photoi_mc%min_weight", phmc_min_weight, &
78 "Minimal photon weight (default: 1.0)")
79 call cfg_add_get(cfg, "photoi_mc%const_dx", phmc_const_dx, &
80 "Whether a constant grid spacing is used for photoionization")
81 call cfg_add_get(cfg, "photoi_mc%min_dx", phmc_min_dx, &
82 "Minimum grid spacing for photoionization")
83 call cfg_add_get(cfg, "photoi_mc%absorp_fac", phmc_absorp_fac, &
84 "At which grid spacing photons are absorbed compared to their mean distance")
85 call cfg_add_get(cfg, "photoi_mc%num_photons", phmc_num_photons, &
86 "Maximum number of discrete photons to use")
87 !< [photoi_mc_parameters]
88
89 if (phmc_absorp_fac <= 0.0_dp) error stop "photoi_mc%absorp_fac <= 0.0"
90 if (phmc_num_photons < 1) error stop "photoi_mc%num_photons < 1"
91
92 if (is_used) then
93 ! Create table for photoionization
94 ix = gas_index("O2")
95 if (ix == -1) then
96 error stop "Photoionization: no oxygen present"
97 else
98 p_o2 = gas_fractions(ix) * gas_pressure
99 end if
100
101 ix = gas_index("H2O")
102 if (ix == -1) then
103 p_h2o = 0.0_dp
104 else
105 p_h2o = gas_fractions(ix) * gas_pressure
106 end if
107
108 select case (source_type)
109 case ("Zheleznyak")
110 phmc_k3 = 0.0_dp
111 ! Standard Zheleznyak model for dry air
112 call phmc_get_table_air(phmc_tbl, p_o2, p_h2o, absfunc_zheleznyak)
113 case ("Naidis_humid")
114 phmc_k3 = (0.26e2_dp / uc_torr_to_bar) * p_h2o
115 ! Naidis model for humid air
116 call phmc_get_table_air(phmc_tbl, p_o2, p_h2o, absfunc_zheleznyak)
117 case ("Aints_humid")
118 phmc_k3 = 0.0_dp
119 ! Aints model for humid air
120 call phmc_get_table_air(phmc_tbl, p_o2, p_h2o, absfunc_aints_humid)
121 case default
122 error stop "Invalid photoionization source_type"
123 end select
124
126
127 if (st_use_dielectric) then
128 if (.not. phmc_const_dx) &
129 error stop "phmc_initialize: with dielectric use const_dx"
130 if (phmc_absorp_fac > 1e-9_dp) &
131 error stop "Use phmc_absorp_fac <= 1e-9 with dielectric"
132 end if
133 end if
134
135 end subroutine phmc_initialize
136
137 !> Print the grid spacing used for the absorption of photons
138 subroutine phmc_print_grid_spacing(dr_base)
139 real(dp), intent(in) :: dr_base
140 real(dp) :: lengthscale, dx
141 integer :: lvl
142
143 ! Get a typical length scale for the absorption of photons
144 lengthscale = lt_get_col(phmc_tbl%tbl, 1, phmc_absorp_fac)
145
146 ! Determine at which level we estimate the photoionization source term. This
147 ! depends on the typical length scale for absorption.
148 lvl = get_lvl_length(dr_base, lengthscale)
149 dx = dr_base * (0.5_dp**(lvl-1))
150
151 if (phmc_const_dx) then
152 write(*, "(A,E12.3)") " Monte-Carlo photoionization spacing: ", dx
153 else
154 write(*, "(A)") " Monte-Carlo photoionization uses adaptive spacing"
155 end if
156 end subroutine phmc_print_grid_spacing
157
158 !> Compute the photonization table for air. If the absorption function is
159 !> f(r), this table contains r as a function of F (the cumulative absorption
160 !> function, between 0 and 1). Later on a distance can be sampled by drawing a
161 !> uniform 0,1 random number and finding the corresponding distance.
162 subroutine phmc_get_table_air(phmc_tbl, p_O2, p_H2O, absfunc)
163 !< The photonization table
164 type(phmc_tbl_t), intent(inout) :: phmc_tbl
165 !> Partial pressure of oxygen (bar)
166 real(dp), intent(in) :: p_o2
167 !> Partial pressure of H2O (bar)
168 real(dp), intent(in) :: p_h2o
169 !> Absorption function to use
170 procedure(absfunc_t) :: absfunc
171
172 ! Large distance, at which all photons have been absorbed
173 real(dp) :: large_distance = 1e3_dp
174 integer :: n
175 integer, parameter :: tbl_size = 500
176 real(dp), allocatable :: fsum(:), dist(:)
177
178 allocate(fsum(tbl_size))
179 allocate(dist(tbl_size))
180
181 call integrate_absfunc(p_o2, p_h2o, absfunc, 1.0_dp, large_distance, &
182 tbl_size, fsum, dist, n)
183
184 ! The full integral can be less than one (e.g. with Naidis_humid).
185 ! Normalize, and store the factor, which will be used when generating
186 ! photons.
187 phmc_tbl%full_integral = fsum(n)
188 fsum(1:n) = fsum(1:n) / fsum(n)
189
190 phmc_tbl%tbl = lt_create(0.0_dp, 1.0_dp, tbl_size, 1)
191 call lt_set_col(phmc_tbl%tbl, 1, fsum(1:n), dist(1:n))
192
193 write(*, "(A,E12.3,A)") " Average photon absorption length", &
194 1e3_dp * sum(dist(1:n-1))/(n-1), " mm"
195 end subroutine phmc_get_table_air
196
197 subroutine integrate_absfunc(p_O2, p_H2O, absfunc, F_max, r_max, n_steps, fsum, r, n)
198 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
199 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
200 procedure(absfunc_t) :: absfunc
201 real(dp), intent(in) :: f_max
202 real(dp), intent(in) :: r_max
203 integer, intent(in) :: n_steps
204 real(dp), intent(inout) :: fsum(n_steps)
205 real(dp), intent(inout) :: r(n_steps)
206 integer, intent(inout) :: n
207 real(dp) :: df, drdf
208
209 df = f_max / (n_steps-1)
210 fsum(1) = 0.0_dp
211 r(1) = 0.0_dp
212
213 do n = 2, n_steps
214 drdf = rk4_drdf(r(n-1), df, p_o2, p_h2o, absfunc)
215 r(n) = r(n-1) + df * drdf
216 fsum(n) = (n-1) * df
217 if (fsum(n) > f_max .or. r(n) > r_max) exit
218 end do
219
220 ! We went past F_max or r_max, or the loop completed fully. Valid data is
221 ! up to index n
222 n = n - 1
223 end subroutine integrate_absfunc
224
225 !> Runge Kutta 4 method to compute dr/dF, where r is the absorption distance
226 !> and F is the cumulative absorption function (between 0 and 1)
227 real(dp) function rk4_drdf(r, df, p_o2, p_h2o, absfunc)
228 real(dp), intent(in) :: r !> Initial point
229 real(dp), intent(in) :: df !> grid size
230 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
231 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
232 procedure(absfunc_t) :: absfunc
233 real(dp) :: drdf
234 real(dp) :: sum_drdf
235 real(dp), parameter :: one_sixth = 1 / 6.0_dp
236
237 ! Step 1 (at initial r)
238 drdf = 1 / absfunc(r, p_o2, p_h2o)
239 sum_drdf = drdf
240
241 ! Step 2 (at initial r + dr/2)
242 drdf = 1 / absfunc(r + 0.5_dp * df * drdf, p_o2, p_h2o)
243 sum_drdf = sum_drdf + 2 * drdf
244
245 ! Step 3 (at initial r + dr/2)
246 drdf = 1 / absfunc(r + 0.5_dp * df * drdf, p_o2, p_h2o)
247 sum_drdf = sum_drdf + 2 * drdf
248
249 ! Step 4 (at initial r + dr)
250 drdf = 1 / absfunc(r + df * drdf, p_o2, p_h2o)
251 sum_drdf = sum_drdf + drdf
252
253 ! Combine r derivatives at steps
254 rk4_drdf = one_sixth * sum_drdf
255 end function rk4_drdf
256
257 !> Absorption function for photons in air according to Zheleznyak's model
258 !>
259 !> Reference: "Photoionization of nitrogen and oxygen mixtures by radiation
260 !> from a gas discharge", Zheleznyak and Mnatsakanian, 1982.
261 real(dp) function absfunc_zheleznyak(dist, p_o2, p_h2o) result(f)
263 real(dp), intent(in) :: dist !< Distance
264 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
265 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
266 real(dp) :: r
267 real(dp), parameter :: k1 = 3.5_dp / uc_torr_to_bar
268 real(dp), parameter :: k2 = 200 / uc_torr_to_bar
269 real(dp), parameter :: eps = epsilon(1.0_dp)
270
271 ! Note: r is here distance times p_O2
272 r = p_o2 * dist
273
274 if (r * k1 < eps) then
275 ! Use Taylor series around r = 0
276 f = (k2 - k1 + 0.5_dp * (k1**2 - k2**2) * r) &
277 * p_o2 / log(k2/k1)
278 else
279 ! Avoid division by zero
280 f = max(1e-100_dp, &
281 (exp(-k1 * r) - exp(-k2 * r)) / (dist * log(k2/k1)))
282 end if
283 end function absfunc_zheleznyak
284
285 !> Absorption function for photons in humid according to Naidis model
286 !>
287 !> Reference: "On photoionization produced by discharges in air", Naidis, 2006
288 real(dp) function absfunc_naidis_humid(dist, p_o2, p_h2o) result(f)
289 real(dp), intent(in) :: dist !< Distance
290 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
291 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
292 real(dp), parameter :: k3 = 0.26e2_dp / uc_torr_to_bar
293
294 f = max(1e-100_dp, &
295 exp(-k3 * p_h2o * dist) * absfunc_zheleznyak(dist, p_o2, p_h2o))
296 end function absfunc_naidis_humid
297
298 !> Absorption function for photons in humid according to Aints model
299 !>
300 !> Reference: "Absorption of Photo-Ionizing Radiation of Corona Discharges in Air",
301 !> Aints et al, 2008. See also Guo et al, PSST 2025.
302 real(dp) function absfunc_aints_humid(dist, p_o2, p_h2o) result(f)
303 real(dp), intent(in) :: dist !< Distance
304 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
305 real(dp), intent(in) :: p_h2o !< Partial pressure of H2O (bar)
306 real(dp), parameter :: k1 = 3.5_dp / uc_torr_to_bar
307 real(dp), parameter :: k2 = 200 / uc_torr_to_bar
308 real(dp), parameter :: k4 = 0.13e2 / uc_torr_to_bar
309 real(dp), parameter :: k5 = 0.57e2 / uc_torr_to_bar
310 real(dp) :: k_a, k_b
311 real(dp), parameter :: eps = epsilon(1.0_dp)
312
313 k_a = k1 * p_o2 + k4 * p_h2o
314 k_b = k2 * p_o2 + k5 * p_h2o
315
316 if (dist * k_a < eps) then
317 ! Use Taylor series around dist = 0
318 f = (k_b - k_a + 0.5_dp * (k_a**2 - k_b**2) * dist) / log(k_b/k_a)
319 else
320 ! Avoid division by zero
321 f = max(1e-100_dp, &
322 (exp(-k_a * dist) - exp(-k_b * dist)) / (dist * log(k_b/k_a)))
323 end if
324 end function absfunc_aints_humid
325
326 !> Determine the lowest level at which the grid spacing is smaller than 'length'.
327 integer function get_lvl_length(dr_base, length)
328 real(dp), intent(in) :: dr_base !< cell spacing at lvl 1
329 real(dp), intent(in) :: length !< Some length
330 real(dp), parameter :: invlog2 = 1 / log(2.0_dp)
331 real(dp) :: ratio
332
333 ratio = dr_base / length
334 if (ratio <= 1) then
335 get_lvl_length = 1
336 else
337 get_lvl_length = 1 + ceiling(log(ratio) * invlog2)
338 end if
339 end function get_lvl_length
340
341 !> As get_lvl_length but with a random choice between lvl and lvl-1
342 integer function get_rlvl_length(dr_base, length, rng)
343 use m_random
344 real(dp), intent(in) :: dr_base !< cell spacing at lvl 1
345 real(dp), intent(in) :: length !< Some length
346 type(rng_t), intent(inout) :: rng !< Random number generator
347 real(dp), parameter :: invlog2 = 1 / log(2.0_dp)
348 real(dp) :: ratio, tmp
349
350 ratio = dr_base / length
351 if (ratio <= 1) then
352 get_rlvl_length = 1
353 else
354 tmp = log(ratio) * invlog2
355 get_rlvl_length = floor(tmp)
356 if (rng%unif_01() < tmp - get_rlvl_length) &
357 get_rlvl_length = get_rlvl_length + 1
358 end if
359 end function get_rlvl_length
360
361 !> Given a list of photon production positions (xyz_in), compute where they
362 !> end up (xyz_out).
363 subroutine phmc_do_absorption(xyz_in, xyz_out, n_dim, n_photons, tbl, prng, k3)
364 use m_lookup_table
365 use m_random
366 use omp_lib
367 !> Input: number of generated photons. Output: number of ionizing photons.
368 !> These can differ when k3 is non-zero.
369 integer, intent(inout) :: n_photons
370 !> Input (x,y,z) values of photon production
371 real(dp), intent(in) :: xyz_in(3, n_photons)
372 !> Output (x,y,z) values of photon absorption
373 real(dp), intent(out) :: xyz_out(3, n_photons)
374 integer, intent(in) :: n_dim !< 2 or 3 dimensional
375 type(lt_t), intent(in) :: tbl !< Lookup table for absorption distance
376 type(prng_t), intent(inout) :: prng !< Random number generator
377 real(dp), intent(in) :: k3 !< non-ionizing absorption coefficient
378 logical, allocatable :: mask(:)
379 integer :: n, proc_id
380 real(dp) :: rr, dist, p_ionization
381
382 if (n_dim < 2 .or. n_dim > 3) error stop "n_dim should be 2 or 3"
383
384 allocate(mask(n_photons))
385
386 !$omp parallel private(n, rr, dist, proc_id, p_ionization)
387 proc_id = 1+omp_get_thread_num()
388
389 !$omp do
390 do n = 1, n_photons
391 rr = prng%rngs(proc_id)%unif_01()
392 ! Sample travel distance
393 dist = lt_get_col(tbl, 1, rr)
394 ! Note that in 2D we ignore the last dimension
395 xyz_out(:, n) = xyz_in(:, n) + prng%rngs(proc_id)%sphere(dist)
396
397 ! Test for non-ionizing absorption
398 if (k3 > 0) then
399 p_ionization = exp(-k3 * dist)
400 mask(n) = (prng%rngs(proc_id)%unif_01() < p_ionization)
401 else
402 mask(n) = .true.
403 end if
404 end do
405 !$omp end parallel
406
407 ! Keep only photons that led to ionization
408 n_photons = 0
409 do n = 1, size(mask)
410 if (mask(n)) then
411 n_photons = n_photons + 1
412 xyz_out(:, n_photons) = xyz_out(:, n)
413 end if
414 end do
415 end subroutine phmc_do_absorption
416
417 ! !> Given a list of photon production positions (xyz_in), compute where they
418 ! !> end up (xyz_out).
419 ! subroutine phe_mc_get_endloc(xyz_in, xyz_out, n_dim, n_photons, prng)
420 ! use m_lookup_table
421 ! use m_random
422 ! use omp_lib
423 ! integer, intent(in) :: n_photons !< Number of photons
424 ! !> Input (x,y,z) values
425 ! real(dp), intent(in) :: xyz_in(3, n_photons)
426 ! !> Output (x,y,z) values
427 ! real(dp), intent(out) :: xyz_out(3, n_photons)
428 ! integer, intent(in) :: n_dim !< 2 or 3 dimensional
429 ! type(PRNG_t), intent(inout) :: prng !< Random number generator
430 ! integer :: n, proc_id
431 ! real(dp) :: rr, dist
432
433 ! !$omp parallel private(n, rr, dist, proc_id)
434 ! proc_id = 1+omp_get_thread_num()
435
436 ! if (n_dim == 2) then
437 ! !$omp do
438 ! do n = 1, n_photons
439 ! rr = prng%rngs(proc_id)%unif_01()
440 ! dist = 0.00001_dp
441 ! ! Pick a random point on a sphere, and ignore the last dimension
442 ! xyz_out(1:3, n) = xyz_in(1:3, n) + &
443 ! prng%rngs(proc_id)%sphere(dist)
444 ! end do
445 ! !$omp end do
446 ! else if (n_dim == 3) then
447 ! !$omp do
448 ! do n = 1, n_photons
449 ! rr = prng%rngs(proc_id)%unif_01()
450 ! dist = 0.00001_dp
451 ! xyz_out(:, n) = xyz_in(:, n) + prng%rngs(proc_id)%sphere(dist)
452 ! end do
453 ! !$omp end do
454 ! else
455 ! print *, "phmc_do_absorption: unknown n_dim", n_dim
456 ! stop
457 ! end if
458 ! !$omp end parallel
459
460 ! end subroutine phe_mc_get_endloc
461
462 !> Set the source term due to photoionization. At most phmc_num_photons
463 !> discrete photons are produced.
464 subroutine phmc_set_src(tree, rng, i_src, i_photo, use_cyl, dt)
465 use m_random
466 use m_lookup_table
467 use m_dielectric
468 use omp_lib
469
470 type(af_t), intent(inout) :: tree !< Tree
471 type(rng_t), intent(inout) :: rng !< Random number generator
472 !> Input variable that contains photon production per cell
473 integer, intent(in) :: i_src
474 !> Input variable that contains photoionization source rate
475 integer, intent(in) :: i_photo
476 !> Use cylindrical coordinates (only in 2D)
477 logical, intent(in) :: use_cyl
478 !> Time step, if present use "physical" photons
479 real(dp), intent(in), optional :: dt
480
481 integer :: lvl, id, nc, min_lvl
482 integer :: ijk, n, n_used
483 integer :: proc_id, n_procs
484 integer :: pho_lvl, max_num_photons
485 real(dp) :: dr(ndim), dt_fac, dist, n_produced
486 real(dp) :: sum_production_rate, pi_lengthscale
487 real(dp), allocatable :: xyz_src(:, :)
488 real(dp), allocatable :: xyz_abs(:, :)
489#if NDIM == 2
490 real(dp) :: tmp, r(3)
491 real(dp), parameter :: pi = acos(-1.0_dp)
492#endif
493 type(prng_t) :: prng
494 type(af_loc_t), allocatable :: ph_loc(:)
495 real(dp), parameter :: small_value = 1e-100_dp
496
497 if (ndim == 3 .and. use_cyl) error stop "phmc_set_src: use_cyl is .true."
498
499 if (st_use_dielectric) then
501 end if
502
503 nc = tree%n_cell
504
505 ! Initialize parallel random numbers
506 n_procs = omp_get_max_threads()
507 call prng%init_parallel(n_procs, rng)
508
509 ! Compute the sum of the photon production rate
510 call af_tree_sum_cc(tree, i_src, sum_production_rate)
511
512 ! dt_fac is used to convert a discrete number of photons to a number of
513 ! photons per unit time. Depending on the arguments dt and phmc_num_photons,
514 ! 'physical' photons with a weight of phmc_min_weight, or super-photons with
515 ! a weight larger (or smaller) than one can be used.
516 if (present(dt)) then
517 ! Create "physical" photons when less than phmc_num_photons are produced,
518 ! otherwise create approximately phmc_num_photons
519 n_produced = dt * sum_production_rate / phmc_min_weight
520
521 if (n_produced < phmc_num_photons) then
522 dt_fac = dt/phmc_min_weight
523 else
524 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
525 end if
526 else
527 ! Create approximately phmc_num_photons by setting dt_fac like this
528 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
529 end if
530
531 ! Over-estimate number of photons because of stochastic production
532 max_num_photons = nint(1.2_dp * dt_fac * sum_production_rate + 1000)
533
534 call phmc_generate_photons(tree, dt_fac, i_src, xyz_src, &
535 max_num_photons, n_used, prng)
536
537 allocate(xyz_abs(3, n_used))
538
539 if (use_cyl) then ! 2D only
540 ! Get location of absorption. On input, xyz is set to (r, z, 0). On
541 ! output, the coordinates thus correspond to (x, z, y)
542 call phmc_do_absorption(xyz_src, xyz_abs, 3, n_used, phmc_tbl%tbl, prng, phmc_k3)
543
544 !$omp do
545 do n = 1, n_used
546 ! Convert back to (r,z) coordinates
547 xyz_abs(1, n) = sqrt(xyz_abs(1, n)**2 + xyz_abs(3, n)**2)
548 xyz_abs(2, n) = xyz_abs(2, n)
549 end do
550 !$omp end do
551
552 if (st_use_dielectric) then
553 ! Handle photons that collide with dielectrics separately
554 call dielectric_photon_absorption(tree, xyz_src, xyz_abs, &
555 2, n_used, 1.0_dp/dt_fac)
556 end if
557 else
558 ! Get location of absorbption
559 call phmc_do_absorption(xyz_src, xyz_abs, ndim, n_used, phmc_tbl%tbl, prng, phmc_k3)
560
561 if (st_use_dielectric) then
562 ! Handle photons that collide with dielectrics separately
563 call dielectric_photon_absorption(tree, xyz_src, xyz_abs, &
564 ndim, n_used, 1.0_dp/dt_fac)
565 end if
566 end if
567
568 call periodify_coordinates(n_used, xyz_abs)
569
570 allocate(ph_loc(n_used))
571
572 if (phmc_const_dx) then
573 ! Get a typical length scale for the absorption of photons
574 pi_lengthscale = lt_get_col(phmc_tbl%tbl, 1, phmc_absorp_fac)
575
576 ! Determine at which level we estimate the photoionization source term. This
577 ! depends on the typical length scale for absorption.
578 pho_lvl = get_lvl_length(maxval(tree%dr_base), pi_lengthscale)
579
580 !$omp parallel do
581 do n = 1, n_used
582 ph_loc(n) = af_get_loc(tree, xyz_abs(1:ndim, n), pho_lvl)
583 end do
584 !$omp end parallel do
585 else
586 !$omp parallel private(n, dist, lvl, proc_id)
587 proc_id = 1+omp_get_thread_num()
588 !$omp do
589 do n = 1, n_used
590 dist = phmc_absorp_fac * norm2(xyz_abs(1:ndim, n) - xyz_src(1:ndim, n))
591 if (dist < phmc_min_dx) dist = phmc_min_dx
592 lvl = get_rlvl_length(maxval(tree%dr_base), dist, prng%rngs(proc_id))
593 ph_loc(n) = af_get_loc(tree, xyz_abs(1:ndim, n), lvl)
594 end do
595 !$omp end do
596 !$omp end parallel
597 end if
598
599 ! Clear variable i_photo, in which we will store the photoionization source term
600 call af_tree_clear_cc(tree, i_photo)
601
602 do n = 1, n_used
603 id = ph_loc(n)%id
604 if (id > af_no_box) then
605#if NDIM == 1
606 i = ph_loc(n)%ix(1)
607 dr = tree%boxes(id)%dr
608
609 tree%boxes(id)%cc(ijk, i_photo) = &
610 tree%boxes(id)%cc(ijk, i_photo) + &
611 phmc_tbl%full_integral/(dt_fac * product(dr))
612#elif NDIM == 2
613 i = ph_loc(n)%ix(1)
614 j = ph_loc(n)%ix(2)
615 dr = tree%boxes(id)%dr
616
617 if (use_cyl) then
618 tmp = dt_fac * 2 * pi
619 r(1:2) = af_r_cc(tree%boxes(id), [i, j])
620 tree%boxes(id)%cc(i, j, i_photo) = &
621 tree%boxes(id)%cc(i, j, i_photo) + &
622 phmc_tbl%full_integral/(tmp * product(dr) * r(1))
623 else
624 tree%boxes(id)%cc(ijk, i_photo) = &
625 tree%boxes(id)%cc(ijk, i_photo) + &
626 phmc_tbl%full_integral/(dt_fac * product(dr))
627 end if
628#elif NDIM == 3
629 i = ph_loc(n)%ix(1)
630 j = ph_loc(n)%ix(2)
631 k = ph_loc(n)%ix(3)
632 dr = tree%boxes(id)%dr
633
634 tree%boxes(id)%cc(ijk, i_photo) = &
635 tree%boxes(id)%cc(ijk, i_photo) + &
636 phmc_tbl%full_integral/(dt_fac * product(dr))
637#endif
638 end if
639 end do
640
641 ! Set ghost cells on highest level with photon source
642 if (phmc_const_dx) then
643 min_lvl = pho_lvl
644 else
645 min_lvl = 1
646 end if
647
648 !$omp parallel private(lvl, i, id)
649 ! Prolong to finer grids
650 do lvl = min_lvl, tree%highest_lvl-1
651 !$omp do
652 do i = 1, size(tree%lvls(lvl)%parents)
653 id = tree%lvls(lvl)%parents(i)
654 call af_gc_box(tree, id, [i_photo])
655 end do
656 !$omp end do
657
658 !$omp do
659 do i = 1, size(tree%lvls(lvl)%parents)
660 id = tree%lvls(lvl)%parents(i)
661 call af_prolong_linear_from(tree%boxes, id, i_photo, add=.true.)
662 end do
663 !$omp end do
664 end do
665 !$omp end parallel
666
667 call prng%update_seed(rng)
668 end subroutine phmc_set_src
669
670 !> Map coordinates to a periodic domain
671 subroutine periodify_coordinates(n_points, coords)
672 integer, intent(in) :: n_points
673 real(dp), intent(inout) :: coords(3, n_points)
674 real(dp) :: xrel
675
676 integer :: i_dim, n
677
678 do i_dim = 1, ndim
679 if (st_periodic(i_dim)) then
680 do n = 1, n_points
681 xrel = coords(i_dim, n) - st_domain_origin(i_dim)
682 coords(i_dim, n) = st_domain_origin(i_dim) + &
683 modulo(xrel, st_domain_len(i_dim))
684 end do
685 end if
686 end do
687 end subroutine periodify_coordinates
688
689! !> Whole photoemission procedure
690! subroutine phe_mc_set_src(tree, rng, i_src, i_photo, use_cyl, dt)
691! use m_random
692! use m_lookup_table
693! use m_dielectric2
694! use omp_lib
695
696! type(af_t), intent(inout) :: tree !< Tree
697! type(RNG_t), intent(inout) :: rng !< Random number generator
698! !> Input variable that contains photon production per cell
699! integer, intent(in) :: i_src
700! !> Input variable that contains photoionization source rate
701! integer, intent(in) :: i_photo
702! !> Use cylindrical coordinates (only in 2D)
703! logical, intent(in) :: use_cyl
704! !> Time step, if present use "physical" photons
705! real(dp), intent(in), optional :: dt
706
707! integer :: nc
708! integer :: n, n_used
709! integer :: proc_id, n_procs
710! integer :: pho_lvl, max_num_photons
711! real(dp) :: tmp, dt_fac, r(3)
712! real(dp) :: sum_production_rate, pi_lengthscale
713! real(dp), allocatable :: xyz_src(:, :)
714! real(dp), allocatable :: xyz_abs(:, :)
715! #if NDIM == 2
716! real(dp), parameter :: pi = acos(-1.0_dp)
717! #endif
718! type(PRNG_t) :: prng
719! type(af_loc_t), allocatable :: ph_loc(:)
720! real(dp), parameter :: small_value = 1e-100_dp
721
722! if (NDIM == 3 .and. use_cyl) error stop "phmc_set_src: use_cyl is .true."
723
724! if (ST_use_dielectric) then
725! call dielectric2_reset_photons()
726! end if
727
728! nc = tree%n_cell
729
730! ! Initialize parallel random numbers
731! n_procs = omp_get_max_threads()
732! call prng%init_parallel(n_procs, rng)
733
734! ! Compute the sum of the photon production rate
735! call af_tree_sum_cc(tree, i_src, sum_production_rate)
736
737! ! dt_fac is used to convert a discrete number of photons to a number of
738! ! photons per unit time. Depending on the arguments dt and phmc_num_photons,
739! ! 'physical' photons with a weight of 1.0, or super-photons with a weight
740! ! larger (or smaller) than one can be used.
741! if (present(dt)) then
742! ! Create "physical" photons when less than phmc_num_photons are produced,
743! ! otherwise create approximately phmc_num_photons
744! dt_fac = min(dt, phmc_num_photons / (sum_production_rate + small_value))
745! else
746! ! Create approximately phmc_num_photons by setting dt_fac like this
747! dt_fac = phmc_num_photons / (sum_production_rate + small_value)
748! end if
749
750! ! Over-estimate number of photons because of stochastic production
751! max_num_photons = nint(1.2_dp * dt_fac * sum_production_rate + 1000)
752
753! call phmc_generate_photons(tree, dt_fac, i_src, xyz_src, &
754! max_num_photons, n_used, prng)
755
756! allocate(xyz_abs(3, n_used))
757
758! if (use_cyl) then ! 2D only
759! ! Get location of absorption. On input, xyz is set to (r, z, 0). On
760! ! output, the coordinates thus correspond to (x, z, y)
761! call phe_mc_get_endloc(xyz_src, xyz_abs, 3, n_used, prng)
762
763! !$omp do
764! do n = 1, n_used
765! ! Convert back to (r,z) coordinates
766! xyz_abs(1, n) = sqrt(xyz_abs(1, n)**2 + xyz_abs(3, n)**2)
767! xyz_abs(2, n) = xyz_abs(2, n)
768! xyz_src(2, n) = xyz_src(3, n)
769! end do
770! !$omp end do
771
772! if (ST_use_dielectric) then
773! ! Handle photons that collide with dielectrics separately
774! call dielectric2_photon_absorption(tree, xyz_src, xyz_abs, &
775! 2, n_used, 1.0_dp/dt_fac)
776! end if
777! else
778! ! Get location of absorbption
779! call phe_mc_get_endloc(xyz_src, xyz_abs, NDIM, n_used, prng)
780
781! if (ST_use_dielectric) then
782! ! Handle photons that collide with dielectrics separately
783! call dielectric2_photon_absorption(tree, xyz_src, xyz_abs, &
784! NDIM, n_used, 1.0_dp/dt_fac)
785! end if
786! end if
787
788! call prng%update_seed(rng)
789! end subroutine phe_mc_set_src
790
791
792 subroutine phmc_generate_photons(tree, dt_fac, i_src, xyz_src, n_max, n_used, prng)
793 use omp_lib
794 use m_random
795 type(af_t), intent(in) :: tree
796 !> Time step to convert photon production rate to number of photons
797 real(dp), intent(in) :: dt_fac
798 !> Index of cell-centered variable containing photon production rate
799 integer, intent(in) :: i_src
800 !> Origin of the produced photons
801 real(dp), allocatable, intent(inout) :: xyz_src(:, :)
802 !> Maximum number of photons that will be produced
803 integer, intent(in) :: n_max
804 !> Number of photons that have been created
805 integer, intent(out) :: n_used
806 !> Parallel random number generator
807 type(prng_t), intent(inout) :: prng
808
809 integer :: proc_id, n_procs
810 integer :: lvl, ix, id, ijk, n, nc, i_ph
811 integer :: n_create, my_num_photons
812 real(dp) :: tmp, dr(ndim), r(3)
813 integer, allocatable :: photons_per_thread(:), photon_thread(:)
814 integer, allocatable :: ix_offset(:)
815 real(dp), allocatable :: xyz_tmp(:, :)
816
817 nc = tree%n_cell
818 n_procs = omp_get_max_threads()
819
820 allocate(xyz_src(3, n_max))
821 allocate(photons_per_thread(n_procs))
822 allocate(photon_thread(n_max))
823
824 ! Loop over all leaves and create photons using random numbers
825 n_used = 0
826
827 !$omp parallel private(lvl, ix, id, IJK, n, r, dr, i_ph, proc_id, &
828 !$omp tmp, n_create, my_num_photons)
829 proc_id = 1+omp_get_thread_num()
830 my_num_photons = 0
831
832 do lvl = 1, tree%highest_lvl
833 dr = af_lvl_dr(tree, lvl)
834 !$omp do
835 do ix = 1, size(tree%lvls(lvl)%leaves)
836 id = tree%lvls(lvl)%leaves(ix)
837
838 do kji_do(1,nc)
839
840#if NDIM == 2
841 if (tree%boxes(id)%coord_t == af_cyl) then
842 tmp = af_cyl_radius_cc(tree%boxes(id), i)
843 tmp = dt_fac * 2 * uc_pi * tmp * &
844 tree%boxes(id)%cc(i, j, i_src) * product(dr)
845 else
846 tmp = dt_fac * tree%boxes(id)%cc(i, j, i_src) * product(dr)
847 end if
848#else
849 tmp = dt_fac * tree%boxes(id)%cc(ijk, i_src) * product(dr)
850#endif
851
852 n_create = floor(tmp)
853
854 if (prng%rngs(proc_id)%unif_01() < tmp - n_create) &
855 n_create = n_create + 1
856
857 if (n_create > 0) then
858 !$omp critical
859 i_ph = n_used
860 n_used = n_used + n_create
861 !$omp end critical
862 my_num_photons = my_num_photons + n_create
863
864 do n = 1, n_create
865 ! Location of production randomly chosen in cell
866 r(1) = prng%rngs(proc_id)%unif_01()
867 r(2) = prng%rngs(proc_id)%unif_01()
868#if NDIM == 2
869 xyz_src(1:ndim, i_ph+n) = af_rr_cc(tree%boxes(id), &
870 [ijk] - 0.5_dp + r(1:ndim))
871 xyz_src(3, i_ph+n) = 0.0_dp
872#elif NDIM == 3
873 r(3) = prng%rngs(proc_id)%unif_01()
874 xyz_src(:, i_ph+n) = af_rr_cc(tree%boxes(id), &
875 [ijk] - 0.5_dp + r(1:ndim))
876#endif
877 photon_thread(i_ph+n) = proc_id
878 end do
879 end if
880 end do; close_do
881 end do
882 !$omp end do nowait
883 end do
884
885 photons_per_thread(proc_id) = my_num_photons
886 !$omp end parallel
887
888 ! Sort the xyz_src array so that runs are deterministic (the order in which
889 ! the threads write is not deterministic)
890 allocate(ix_offset(n_procs))
891 allocate(xyz_tmp(3, n_used))
892
893 ix_offset(1) = 0
894 do n = 2, n_procs
895 ix_offset(n) = ix_offset(n-1) + photons_per_thread(n-1)
896 end do
897
898 photons_per_thread = 0
899 do n = 1, n_used
900 i = photon_thread(n)
901 photons_per_thread(i) = photons_per_thread(i) + 1
902 ix = ix_offset(i) + photons_per_thread(i)
903 xyz_tmp(:, ix) = xyz_src(:, n)
904 end do
905 xyz_src(:, 1:n_used) = xyz_tmp(:, 1:n_used)
906
907 end subroutine phmc_generate_photons
908
909end module m_photoi_mc
Module with settings and routines to handle dielectrics.
subroutine, public dielectric_reset_photons()
subroutine, public dielectric_photon_absorption(tree, xyz_start, xyz_end, n_dim, n_photons, photon_weight)
Determine whether and where photons hit a dielectric, and change their absorption location to the fir...
Module that stores parameters related to the gas.
Definition m_gas.f90:2
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), dimension(:), allocatable, public, protected gas_fractions
Definition m_gas.f90:24
real(dp), public, protected gas_pressure
Definition m_gas.f90:18
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_do_absorption(xyz_in, xyz_out, n_dim, n_photons, tbl, prng, k3)
Given a list of photon production positions (xyz_in), compute where they end up (xyz_out).
subroutine, public phmc_get_table_air(phmc_tbl, p_o2, p_h2o, absfunc)
Compute the photonization table for air. If the absorption function is f(r), this table contains r as...
type(phmc_tbl_t), public, protected phmc_tbl
Table for photoionization.
subroutine, public phmc_initialize(cfg, is_used, source_type)
Initialize photoionization parameters.
subroutine, public phmc_print_grid_spacing(dr_base)
Print the grid spacing used for the absorption of photons.
logical, public, protected phmc_physical_photons
Whether physical photons are used.
This module contains several pre-defined variables like:
Definition m_streamer.f90:6
integer, public, protected st_box_size
The size of the boxes that we use to construct our mesh.
logical, public, protected st_use_dielectric
Whether a dielectric is used.
real(dp), dimension(ndim), public, protected st_domain_len
Domain length per dimension.
Module that contains physical and numerical constants.
real(dp), parameter uc_pi
real(dp), parameter uc_torr_to_bar
Type to quickly look up absorption lengths from a table.