4 program simple_streamer_2d
19 character(len=100) :: fname
22 type(ref_info_t) :: refine_info
25 integer,
parameter :: n_var_cell = 7
26 integer,
parameter :: i_elec = 1
27 integer,
parameter :: i_pion = 2
28 integer,
parameter :: i_elec_old = 3
29 integer,
parameter :: i_pion_old = 4
30 integer,
parameter :: i_phi = 5
31 integer,
parameter :: i_fld = 6
32 integer,
parameter :: i_rhs = 7
35 integer,
parameter :: n_var_face = 2
36 integer,
parameter :: f_elec = 1
37 integer,
parameter :: f_fld = 2
40 character(len=10) :: st_cc_names(n_var_cell) = &
41 [character(len=10) ::
"elec",
"pion",
"elec_old", &
42 "pion_old",
"phi",
"fld",
"rhs"]
45 real(dp),
parameter :: end_time = 8e-9_dp
46 real(dp),
parameter :: dt_output = 20e-11_dp
47 real(dp),
parameter :: dt_max = 20e-11_dp
48 integer,
parameter :: ref_per_steps = 2
49 integer,
parameter :: box_size = 8
52 real(dp),
parameter :: applied_field = -0.8e7_dp
53 real(dp),
parameter :: mobility = 0.03_dp
54 real(dp),
parameter :: diffusion_c = 0.2_dp
57 real(dp),
parameter :: domain_length = 10e-3_dp
58 real(dp),
parameter :: refine_max_dx = 1e-3_dp
59 real(dp),
parameter :: refine_min_dx = 1e-9_dp
62 real(dp),
parameter :: init_density = 1e15_dp
63 real(dp),
parameter :: init_y_min = 8.0e-3_dp
64 real(dp),
parameter :: init_y_max = 9.0e-3_dp
69 integer :: output_count
72 real(dp) :: sum_elec, sum_pion
83 mg%sides_bc => sides_bc_pot
88 call a2_set_cc_methods(tree, i_elec, a2_bc_dirichlet_zero, &
89 a2_gc_interp_lim, a2_prolong_limit)
90 call a2_set_cc_methods(tree, i_pion, a2_bc_dirichlet_zero, &
91 a2_gc_interp_lim, a2_prolong_limit)
92 call a2_set_cc_methods(tree, i_phi, mg%sides_bc, mg%sides_rb)
100 call a2_loop_box(tree, set_initial_condition)
105 call compute_fld(tree, .false.)
112 call a2_adjust_refinement(tree, &
113 refinement_routine, &
117 if (refine_info%n_add == 0 .and. refine_info%n_rm == 0)
exit 120 call a2_print_info(tree)
124 call a2_reduction(tree, &
131 print *,
"dt getting too small, instability?", dt
132 time = end_time + 1.0_dp
136 if (output_count * dt_output <= time)
then 137 output_count = output_count + 1
138 write(fname,
"(A,I6.6)")
"simple_streamer_2d_", output_count
142 call a2_write_silo(tree, fname, output_count, time, dir=
"output")
144 call a2_tree_sum_cc(tree, i_elec, sum_elec)
145 call a2_tree_sum_cc(tree, i_pion, sum_pion)
146 print *,
"sum(pion-elec)/sum(pion)", (sum_pion - sum_elec)/sum_pion
149 if (time > end_time)
exit 152 do n = 1, ref_per_steps
156 call a2_tree_copy_cc(tree, i_elec, i_elec_old)
157 call a2_tree_copy_cc(tree, i_pion, i_pion_old)
162 call a2_loop_boxes(tree, fluxes_koren, leaves_only=.true.)
163 call a2_consistent_fluxes(tree, [f_elec])
166 call a2_loop_box_arg(tree, update_solution, [dt], &
170 if (i == 1)
call compute_fld(tree, .true.)
174 call a2_loop_box(tree, average_dens)
177 call compute_fld(tree, .true.)
181 call a2_gc_tree(tree, i_elec)
182 call a2_gc_tree(tree, i_pion)
190 call a2_adjust_refinement(tree, &
191 refinement_routine, &
195 if (refine_info%n_add > 0 .or. refine_info%n_rm > 0)
then 197 call compute_fld(tree, .true.)
204 subroutine init_tree(tree)
205 type(
a2_t),
intent(inout) :: tree
210 integer :: ix_list(2, 1)
211 integer :: nb_list(4, 1)
212 integer :: n_boxes_init = 1000
214 dr = domain_length / box_size
223 n_boxes = n_boxes_init, &
224 cc_names=st_cc_names)
228 ix_list(:, id) = [1,1]
230 nb_list(a2_neighb_lowy, id) = -1
231 nb_list(a2_neighb_highy, id) = -1
232 nb_list(a2_neighb_lowx, id) = id
233 nb_list(a2_neighb_highx, id) = id
236 call a2_set_base(tree, 1, ix_list, nb_list)
238 end subroutine init_tree
241 subroutine refinement_routine(box, cell_flags)
242 type(
box2_t),
intent(in) :: box
243 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
245 real(dp) :: dx, dens, fld, adx, xy(2)
252 xy = a2_r_cc(box, [i,j])
253 dens = box%cc(i, j, i_elec)
254 fld = box%cc(i, j, i_fld)
255 adx = get_alpha(fld) * dx
257 if (dens > 1e0_dp .and. adx > 0.8_dp)
then 258 cell_flags(i, j) = af_do_ref
259 else if (dx < 1.25e-5_dp .and. adx < 0.1_dp)
then 260 cell_flags(i, j) = af_rm_ref
262 cell_flags(i, j) = af_keep_ref
268 if (dx > refine_max_dx)
then 269 cell_flags = af_do_ref
270 else if (dx < 2 * refine_min_dx)
then 271 where(cell_flags == af_do_ref) cell_flags = af_keep_ref
272 else if (dx > 0.5_dp * refine_max_dx)
then 273 where(cell_flags == af_rm_ref) cell_flags = af_keep_ref
276 end subroutine refinement_routine
279 subroutine set_initial_condition(box)
280 type(
box2_t),
intent(inout) :: box
282 real(dp) :: xy(2), normal_rands(2)
288 xy = a2_r_cc(box, [i,j])
290 if (xy(2) > init_y_min .and. xy(2) < init_y_max)
then 292 normal_rands = two_normals(box%dr**3 * init_density, &
293 sqrt(box%dr**3 * init_density))
295 box%cc(i, j, i_elec) = abs(normal_rands(1)) / box%dr**3
297 box%cc(i, j, i_elec) = 0
302 box%cc(:, :, i_pion) = box%cc(:, :, i_elec)
303 box%cc(:, :, i_phi) = 0
305 end subroutine set_initial_condition
309 function two_normals(mean, sigma)
result(rands)
310 real(dp),
intent(in) :: mean, sigma
311 real(dp) :: rands(2), sum_sq
314 call random_number(rands)
315 rands = 2 * rands - 1
316 sum_sq = sum(rands**2)
317 if (sum_sq < 1.0_dp .and. sum_sq > 0.0_dp)
exit 319 rands = rands * sqrt(-2 * log(sum_sq) / sum_sq)
320 rands = mean + rands * sigma
321 end function two_normals
324 real(dp) function get_min(a, b)
325 real(dp),
intent(in) :: a, b
330 real(dp) function max_dt(box)
331 type(
box2_t),
intent(in) :: box
332 real(dp),
parameter :: uc_eps0 = 8.8541878176d-12
333 real(dp),
parameter :: uc_elem_charge = 1.6022d-19
336 real(dp) :: dt_cfl, dt_dif, dt_drt
343 fld(1) = 0.5_dp * (box%fc(i, j, 1, f_fld) + &
344 box%fc(i+1, j, 1, f_fld))
345 fld(2) = 0.5_dp * (box%fc(i, j, 2, f_fld) + &
346 box%fc(i, j+1, 2, f_fld))
349 dt_cfl = min(dt_cfl, 0.5_dp / sum(abs(fld * mobility) / box%dr))
354 dt_drt = uc_eps0 / (uc_elem_charge * mobility * &
355 maxval(box%cc(1:nc, 1:nc, i_elec)) + epsilon(1.0_dp))
358 dt_dif = 0.25_dp * box%dr**2 / max(diffusion_c, epsilon(1.0_dp))
360 max_dt = min(dt_cfl, dt_drt, dt_dif)
369 elemental function get_alpha(fld)
result(alpha)
370 real(dp),
intent(in) :: fld
372 real(dp),
parameter :: c0 = 1.04e1_dp
373 real(dp),
parameter :: c1 = 0.601_dp
374 real(dp),
parameter :: c2 = 1.86e7_dp
376 alpha = exp(c0) * (abs(fld) * 1e-5_dp)**c1 * exp(-c2 / abs(fld))
377 end function get_alpha
381 subroutine compute_fld(tree, have_guess)
382 type(
a2_t),
intent(inout) :: tree
383 logical,
intent(in) :: have_guess
384 real(dp),
parameter :: fac = 1.6021766208e-19_dp / 8.8541878176e-12_dp
385 integer :: lvl, i, id, nc
391 do lvl = 1, tree%highest_lvl
393 do i = 1,
size(tree%lvls(lvl)%leaves)
394 id = tree%lvls(lvl)%leaves(i)
395 tree%boxes(id)%cc(:, :, i_rhs) = fac * (&
396 tree%boxes(id)%cc(:, :, i_elec) - &
397 tree%boxes(id)%cc(:, :, i_pion))
404 call mg2_fas_fmg(tree, mg, .false., have_guess)
407 call a2_loop_box(tree, fld_from_pot)
410 call a2_gc_tree(tree, i_fld, a2_gc_interp, a2_bc_neumann_zero)
411 end subroutine compute_fld
414 subroutine fld_from_pot(box)
415 type(
box2_t),
intent(inout) :: box
422 box%fc(1:nc+1, 1:nc, 1, f_fld) = inv_dr * &
423 (box%cc(0:nc, 1:nc, i_phi) - box%cc(1:nc+1, 1:nc, i_phi))
424 box%fc(1:nc, 1:nc+1, 2, f_fld) = inv_dr * &
425 (box%cc(1:nc, 0:nc, i_phi) - box%cc(1:nc, 1:nc+1, i_phi))
427 box%cc(1:nc, 1:nc, i_fld) = sqrt(&
428 0.25_dp * (box%fc(1:nc, 1:nc, 1, f_fld) + &
429 box%fc(2:nc+1, 1:nc, 1, f_fld))**2 + &
430 0.25_dp * (box%fc(1:nc, 1:nc, 2, f_fld) + &
431 box%fc(1:nc, 2:nc+1, 2, f_fld))**2)
432 end subroutine fld_from_pot
435 subroutine fluxes_koren(boxes, id)
437 type(
box2_t),
intent(inout) :: boxes(:)
438 integer,
intent(in) :: id
440 real(dp) :: cc(-1:boxes(id)%n_cell+2, -1:boxes(id)%n_cell+2)
441 real(dp),
allocatable :: v(:, :, :), dc(:, :, :)
444 nc = boxes(id)%n_cell
445 inv_dr = 1/boxes(id)%dr
447 allocate(v(1:nc+1, 1:nc+1, 2))
448 allocate(dc(1:nc+1, 1:nc+1, 2))
450 call a2_gc_box(boxes, id, i_elec, a2_gc_interp_lim, &
451 a2_bc_dirichlet_zero)
452 call a2_gc2_box(boxes, id, i_elec, a2_gc2_prolong_linear, &
453 a2_bc2_dirichlet_zero, cc, nc)
455 v = -mobility * boxes(id)%fc(:, :, :, f_fld)
458 call flux_koren_2d(cc, v, nc, 2)
459 call flux_diff_2d(cc, dc, inv_dr, nc, 2)
461 boxes(id)%fc(:, :, :, f_elec) = v + dc
462 end subroutine fluxes_koren
465 subroutine average_dens(box)
466 type(
box2_t),
intent(inout) :: box
467 box%cc(:, :, i_elec) = 0.5_dp * (box%cc(:, :, i_elec) + box%cc(:, :, i_elec_old))
468 box%cc(:, :, i_pion) = 0.5_dp * (box%cc(:, :, i_pion) + box%cc(:, :, i_pion_old))
469 end subroutine average_dens
472 subroutine update_solution(box, dt)
473 type(
box2_t),
intent(inout) :: box
474 real(dp),
intent(in) :: dt(:)
475 real(dp) :: inv_dr, src, sflux, fld
483 fld = box%cc(i,j, i_fld)
484 alpha = get_alpha(fld)
485 src = box%cc(i, j, i_elec) * mobility * abs(fld) * alpha
487 sflux = (sum(box%fc(i, j, :, f_elec)) - &
488 box%fc(i+1, j, 1, f_elec) - &
489 box%fc(i, j+1, 2, f_elec)) * inv_dr
491 box%cc(i, j, i_elec) = box%cc(i, j, i_elec) + (src + sflux) * dt(1)
492 box%cc(i, j, i_pion) = box%cc(i, j, i_pion) + src * dt(1)
496 end subroutine update_solution
499 subroutine sides_bc_pot(box, nb, iv, bc_type)
500 type(
box2_t),
intent(inout) :: box
501 integer,
intent(in) :: nb
502 integer,
intent(in) :: iv
503 integer,
intent(out) :: bc_type
509 case (a2_neighb_lowy)
510 bc_type = af_bc_neumann
511 box%cc(1:nc, 0, iv) = applied_field
512 case (a2_neighb_highy)
513 bc_type = af_bc_dirichlet
514 box%cc(1:nc, nc+1, iv) = 0
516 stop
"sides_bc_pot: unspecified boundary" 518 end subroutine sides_bc_pot
520 end program simple_streamer_2d