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 real(dp) :: frac_in_tbl !< Fraction photons in table
16 end type phmc_tbl_t
17
18 !> Whether physical photons are used
19 logical, public, protected :: phmc_physical_photons = .true.
20
21 !> Whether a constant or adaptive grid spacing is used
22 logical, protected :: phmc_const_dx = .true.
23
24 !> Minimum grid spacing for photoionization
25 real(dp), protected :: phmc_min_dx = 1e-9_dp
26
27 !> At which grid spacing photons are absorbed compared to their mean distance
28 real(dp), protected :: phmc_absorp_fac = 0.25_dp
29
30 !> Number of photons to use
31 integer, protected :: phmc_num_photons = 5000*1000
32
33 !> Minimal photon weight
34 real(dp), protected :: phmc_min_weight = 1.0_dp
35
36 !> Table for photoionization
37 type(phmc_tbl_t), public, protected :: phmc_tbl
38
39 ! Public methods
40 public :: phmc_initialize
42 public :: phmc_get_table_air
43 public :: phmc_do_absorption
45 public :: phmc_set_src
46 ! public :: phe_mc_set_src
47
48contains
49
50 !> Initialize photoionization parameters
51 subroutine phmc_initialize(cfg, is_used)
52 use m_config
53 use m_streamer
54 use m_gas
55 type(cfg_t), intent(inout) :: cfg !< The configuration for the simulation
56 logical, intent(in) :: is_used !< Whether Monte Carlo photoionization is used
57 integer :: ix
58
59 !< [photoi_mc_parameters]
60 call cfg_add_get(cfg, "photoi_mc%physical_photons", phmc_physical_photons, &
61 "Whether physical photons are used")
62 call cfg_add_get(cfg, "photoi_mc%min_weight", phmc_min_weight, &
63 "Minimal photon weight (default: 1.0)")
64 call cfg_add_get(cfg, "photoi_mc%const_dx", phmc_const_dx, &
65 "Whether a constant grid spacing is used for photoionization")
66 call cfg_add_get(cfg, "photoi_mc%min_dx", phmc_min_dx, &
67 "Minimum grid spacing for photoionization")
68 call cfg_add_get(cfg, "photoi_mc%absorp_fac", phmc_absorp_fac, &
69 "At which grid spacing photons are absorbed compared to their mean distance")
70 call cfg_add_get(cfg, "photoi_mc%num_photons", phmc_num_photons, &
71 "Maximum number of discrete photons to use")
72 !< [photoi_mc_parameters]
73
74 if (phmc_absorp_fac <= 0.0_dp) error stop "photoi_mc%absorp_fac <= 0.0"
75 if (phmc_num_photons < 1) error stop "photoi_mc%num_photons < 1"
76
77 if (is_used) then
78 ! Create table for photoionization
79 ix = gas_index("O2")
80 if (ix == -1) error stop "Photoionization: no oxygen present"
82 2 * maxval(st_domain_len))
83
85
86 if (st_use_dielectric) then
87 if (.not. phmc_const_dx) &
88 error stop "phmc_initialize: with dielectric use const_dx"
89 if (phmc_absorp_fac > 1e-9_dp) &
90 error stop "Use phmc_absorp_fac <= 1e-9 with dielectric"
91 end if
92 end if
93
94 end subroutine phmc_initialize
95
96 !> Print the grid spacing used for the absorption of photons
97 subroutine phmc_print_grid_spacing(dr_base)
98 real(dp), intent(in) :: dr_base
99 real(dp) :: lengthscale, dx
100 integer :: lvl
101
102 ! Get a typical length scale for the absorption of photons
103 lengthscale = lt_get_col(phmc_tbl%tbl, 1, phmc_absorp_fac)
104
105 ! Determine at which level we estimate the photoionization source term. This
106 ! depends on the typical length scale for absorption.
107 lvl = get_lvl_length(dr_base, lengthscale)
108 dx = dr_base * (0.5_dp**(lvl-1))
109
110 if (phmc_const_dx) then
111 write(*, "(A,E12.3)") " Monte-Carlo photoionization spacing: ", dx
112 else
113 write(*, "(A)") " Monte-Carlo photoionization uses adaptive spacing"
114 end if
115 end subroutine phmc_print_grid_spacing
116
117 !> Compute the photonization table for air. If the absorption function is
118 !> f(r), this table contains r as a function of F (the cumulative absorption
119 !> function, between 0 and 1). Later on a distance can be sampled by drawing a
120 !> uniform 0,1 random number and finding the corresponding distance. The table
121 !> constructed up to max_dist; we can ignore photons that fly very far.
122 subroutine phmc_get_table_air(phmc_tbl, p_O2, max_dist)
123 !< The photonization table
124 type(phmc_tbl_t), intent(inout) :: phmc_tbl
125 !< Partial pressure of oxygen (bar)
126 real(dp), intent(in) :: p_o2
127 !< Maximum distance in lookup table
128 real(dp), intent(in) :: max_dist
129
130 integer :: n
131 integer, parameter :: tbl_size = 500
132 real(dp), allocatable :: fsum(:), dist(:)
133 real(dp) :: df, drdf, r, f, fmax_guess
134
135 ! First estimate which fraction of photons are within max_dist (start with
136 ! the upper bound of 1.0)
137 fmax_guess = 1.0_dp
138
139 ! 5 loops should be enough for a good guess
140 do n = 1, 5
141 ! When Fmax_guess decreases so will dF, since we keep the number of
142 ! points the same
143 df = fmax_guess / (tbl_size-1)
144 r = 0
145 f = 0
146
147 do
148 drdf = rk4_drdf(r, df, p_o2)
149 r = r + df * drdf
150 f = f + df
151
152 if (r > max_dist) then
153 ! A better estimate for the upper bound
154 fmax_guess = f
155 exit
156 end if
157 end do
158 end do
159
160 ! Make arrays larger so that we surely are in range of maxdist
161 allocate(fsum(2 * tbl_size))
162 allocate(dist(2 * tbl_size))
163
164 ! Now create table
165 df = fmax_guess / (tbl_size-1)
166 dist(1) = 0
167 fsum(1) = 0
168
169 ! Compute r(F) for F = dF, 2 dF, 3 dF, ...
170 do n = 2, 2 * tbl_size
171 drdf = rk4_drdf(dist(n-1), df, p_o2)
172 fsum(n) = fsum(n-1) + df
173 dist(n) = dist(n-1) + df * drdf
174 if (dist(n) > max_dist) exit
175 end do
176
177 if (n > tbl_size + 10) &
178 stop "phmc_get_table_air: integration accuracy fail"
179
180 if (st_use_dielectric) then
181 ! Photons that fly outside the domain can be absorbed by the dielectric,
182 ! so don't apply any tricks
183 phmc_tbl%frac_in_tbl = 1.0_dp
184 else
185 ! Scale table to lie between 0 and 1 (later on we have to correct for this)
186 phmc_tbl%frac_in_tbl = fsum(n-1)
187 fsum(1:n-1) = fsum(1:n-1) / fsum(n-1)
188 end if
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-1), dist(1:n-1))
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 !> Runge Kutta 4 method to compute dr/dF, where r is the absorption distance
198 !> and F is the cumulative absorption function (between 0 and 1)
199 real(dp) function rk4_drdf(r, df, p_o2)
200 real(dp), intent(in) :: r !> Initial point
201 real(dp), intent(in) :: df !> grid size
202 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
203 real(dp) :: drdf
204 real(dp) :: sum_drdf
205 real(dp), parameter :: one_sixth = 1 / 6.0_dp
206
207 ! Step 1 (at initial r)
208 drdf = 1 / phmc_absorption_func_air(r, p_o2)
209 sum_drdf = drdf
210
211 ! Step 2 (at initial r + dr/2)
212 drdf = 1 / phmc_absorption_func_air(r + 0.5_dp * df * drdf, p_o2)
213 sum_drdf = sum_drdf + 2 * drdf
214
215 ! Step 3 (at initial r + dr/2)
216 drdf = 1 / phmc_absorption_func_air(r + 0.5_dp * df * drdf, p_o2)
217 sum_drdf = sum_drdf + 2 * drdf
218
219 ! Step 4 (at initial r + dr)
220 drdf = 1 / phmc_absorption_func_air(r + df * drdf, p_o2)
221 sum_drdf = sum_drdf + drdf
222
223 ! Combine r derivatives at steps
224 rk4_drdf = one_sixth * sum_drdf
225 end function rk4_drdf
226
227 !> The absorption function for photons in air according to Zheleznyak's model
228 real(dp) function phmc_absorption_func_air(dist, p_o2)
230 real(dp), intent(in) :: dist !< Distance
231 real(dp), intent(in) :: p_o2 !< Partial pressure of oxygen (bar)
232 real(dp) :: r
233 real(dp), parameter :: c0 = 3.5_dp / uc_torr_to_bar
234 real(dp), parameter :: c1 = 200 / uc_torr_to_bar
235 real(dp), parameter :: eps = epsilon(1.0_dp)
236
237 r = p_o2 * dist
238 if (r * (c0 + c1) < eps) then
239 ! Use limit to prevent over/underflow
240 phmc_absorption_func_air = (c1 - c0 + 0.5_dp * (c0**2 - c1**2) * r) &
241 * p_o2 / log(c1/c0)
242 else if (r * c0 > -log(eps)) then
243 ! Use limit to prevent over/underflow
245 else
246 phmc_absorption_func_air = (exp(-c0 * r) - exp(-c1 * r)) / (dist * log(c1/c0))
247 end if
248 end function phmc_absorption_func_air
249
250 !> Determine the lowest level at which the grid spacing is smaller than 'length'.
251 integer function get_lvl_length(dr_base, length)
252 real(dp), intent(in) :: dr_base !< cell spacing at lvl 1
253 real(dp), intent(in) :: length !< Some length
254 real(dp), parameter :: invlog2 = 1 / log(2.0_dp)
255 real(dp) :: ratio
256
257 ratio = dr_base / length
258 if (ratio <= 1) then
259 get_lvl_length = 1
260 else
261 get_lvl_length = 1 + ceiling(log(ratio) * invlog2)
262 end if
263 end function get_lvl_length
264
265 !> As get_lvl_length but with a random choice between lvl and lvl-1
266 integer function get_rlvl_length(dr_base, length, rng)
267 use m_random
268 real(dp), intent(in) :: dr_base !< cell spacing at lvl 1
269 real(dp), intent(in) :: length !< Some length
270 type(rng_t), intent(inout) :: rng !< Random number generator
271 real(dp), parameter :: invlog2 = 1 / log(2.0_dp)
272 real(dp) :: ratio, tmp
273
274 ratio = dr_base / length
275 if (ratio <= 1) then
276 get_rlvl_length = 1
277 else
278 tmp = log(ratio) * invlog2
279 get_rlvl_length = floor(tmp)
280 if (rng%unif_01() < tmp - get_rlvl_length) &
281 get_rlvl_length = get_rlvl_length + 1
282 end if
283 end function get_rlvl_length
284
285 !> Given a list of photon production positions (xyz_in), compute where they
286 !> end up (xyz_out).
287 subroutine phmc_do_absorption(xyz_in, xyz_out, n_dim, n_photons, tbl, prng)
288 use m_lookup_table
289 use m_random
290 use omp_lib
291 integer, intent(in) :: n_photons !< Number of photons
292 !> Input (x,y,z) values
293 real(dp), intent(in) :: xyz_in(3, n_photons)
294 !> Output (x,y,z) values
295 real(dp), intent(out) :: xyz_out(3, n_photons)
296 integer, intent(in) :: n_dim !< 2 or 3 dimensional
297 !< Lookup table
298 type(lt_t), intent(in) :: tbl
299 type(prng_t), intent(inout) :: prng !< Random number generator
300 integer :: n, proc_id
301 real(dp) :: rr, dist
302
303 !$omp parallel private(n, rr, dist, proc_id)
304 proc_id = 1+omp_get_thread_num()
305
306 if (n_dim == 2) then
307 !$omp do
308 do n = 1, n_photons
309 rr = prng%rngs(proc_id)%unif_01()
310 dist = lt_get_col(tbl, 1, rr)
311 ! Pick a random point on a sphere, and ignore the last dimension
312 xyz_out(1:3, n) = xyz_in(1:3, n) + &
313 prng%rngs(proc_id)%sphere(dist)
314 end do
315 !$omp end do
316 else if (n_dim == 3) then
317 !$omp do
318 do n = 1, n_photons
319 rr = prng%rngs(proc_id)%unif_01()
320 dist = lt_get_col(tbl, 1, rr)
321 xyz_out(:, n) = xyz_in(:, n) + prng%rngs(proc_id)%sphere(dist)
322 end do
323 !$omp end do
324 else
325 print *, "phmc_do_absorption: unknown n_dim", n_dim
326 stop
327 end if
328 !$omp end parallel
329
330 end subroutine phmc_do_absorption
331
332 ! !> Given a list of photon production positions (xyz_in), compute where they
333 ! !> end up (xyz_out).
334 ! subroutine phe_mc_get_endloc(xyz_in, xyz_out, n_dim, n_photons, prng)
335 ! use m_lookup_table
336 ! use m_random
337 ! use omp_lib
338 ! integer, intent(in) :: n_photons !< Number of photons
339 ! !> Input (x,y,z) values
340 ! real(dp), intent(in) :: xyz_in(3, n_photons)
341 ! !> Output (x,y,z) values
342 ! real(dp), intent(out) :: xyz_out(3, n_photons)
343 ! integer, intent(in) :: n_dim !< 2 or 3 dimensional
344 ! type(PRNG_t), intent(inout) :: prng !< Random number generator
345 ! integer :: n, proc_id
346 ! real(dp) :: rr, dist
347
348 ! !$omp parallel private(n, rr, dist, proc_id)
349 ! proc_id = 1+omp_get_thread_num()
350
351 ! if (n_dim == 2) then
352 ! !$omp do
353 ! do n = 1, n_photons
354 ! rr = prng%rngs(proc_id)%unif_01()
355 ! dist = 0.00001_dp
356 ! ! Pick a random point on a sphere, and ignore the last dimension
357 ! xyz_out(1:3, n) = xyz_in(1:3, n) + &
358 ! prng%rngs(proc_id)%sphere(dist)
359 ! end do
360 ! !$omp end do
361 ! else if (n_dim == 3) then
362 ! !$omp do
363 ! do n = 1, n_photons
364 ! rr = prng%rngs(proc_id)%unif_01()
365 ! dist = 0.00001_dp
366 ! xyz_out(:, n) = xyz_in(:, n) + prng%rngs(proc_id)%sphere(dist)
367 ! end do
368 ! !$omp end do
369 ! else
370 ! print *, "phmc_do_absorption: unknown n_dim", n_dim
371 ! stop
372 ! end if
373 ! !$omp end parallel
374
375 ! end subroutine phe_mc_get_endloc
376
377 !> Set the source term due to photoionization. At most phmc_num_photons
378 !> discrete photons are produced.
379 subroutine phmc_set_src(tree, rng, i_src, i_photo, use_cyl, dt)
380 use m_random
381 use m_lookup_table
382 use m_dielectric
383 use omp_lib
384
385 type(af_t), intent(inout) :: tree !< Tree
386 type(rng_t), intent(inout) :: rng !< Random number generator
387 !> Input variable that contains photon production per cell
388 integer, intent(in) :: i_src
389 !> Input variable that contains photoionization source rate
390 integer, intent(in) :: i_photo
391 !> Use cylindrical coordinates (only in 2D)
392 logical, intent(in) :: use_cyl
393 !> Time step, if present use "physical" photons
394 real(dp), intent(in), optional :: dt
395
396 integer :: lvl, id, nc, min_lvl
397 integer :: ijk, n, n_used
398 integer :: proc_id, n_procs
399 integer :: pho_lvl, max_num_photons
400 real(dp) :: dr(ndim), dt_fac, dist, n_produced
401 real(dp) :: sum_production_rate, pi_lengthscale
402 real(dp), allocatable :: xyz_src(:, :)
403 real(dp), allocatable :: xyz_abs(:, :)
404#if NDIM == 2
405 real(dp) :: tmp, r(3)
406 real(dp), parameter :: pi = acos(-1.0_dp)
407#endif
408 type(prng_t) :: prng
409 type(af_loc_t), allocatable :: ph_loc(:)
410 real(dp), parameter :: small_value = 1e-100_dp
411
412 if (ndim == 3 .and. use_cyl) error stop "phmc_set_src: use_cyl is .true."
413
414 if (st_use_dielectric) then
416 end if
417
418 nc = tree%n_cell
419
420 ! Initialize parallel random numbers
421 n_procs = omp_get_max_threads()
422 call prng%init_parallel(n_procs, rng)
423
424 ! Compute the sum of the photon production rate
425 call af_tree_sum_cc(tree, i_src, sum_production_rate)
426
427 ! dt_fac is used to convert a discrete number of photons to a number of
428 ! photons per unit time. Depending on the arguments dt and phmc_num_photons,
429 ! 'physical' photons with a weight of phmc_min_weight, or super-photons with
430 ! a weight larger (or smaller) than one can be used.
431 if (present(dt)) then
432 ! Create "physical" photons when less than phmc_num_photons are produced,
433 ! otherwise create approximately phmc_num_photons
434 n_produced = dt * sum_production_rate / phmc_min_weight
435
436 if (n_produced < phmc_num_photons) then
437 dt_fac = dt/phmc_min_weight
438 else
439 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
440 end if
441 else
442 ! Create approximately phmc_num_photons by setting dt_fac like this
443 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
444 end if
445
446 ! Over-estimate number of photons because of stochastic production
447 max_num_photons = nint(1.2_dp * dt_fac * sum_production_rate + 1000)
448
449 call phmc_generate_photons(tree, dt_fac, i_src, xyz_src, &
450 max_num_photons, n_used, prng)
451
452 allocate(xyz_abs(3, n_used))
453
454 if (use_cyl) then ! 2D only
455 ! Get location of absorption. On input, xyz is set to (r, z, 0). On
456 ! output, the coordinates thus correspond to (x, z, y)
457 call phmc_do_absorption(xyz_src, xyz_abs, 3, n_used, phmc_tbl%tbl, prng)
458
459 !$omp do
460 do n = 1, n_used
461 ! Convert back to (r,z) coordinates
462 xyz_abs(1, n) = sqrt(xyz_abs(1, n)**2 + xyz_abs(3, n)**2)
463 xyz_abs(2, n) = xyz_abs(2, n)
464 end do
465 !$omp end do
466
467 if (st_use_dielectric) then
468 ! Handle photons that collide with dielectrics separately
469 call dielectric_photon_absorption(tree, xyz_src, xyz_abs, &
470 2, n_used, 1.0_dp/dt_fac)
471 end if
472 else
473 ! Get location of absorbption
474 call phmc_do_absorption(xyz_src, xyz_abs, ndim, n_used, phmc_tbl%tbl, prng)
475
476 if (st_use_dielectric) then
477 ! Handle photons that collide with dielectrics separately
478 call dielectric_photon_absorption(tree, xyz_src, xyz_abs, &
479 ndim, n_used, 1.0_dp/dt_fac)
480 end if
481 end if
482
483 allocate(ph_loc(n_used))
484
485 if (phmc_const_dx) then
486 ! Get a typical length scale for the absorption of photons
487 pi_lengthscale = lt_get_col(phmc_tbl%tbl, 1, phmc_absorp_fac)
488
489 ! Determine at which level we estimate the photoionization source term. This
490 ! depends on the typical length scale for absorption.
491 pho_lvl = get_lvl_length(maxval(tree%dr_base), pi_lengthscale)
492
493 !$omp parallel do
494 do n = 1, n_used
495 ph_loc(n) = af_get_loc(tree, xyz_abs(1:ndim, n), pho_lvl)
496 end do
497 !$omp end parallel do
498 else
499 !$omp parallel private(n, dist, lvl, proc_id)
500 proc_id = 1+omp_get_thread_num()
501 !$omp do
502 do n = 1, n_used
503 dist = phmc_absorp_fac * norm2(xyz_abs(1:ndim, n) - xyz_src(1:ndim, n))
504 if (dist < phmc_min_dx) dist = phmc_min_dx
505 lvl = get_rlvl_length(maxval(tree%dr_base), dist, prng%rngs(proc_id))
506 ph_loc(n) = af_get_loc(tree, xyz_abs(1:ndim, n), lvl)
507 end do
508 !$omp end do
509 !$omp end parallel
510 end if
511
512 ! Clear variable i_photo, in which we will store the photoionization source term
513 call af_tree_clear_cc(tree, i_photo)
514
515 do n = 1, n_used
516 id = ph_loc(n)%id
517 if (id > af_no_box) then
518#if NDIM == 1
519 i = ph_loc(n)%ix(1)
520 dr = tree%boxes(id)%dr
521
522 tree%boxes(id)%cc(ijk, i_photo) = &
523 tree%boxes(id)%cc(ijk, i_photo) + &
524 phmc_tbl%frac_in_tbl/(dt_fac * product(dr))
525#elif NDIM == 2
526 i = ph_loc(n)%ix(1)
527 j = ph_loc(n)%ix(2)
528 dr = tree%boxes(id)%dr
529
530 if (use_cyl) then
531 tmp = dt_fac * 2 * pi
532 r(1:2) = af_r_cc(tree%boxes(id), [i, j])
533 tree%boxes(id)%cc(i, j, i_photo) = &
534 tree%boxes(id)%cc(i, j, i_photo) + &
535 phmc_tbl%frac_in_tbl/(tmp * product(dr) * r(1))
536 else
537 tree%boxes(id)%cc(ijk, i_photo) = &
538 tree%boxes(id)%cc(ijk, i_photo) + &
539 phmc_tbl%frac_in_tbl/(dt_fac * product(dr))
540 end if
541#elif NDIM == 3
542 i = ph_loc(n)%ix(1)
543 j = ph_loc(n)%ix(2)
544 k = ph_loc(n)%ix(3)
545 dr = tree%boxes(id)%dr
546
547 tree%boxes(id)%cc(ijk, i_photo) = &
548 tree%boxes(id)%cc(ijk, i_photo) + &
549 phmc_tbl%frac_in_tbl/(dt_fac * product(dr))
550#endif
551 end if
552 end do
553
554 ! Set ghost cells on highest level with photon source
555 if (phmc_const_dx) then
556 min_lvl = pho_lvl
557 else
558 min_lvl = 1
559 end if
560
561 !$omp parallel private(lvl, i, id)
562 ! Prolong to finer grids
563 do lvl = min_lvl, tree%highest_lvl-1
564 !$omp do
565 do i = 1, size(tree%lvls(lvl)%parents)
566 id = tree%lvls(lvl)%parents(i)
567 call af_gc_box(tree, id, [i_photo])
568 end do
569 !$omp end do
570
571 !$omp do
572 do i = 1, size(tree%lvls(lvl)%parents)
573 id = tree%lvls(lvl)%parents(i)
574 call af_prolong_linear_from(tree%boxes, id, i_photo, add=.true.)
575 end do
576 !$omp end do
577 end do
578 !$omp end parallel
579
580 call prng%update_seed(rng)
581 end subroutine phmc_set_src
582
583! !> Whole photoemission procedure
584! subroutine phe_mc_set_src(tree, rng, i_src, i_photo, use_cyl, dt)
585! use m_random
586! use m_lookup_table
587! use m_dielectric2
588! use omp_lib
589
590! type(af_t), intent(inout) :: tree !< Tree
591! type(RNG_t), intent(inout) :: rng !< Random number generator
592! !> Input variable that contains photon production per cell
593! integer, intent(in) :: i_src
594! !> Input variable that contains photoionization source rate
595! integer, intent(in) :: i_photo
596! !> Use cylindrical coordinates (only in 2D)
597! logical, intent(in) :: use_cyl
598! !> Time step, if present use "physical" photons
599! real(dp), intent(in), optional :: dt
600
601! integer :: nc
602! integer :: n, n_used
603! integer :: proc_id, n_procs
604! integer :: pho_lvl, max_num_photons
605! real(dp) :: tmp, dt_fac, r(3)
606! real(dp) :: sum_production_rate, pi_lengthscale
607! real(dp), allocatable :: xyz_src(:, :)
608! real(dp), allocatable :: xyz_abs(:, :)
609! #if NDIM == 2
610! real(dp), parameter :: pi = acos(-1.0_dp)
611! #endif
612! type(PRNG_t) :: prng
613! type(af_loc_t), allocatable :: ph_loc(:)
614! real(dp), parameter :: small_value = 1e-100_dp
615
616! if (NDIM == 3 .and. use_cyl) error stop "phmc_set_src: use_cyl is .true."
617
618! if (ST_use_dielectric) then
619! call dielectric2_reset_photons()
620! end if
621
622! nc = tree%n_cell
623
624! ! Initialize parallel random numbers
625! n_procs = omp_get_max_threads()
626! call prng%init_parallel(n_procs, rng)
627
628! ! Compute the sum of the photon production rate
629! call af_tree_sum_cc(tree, i_src, sum_production_rate)
630
631! ! dt_fac is used to convert a discrete number of photons to a number of
632! ! photons per unit time. Depending on the arguments dt and phmc_num_photons,
633! ! 'physical' photons with a weight of 1.0, or super-photons with a weight
634! ! larger (or smaller) than one can be used.
635! if (present(dt)) then
636! ! Create "physical" photons when less than phmc_num_photons are produced,
637! ! otherwise create approximately phmc_num_photons
638! dt_fac = min(dt, phmc_num_photons / (sum_production_rate + small_value))
639! else
640! ! Create approximately phmc_num_photons by setting dt_fac like this
641! dt_fac = phmc_num_photons / (sum_production_rate + small_value)
642! end if
643
644! ! Over-estimate number of photons because of stochastic production
645! max_num_photons = nint(1.2_dp * dt_fac * sum_production_rate + 1000)
646
647! call phmc_generate_photons(tree, dt_fac, i_src, xyz_src, &
648! max_num_photons, n_used, prng)
649
650! allocate(xyz_abs(3, n_used))
651
652! if (use_cyl) then ! 2D only
653! ! Get location of absorption. On input, xyz is set to (r, z, 0). On
654! ! output, the coordinates thus correspond to (x, z, y)
655! call phe_mc_get_endloc(xyz_src, xyz_abs, 3, n_used, prng)
656
657! !$omp do
658! do n = 1, n_used
659! ! Convert back to (r,z) coordinates
660! xyz_abs(1, n) = sqrt(xyz_abs(1, n)**2 + xyz_abs(3, n)**2)
661! xyz_abs(2, n) = xyz_abs(2, n)
662! xyz_src(2, n) = xyz_src(3, n)
663! end do
664! !$omp end do
665
666! if (ST_use_dielectric) then
667! ! Handle photons that collide with dielectrics separately
668! call dielectric2_photon_absorption(tree, xyz_src, xyz_abs, &
669! 2, n_used, 1.0_dp/dt_fac)
670! end if
671! else
672! ! Get location of absorbption
673! call phe_mc_get_endloc(xyz_src, xyz_abs, NDIM, n_used, prng)
674
675! if (ST_use_dielectric) then
676! ! Handle photons that collide with dielectrics separately
677! call dielectric2_photon_absorption(tree, xyz_src, xyz_abs, &
678! NDIM, n_used, 1.0_dp/dt_fac)
679! end if
680! end if
681
682! call prng%update_seed(rng)
683! end subroutine phe_mc_set_src
684
685
686 subroutine phmc_generate_photons(tree, dt_fac, i_src, xyz_src, n_max, n_used, prng)
687 use omp_lib
688 use m_random
689 type(af_t), intent(in) :: tree
690 !> Time step to convert photon production rate to number of photons
691 real(dp), intent(in) :: dt_fac
692 !> Index of cell-centered variable containing photon production rate
693 integer, intent(in) :: i_src
694 !> Origin of the produced photons
695 real(dp), allocatable, intent(inout) :: xyz_src(:, :)
696 !> Maximum number of photons that will be produced
697 integer, intent(in) :: n_max
698 !> Number of photons that have been created
699 integer, intent(out) :: n_used
700 !> Parallel random number generator
701 type(prng_t), intent(inout) :: prng
702
703 integer :: proc_id, n_procs
704 integer :: lvl, ix, id, ijk, n, nc, i_ph
705 integer :: n_create, my_num_photons
706 real(dp) :: tmp, dr(ndim), r(3)
707 integer, allocatable :: photons_per_thread(:), photon_thread(:)
708 integer, allocatable :: ix_offset(:)
709 real(dp), allocatable :: xyz_tmp(:, :)
710
711 nc = tree%n_cell
712 n_procs = omp_get_max_threads()
713
714 allocate(xyz_src(3, n_max))
715 allocate(photons_per_thread(n_procs))
716 allocate(photon_thread(n_max))
717
718 ! Loop over all leaves and create photons using random numbers
719 n_used = 0
720
721 !$omp parallel private(lvl, ix, id, IJK, n, r, dr, i_ph, proc_id, &
722 !$omp tmp, n_create, my_num_photons)
723 proc_id = 1+omp_get_thread_num()
724 my_num_photons = 0
725
726 do lvl = 1, tree%highest_lvl
727 dr = af_lvl_dr(tree, lvl)
728 !$omp do
729 do ix = 1, size(tree%lvls(lvl)%leaves)
730 id = tree%lvls(lvl)%leaves(ix)
731
732 do kji_do(1,nc)
733
734#if NDIM == 2
735 if (tree%boxes(id)%coord_t == af_cyl) then
736 tmp = af_cyl_radius_cc(tree%boxes(id), i)
737 tmp = dt_fac * 2 * uc_pi * tmp * &
738 tree%boxes(id)%cc(i, j, i_src) * product(dr)
739 else
740 tmp = dt_fac * tree%boxes(id)%cc(i, j, i_src) * product(dr)
741 end if
742#else
743 tmp = dt_fac * tree%boxes(id)%cc(ijk, i_src) * product(dr)
744#endif
745
746 n_create = floor(tmp)
747
748 if (prng%rngs(proc_id)%unif_01() < tmp - n_create) &
749 n_create = n_create + 1
750
751 if (n_create > 0) then
752 !$omp critical
753 i_ph = n_used
754 n_used = n_used + n_create
755 !$omp end critical
756 my_num_photons = my_num_photons + n_create
757
758 do n = 1, n_create
759 ! Location of production randomly chosen in cell
760 r(1) = prng%rngs(proc_id)%unif_01()
761 r(2) = prng%rngs(proc_id)%unif_01()
762#if NDIM == 2
763 xyz_src(1:ndim, i_ph+n) = af_rr_cc(tree%boxes(id), &
764 [ijk] - 0.5_dp + r(1:ndim))
765 xyz_src(3, i_ph+n) = 0.0_dp
766#elif NDIM == 3
767 r(3) = prng%rngs(proc_id)%unif_01()
768 xyz_src(:, i_ph+n) = af_rr_cc(tree%boxes(id), &
769 [ijk] - 0.5_dp + r(1:ndim))
770#endif
771 photon_thread(i_ph+n) = proc_id
772 end do
773 end if
774 end do; close_do
775 end do
776 !$omp end do nowait
777 end do
778
779 photons_per_thread(proc_id) = my_num_photons
780 !$omp end parallel
781
782 ! Sort the xyz_src array so that runs are deterministic (the order in which
783 ! the threads write is not deterministic)
784 allocate(ix_offset(n_procs))
785 allocate(xyz_tmp(3, n_used))
786
787 ix_offset(1) = 0
788 do n = 2, n_procs
789 ix_offset(n) = ix_offset(n-1) + photons_per_thread(n-1)
790 end do
791
792 photons_per_thread = 0
793 do n = 1, n_used
794 i = photon_thread(n)
795 photons_per_thread(i) = photons_per_thread(i) + 1
796 ix = ix_offset(i) + photons_per_thread(i)
797 xyz_tmp(:, ix) = xyz_src(:, n)
798 end do
799 xyz_src(:, 1:n_used) = xyz_tmp(:, 1:n_used)
800
801 end subroutine phmc_generate_photons
802
803end 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.
real(dp) function, public phmc_absorption_func_air(dist, p_o2)
The absorption function for photons in air according to Zheleznyak's model.
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.
type(phmc_tbl_t), public, protected phmc_tbl
Table for photoionization.
subroutine, public phmc_do_absorption(xyz_in, xyz_out, n_dim, n_photons, tbl, prng)
Given a list of photon production positions (xyz_in), compute where they end up (xyz_out).
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.
subroutine, public phmc_get_table_air(phmc_tbl, p_o2, max_dist)
Compute the photonization table for air. If the absorption function is f(r), this table contains r as...
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.