3#include "../afivo/src/cpp_macros.h"
17 real(dp) :: full_integral
24 logical,
protected :: phmc_const_dx = .true.
27 real(dp),
protected :: phmc_min_dx = 1e-9_dp
30 real(dp),
protected :: phmc_absorp_fac = 0.25_dp
33 integer,
protected :: phmc_num_photons = 5000*1000
36 real(dp),
protected :: phmc_min_weight = 1.0_dp
39 real(dp),
protected :: phmc_k3 = 0.0_dp
45 real(dp) function absfunc_t(dist, p_o2, p_h2o)
47 real(dp),
intent(in) :: dist
48 real(dp),
intent(in) :: p_o2
49 real(dp),
intent(in) :: p_h2o
50 end function absfunc_t
68 type(cfg_t),
intent(inout) :: cfg
69 logical,
intent(in) :: is_used
70 character(len=*),
intent(in) :: source_type
72 real(dp) :: p_o2, p_h2o
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")
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"
96 error stop
"Photoionization: no oxygen present"
108 select case (source_type)
113 case (
"Naidis_humid")
114 phmc_k3 = (0.26e2_dp / uc_torr_to_bar) * p_h2o
122 error stop
"Invalid photoionization source_type"
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"
139 real(dp),
intent(in) :: dr_base
140 real(dp) :: lengthscale, dx
144 lengthscale = lt_get_col(
phmc_tbl%tbl, 1, phmc_absorp_fac)
148 lvl = get_lvl_length(dr_base, lengthscale)
149 dx = dr_base * (0.5_dp**(lvl-1))
151 if (phmc_const_dx)
then
152 write(*,
"(A,E12.3)")
" Monte-Carlo photoionization spacing: ", dx
154 write(*,
"(A)")
" Monte-Carlo photoionization uses adaptive spacing"
166 real(dp),
intent(in) :: p_o2
168 real(dp),
intent(in) :: p_h2o
170 procedure(absfunc_t) :: absfunc
173 real(dp) :: large_distance = 1e3_dp
175 integer,
parameter :: tbl_size = 500
176 real(dp),
allocatable :: fsum(:), dist(:)
178 allocate(fsum(tbl_size))
179 allocate(dist(tbl_size))
181 call integrate_absfunc(p_o2, p_h2o, absfunc, 1.0_dp, large_distance, &
182 tbl_size, fsum, dist, n)
188 fsum(1:n) = fsum(1:n) / fsum(n)
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))
193 write(*,
"(A,E12.3,A)")
" Average photon absorption length", &
194 1e3_dp * sum(dist(1:n-1))/(n-1),
" mm"
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
199 real(dp),
intent(in) :: p_h2o
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
209 df = f_max / (n_steps-1)
214 drdf = rk4_drdf(r(n-1), df, p_o2, p_h2o, absfunc)
215 r(n) = r(n-1) + df * drdf
217 if (fsum(n) > f_max .or. r(n) > r_max)
exit
223 end subroutine integrate_absfunc
227 real(dp) function rk4_drdf(r, df, p_o2, p_h2o, absfunc)
228 real(dp),
intent(in) :: r
229 real(dp),
intent(in) :: df
230 real(dp),
intent(in) :: p_o2
231 real(dp),
intent(in) :: p_h2o
232 procedure(absfunc_t) :: absfunc
235 real(dp),
parameter :: one_sixth = 1 / 6.0_dp
238 drdf = 1 / absfunc(r, p_o2, p_h2o)
242 drdf = 1 / absfunc(r + 0.5_dp * df * drdf, p_o2, p_h2o)
243 sum_drdf = sum_drdf + 2 * drdf
246 drdf = 1 / absfunc(r + 0.5_dp * df * drdf, p_o2, p_h2o)
247 sum_drdf = sum_drdf + 2 * drdf
250 drdf = 1 / absfunc(r + df * drdf, p_o2, p_h2o)
251 sum_drdf = sum_drdf + drdf
254 rk4_drdf = one_sixth * sum_drdf
255 end function rk4_drdf
261 real(dp) function absfunc_zheleznyak(dist, p_o2, p_h2o) result(f)
263 real(dp),
intent(in) :: dist
264 real(dp),
intent(in) :: p_o2
265 real(dp),
intent(in) :: p_h2o
269 real(dp),
parameter :: eps = epsilon(1.0_dp)
274 if (r * k1 < eps)
then
276 f = (k2 - k1 + 0.5_dp * (k1**2 - k2**2) * r) &
281 (exp(-k1 * r) - exp(-k2 * r)) / (dist * log(k2/k1)))
283 end function absfunc_zheleznyak
288 real(dp) function absfunc_naidis_humid(dist, p_o2, p_h2o) result(f)
289 real(dp),
intent(in) :: dist
290 real(dp),
intent(in) :: p_o2
291 real(dp),
intent(in) :: p_h2o
295 exp(-k3 * p_h2o * dist) * absfunc_zheleznyak(dist, p_o2, p_h2o))
296 end function absfunc_naidis_humid
302 real(dp) function absfunc_aints_humid(dist, p_o2, p_h2o) result(f)
303 real(dp),
intent(in) :: dist
304 real(dp),
intent(in) :: p_o2
305 real(dp),
intent(in) :: p_h2o
311 real(dp),
parameter :: eps = epsilon(1.0_dp)
313 k_a = k1 * p_o2 + k4 * p_h2o
314 k_b = k2 * p_o2 + k5 * p_h2o
316 if (dist * k_a < eps)
then
318 f = (k_b - k_a + 0.5_dp * (k_a**2 - k_b**2) * dist) / log(k_b/k_a)
322 (exp(-k_a * dist) - exp(-k_b * dist)) / (dist * log(k_b/k_a)))
324 end function absfunc_aints_humid
327 integer function get_lvl_length(dr_base, length)
328 real(dp),
intent(in) :: dr_base
329 real(dp),
intent(in) :: length
330 real(dp),
parameter :: invlog2 = 1 / log(2.0_dp)
333 ratio = dr_base / length
337 get_lvl_length = 1 + ceiling(log(ratio) * invlog2)
339 end function get_lvl_length
342 integer function get_rlvl_length(dr_base, length, rng)
344 real(dp),
intent(in) :: dr_base
345 real(dp),
intent(in) :: length
346 type(rng_t),
intent(inout) :: rng
347 real(dp),
parameter :: invlog2 = 1 / log(2.0_dp)
348 real(dp) :: ratio, tmp
350 ratio = dr_base / length
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
359 end function get_rlvl_length
369 integer,
intent(inout) :: n_photons
371 real(dp),
intent(in) :: xyz_in(3, n_photons)
373 real(dp),
intent(out) :: xyz_out(3, n_photons)
374 integer,
intent(in) :: n_dim
375 type(lt_t),
intent(in) :: tbl
376 type(prng_t),
intent(inout) :: prng
377 real(dp),
intent(in) :: k3
378 logical,
allocatable :: mask(:)
379 integer :: n, proc_id
380 real(dp) :: rr, dist, p_ionization
382 if (n_dim < 2 .or. n_dim > 3) error stop
"n_dim should be 2 or 3"
384 allocate(mask(n_photons))
387 proc_id = 1+omp_get_thread_num()
391 rr = prng%rngs(proc_id)%unif_01()
393 dist = lt_get_col(tbl, 1, rr)
395 xyz_out(:, n) = xyz_in(:, n) + prng%rngs(proc_id)%sphere(dist)
399 p_ionization = exp(-k3 * dist)
400 mask(n) = (prng%rngs(proc_id)%unif_01() < p_ionization)
411 n_photons = n_photons + 1
412 xyz_out(:, n_photons) = xyz_out(:, n)
470 type(af_t),
intent(inout) :: tree
471 type(rng_t),
intent(inout) :: rng
473 integer,
intent(in) :: i_src
475 integer,
intent(in) :: i_photo
477 logical,
intent(in) :: use_cyl
479 real(dp),
intent(in),
optional :: dt
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(:, :)
490 real(dp) :: tmp, r(3)
491 real(dp),
parameter :: pi = acos(-1.0_dp)
494 type(af_loc_t),
allocatable :: ph_loc(:)
495 real(dp),
parameter :: small_value = 1e-100_dp
497 if (ndim == 3 .and. use_cyl) error stop
"phmc_set_src: use_cyl is .true."
499 if (st_use_dielectric)
then
506 n_procs = omp_get_max_threads()
507 call prng%init_parallel(n_procs, rng)
510 call af_tree_sum_cc(tree, i_src, sum_production_rate)
516 if (
present(dt))
then
519 n_produced = dt * sum_production_rate / phmc_min_weight
521 if (n_produced < phmc_num_photons)
then
522 dt_fac = dt/phmc_min_weight
524 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
528 dt_fac = phmc_num_photons / (sum_production_rate + small_value)
532 max_num_photons = nint(1.2_dp * dt_fac * sum_production_rate + 1000)
534 call phmc_generate_photons(tree, dt_fac, i_src, xyz_src, &
535 max_num_photons, n_used, prng)
537 allocate(xyz_abs(3, n_used))
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)
552 if (st_use_dielectric)
then
555 2, n_used, 1.0_dp/dt_fac)
561 if (st_use_dielectric)
then
564 ndim, n_used, 1.0_dp/dt_fac)
568 call periodify_coordinates(n_used, xyz_abs)
570 allocate(ph_loc(n_used))
572 if (phmc_const_dx)
then
574 pi_lengthscale = lt_get_col(
phmc_tbl%tbl, 1, phmc_absorp_fac)
578 pho_lvl = get_lvl_length(maxval(tree%dr_base), pi_lengthscale)
582 ph_loc(n) = af_get_loc(tree, xyz_abs(1:ndim, n), pho_lvl)
587 proc_id = 1+omp_get_thread_num()
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)
600 call af_tree_clear_cc(tree, i_photo)
604 if (id > af_no_box)
then
607 dr = tree%boxes(id)%dr
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))
615 dr = tree%boxes(id)%dr
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))
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))
632 dr = tree%boxes(id)%dr
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))
642 if (phmc_const_dx)
then
650 do lvl = min_lvl, tree%highest_lvl-1
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])
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.)
667 call prng%update_seed(rng)
671 subroutine periodify_coordinates(n_points, coords)
672 integer,
intent(in) :: n_points
673 real(dp),
intent(inout) :: coords(3, n_points)
679 if (st_periodic(i_dim))
then
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))
687 end subroutine periodify_coordinates
792 subroutine phmc_generate_photons(tree, dt_fac, i_src, xyz_src, n_max, n_used, prng)
795 type(af_t),
intent(in) :: tree
797 real(dp),
intent(in) :: dt_fac
799 integer,
intent(in) :: i_src
801 real(dp),
allocatable,
intent(inout) :: xyz_src(:, :)
803 integer,
intent(in) :: n_max
805 integer,
intent(out) :: n_used
807 type(prng_t),
intent(inout) :: prng
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(:, :)
818 n_procs = omp_get_max_threads()
820 allocate(xyz_src(3, n_max))
821 allocate(photons_per_thread(n_procs))
822 allocate(photon_thread(n_max))
829 proc_id = 1+omp_get_thread_num()
832 do lvl = 1, tree%highest_lvl
833 dr = af_lvl_dr(tree, lvl)
835 do ix = 1,
size(tree%lvls(lvl)%leaves)
836 id = tree%lvls(lvl)%leaves(ix)
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)
846 tmp = dt_fac * tree%boxes(id)%cc(i, j, i_src) * product(dr)
849 tmp = dt_fac * tree%boxes(id)%cc(ijk, i_src) * product(dr)
852 n_create = floor(tmp)
854 if (prng%rngs(proc_id)%unif_01() < tmp - n_create) &
855 n_create = n_create + 1
857 if (n_create > 0)
then
860 n_used = n_used + n_create
862 my_num_photons = my_num_photons + n_create
866 r(1) = prng%rngs(proc_id)%unif_01()
867 r(2) = prng%rngs(proc_id)%unif_01()
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
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))
877 photon_thread(i_ph+n) = proc_id
885 photons_per_thread(proc_id) = my_num_photons
890 allocate(ix_offset(n_procs))
891 allocate(xyz_tmp(3, n_used))
895 ix_offset(n) = ix_offset(n-1) + photons_per_thread(n-1)
898 photons_per_thread = 0
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)
905 xyz_src(:, 1:n_used) = xyz_tmp(:, 1:n_used)
907 end subroutine phmc_generate_photons
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.
elemental integer function, public gas_index(name)
Find index of a gas component, return -1 if not found.
real(dp), dimension(:), allocatable, public, protected gas_fractions
real(dp), public, protected gas_pressure
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:
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.