3#include "../afivo/src/cpp_macros.h"
4 use iso_fortran_env,
only: error_unit
28 integer,
parameter :: int8 = selected_int_kind(18)
29 integer,
parameter :: max_attemps_per_time_step = 10
30 integer,
parameter :: datfile_version = 30
31 integer(int8) :: t_start, t_current, count_rate
32 real(dp) :: wc_time = 0.0_dp, inv_count_rate
33 real(dp) :: time_last_print, time_last_output
34 integer :: it_last_output
35 integer :: i, it, n, coord_type, box_bytes
36 integer :: n_steps_rejected
37 integer,
allocatable :: ref_links(:, :)
39 real(dp) :: time, dt, dt_lim, photoi_prev_time
40 real(dp) :: dt_gas_lim, dt_lim_step
41 real(dp) :: fraction_steps_rejected
42 real(dp) :: memory_limit_gb
43 type(af_t) :: tree_copy
44 type(ref_info_t) :: ref_info
45 integer :: output_cnt = 0
46 character(len=string_len) :: restart_from_file =
undefined_str
48 type(af_loc_t) :: loc_field, loc_field_t0
49 real(dp) :: pos_emax(ndim), pos_emax_t0(ndim)
50 real(dp) :: breakdown_field_td, current_output_dt
51 real(dp) :: time_until_next_pulse
52 real(dp) :: field_energy, field_energy_prev
53 real(dp) :: tmp, field_energy_prev_time
54 logical :: step_accepted, start_of_new_pulse
57 real(dp) :: t1, t2, t3
64 call print_program_name()
67 call cfg_update_from_arguments(cfg)
69 call cfg_add_get(cfg,
"restart_from_file", restart_from_file, &
70 "If set, restart simulation from a previous .dat file")
72 memory_limit_gb = 4.0_dp**(ndim-1)
73 call cfg_add_get(cfg,
"memory_limit_GB", memory_limit_gb, &
78 call cfg_write(cfg, trim(
output_name) //
"_out.cfg", custom_first=.true.)
81 write(*,
'(A,E12.4)')
" Estimated breakdown field (Td): ", breakdown_field_td
97 funcval=set_gas_density_from_user_function)
101 do i = 1, tree%n_var_cell
102 if (tree%cc_write_output(i) .and. .not. &
103 (tree%has_cc_method(i) .or. i ==
i_phi))
then
104 call af_set_cc_methods(tree, i, af_bc_neumann_zero, &
113 photoi_prev_time = time
116 fraction_steps_rejected = 0.0_dp
122 call af_read_tree(tree, restart_from_file, read_sim_data)
125 tree%cc_write_output = tree_copy%cc_write_output
126 tree%cc_write_binary = tree_copy%cc_write_binary
127 tree%fc_write_binary = tree_copy%fc_write_binary
129 box_bytes = af_box_bytes(tree%n_cell, tree%n_var_cell, tree%n_var_face)
130 tree%box_limit = nint(memory_limit_gb * 2.0_dp**30 / box_bytes)
133 error stop
"restart_from_file: incompatible box size"
135 if (tree%n_var_cell /= tree_copy%n_var_cell)
then
136 write(error_unit, *)
"n_var_cell here:", tree_copy%n_var_cell
137 write(error_unit, *)
"n_var_cell file:", tree%n_var_cell
138 error stop
"restart_from_file: incompatible variable list"
146 call mg_init(tree,
mg)
159 call set_initial_conditions(tree,
mg)
163 call output_write(tree, output_cnt, 0.0_dp, write_sim_data)
167 print *,
"Number of threads: ", af_get_max_threads()
168 call af_print_info(tree)
169 call af_stencil_print_info(tree)
172 call system_clock(t_start, count_rate)
173 inv_count_rate = 1.0_dp / count_rate
174 time_last_print = -1e10_dp
175 time_last_output = time
179 field_energy_prev_time = time
191 call af_tree_max_cc(tree,
i_electric_fld, max_field, loc_field_t0)
192 pos_emax_t0 = af_r_loc(tree, loc_field_t0)
198 pos_emax = af_r_loc(tree, loc_field)
201 print *,
"Streamer reached its desired length"
207 call system_clock(t_current)
208 wc_time = (t_current - t_start) * inv_count_rate
213 time_last_print = wc_time
228 write_out = (time + dt >= time_last_output + current_output_dt .or. &
233 dt = min(dt, max(0.0_dp, time_last_output + current_output_dt - time))
237 start_of_new_pulse = (dt >= time_until_next_pulse)
238 if (start_of_new_pulse)
then
239 dt = max(time_until_next_pulse,
dt_min)
247 photoi_prev_time = time
251 call set_electrode_densities(tree)
256 step_accepted = .false.
257 do n = 1, max_attemps_per_time_step
259 call copy_current_state()
263 call af_advance(tree, dt, dt_lim_step, time,
all_densities, &
267 dt_lim = min(dt_lim, dt_lim_step)
270 step_accepted = (dt <= dt_lim_step)
272 if (step_accepted)
then
275 n_steps_rejected = n_steps_rejected + 1
276 write (*,
"(I0,A,I0,A,2E12.4,A,F6.2)") it,
" Step rejected (#", &
277 n_steps_rejected,
") (dt, dt_lim) = ", dt, dt_lim, &
278 " frac", fraction_steps_rejected
285 call restore_previous_state()
290 fraction_steps_rejected = 0.99_dp * fraction_steps_rejected
291 if (n > 1) fraction_steps_rejected = fraction_steps_rejected + 0.01_dp
293 if (n == max_attemps_per_time_step+1) &
294 error stop
"All time steps were rejected"
311 tmp = (field_energy - field_energy_prev)/(time - field_energy_prev_time)
312 field_energy_prev = field_energy
313 field_energy_prev_time = time
337 call af_advance(tree, dt, dt_gas_lim, time, &
351 if (fraction_steps_rejected > 0.1_dp) tmp = 1.0_dp
355 if (start_of_new_pulse)
then
368 write(error_unit,
"(A,E12.4,A)")
" Time step (dt =",
global_dt, &
369 ") getting too small"
370 write(error_unit, *)
"See the documentation on time integration"
377 output_cnt = output_cnt + 1
380 call output_write(tree, output_cnt, wc_time, write_sim_data)
383 [
"photon_flux",
"surf_dens "],
output_name, output_cnt)
392 if (axisymmetric_is_branching(tree))
then
393 error stop
"Branching detected, abort"
403 call af_restrict_tree(tree,
gas_vars)
409 call surface_get_refinement_links(
diel, ref_links)
412 call surface_update_after_refinement(tree,
diel, ref_info)
418 if (ref_info%n_add > 0 .or. ref_info%n_rm > 0)
then
425 photoi_prev_time = time
436 write(*,
"(A)")
"Computational cost breakdown (%)"
437 write(*,
"(7(A10))")
"flux",
"source",
"copy",
"field",
"output", &
451 type(cfg_t),
intent(inout) :: cfg
452 type(af_t),
intent(inout) :: tree
453 type(mg_t),
intent(inout) :: mg
454 logical,
intent(in) :: restart
476 error stop
"Stochastic species not supported in 1D/2D Cartesian"
480 subroutine set_initial_conditions(tree, mg)
481 type(af_t),
intent(inout) :: tree
482 type(mg_t),
intent(inout) :: mg
483 integer :: n, lvl, i, id, n_surface_variables
486 do lvl = 1, af_max_lvl-1
489 call af_refine_up_to_lvl(tree, lvl)
496 error stop
"use_dielectric requires user_initial_conditions to be set"
501 call mg_init(tree, mg)
507 call surface_initialize(tree,
i_eps,
diel, n_surface_variables)
515 call surface_get_refinement_links(
diel, ref_links)
518 call surface_update_after_refinement(tree,
diel, ref_info)
525 do lvl = 1,
size(ref_info%lvls)
527 do i = 1,
size(ref_info%lvls(lvl)%add)
528 id = ref_info%lvls(lvl)%add(i)
536 if (ref_info%n_add == 0)
exit
539 end subroutine set_initial_conditions
541 subroutine write_sim_data(my_unit)
542 integer,
intent(in) :: my_unit
544 write(my_unit) datfile_version
547 write(my_unit) output_cnt
550 write(my_unit) photoi_prev_time
555 write(my_unit) fraction_steps_rejected
556 end subroutine write_sim_data
558 subroutine read_sim_data(my_unit)
559 integer,
intent(in) :: my_unit
562 read(my_unit) version
563 if (version /= datfile_version) error stop
"Different datfile version"
566 read(my_unit) output_cnt
569 read(my_unit) photoi_prev_time
574 read(my_unit) fraction_steps_rejected
577 end subroutine read_sim_data
579 subroutine print_program_name()
580 print *,
" __ _ _____ _ "
581 print *,
" /\ / _(_) / ____| | "
582 print *,
" / \ | |_ ___ _____ _____| (___ | |_ _ __ ___ __ _ _ __ ___ ___ _ __ "
583 print *,
" / /\ \ | _| \ \ / / _ \______\___ \| __| '__/ _ \/ _` | '_ ` _ \ / _ \ '__|"
584 print *,
" / ____ \| | | |\ V / (_) | ____) | |_| | | __/ (_| | | | | | | __/ | "
585 print *,
" /_/ \_\_| |_| \_/ \___/ |_____/ \__|_| \___|\__,_|_| |_| |_|\___|_| "
587 end subroutine print_program_name
590 subroutine set_electrode_densities(tree)
591 type(af_t),
intent(inout) :: tree
593 call af_loop_box(tree, electrode_species_bc)
594 end subroutine set_electrode_densities
598 subroutine electrode_species_bc(box)
599 type(box_t),
intent(inout) :: box
600 real(dp) :: lsf_nb(2*NDIM), dens_nb(2*NDIM)
603 if (iand(box%tag, mg_lsf_box) == 0)
return
607 if (box%cc(ijk,
i_lsf) < 0)
then
612 lsf_nb = [box%cc(i-1,
i_lsf), &
615 lsf_nb = [box%cc(i-1, j,
i_lsf), &
616 box%cc(i+1, j,
i_lsf), &
617 box%cc(i, j-1,
i_lsf), &
618 box%cc(i, j+1,
i_lsf)]
620 lsf_nb = [box%cc(i-1, j, k,
i_lsf), &
621 box%cc(i+1, j, k,
i_lsf), &
622 box%cc(i, j-1, k,
i_lsf), &
623 box%cc(i, j+1, k,
i_lsf), &
624 box%cc(i, j, k-1,
i_lsf), &
625 box%cc(i, j, k+1,
i_lsf)]
628 if (any(lsf_nb > 0) .and. &
629 associated(
bc_species, af_bc_neumann_zero))
then
649 box%cc(ijk,
i_electron) = sum(dens_nb, mask=(lsf_nb > 0)) / &
656 end subroutine electrode_species_bc
659 subroutine copy_current_state()
662 n_states = af_advance_num_steps(time_integrator)
663 call af_tree_copy_ccs(tree, all_densities, all_densities+n_states)
666 call af_tree_copy_cc(tree, i_phi, i_phi+1)
668 if (st_use_dielectric)
then
669 call surface_copy_variable(diel, i_surf_dens, i_surf_dens+n_states)
671 end subroutine copy_current_state
674 subroutine restore_previous_state()
677 n_states = af_advance_num_steps(time_integrator)
679 call af_tree_copy_ccs(tree, all_densities+n_states, all_densities)
682 call af_tree_copy_cc(tree, i_phi+1, i_phi)
683 call field_from_potential(tree, mg)
685 if (st_use_dielectric)
then
686 call surface_copy_variable(diel, i_surf_dens+n_states, i_surf_dens )
688 end subroutine restore_previous_state
692 subroutine set_gas_density_from_user_function(box, iv)
693 type(box_t),
intent(inout) :: box
694 integer,
intent(in) :: iv
699 box%cc(ijk, iv) = user_gas_density(box, ijk)
701 end subroutine set_gas_density_from_user_function
703 logical function axisymmetric_is_branching(tree)
704 type(af_t),
intent(in) :: tree
705 real(dp) :: max_field, axis_field, rmin(NDIM), rmax(NDIM)
706 type(af_loc_t) :: loc_field
709 call af_tree_max_cc(tree, i_electric_fld, max_field, loc_field)
713 rmax(1) = 1e-2_dp * st_domain_len(1)
714 rmax(2:) = st_domain_len(2:)
715 call analysis_max_var_region(tree, i_electric_fld, rmin, &
716 rmax, axis_field, loc_field)
718 axisymmetric_is_branching = (max_field > 1.1_dp * axis_field)
720 if (axisymmetric_is_branching) print *, max_field, axis_field
721 end function axisymmetric_is_branching
Module with routines to help analyze simulations.
Module for handling chemical reactions.
subroutine, public chemistry_initialize(tree, cfg)
Initialize module and load chemical reactions.
integer, public, protected n_stochastic_species
Number of stochastic species.
subroutine, public chemistry_convert_stochastic_species(tree, rng)
Convert stochastic species into regular species, using random numbers.
subroutine, public chemistry_get_breakdown_field(field_td, min_growth_rate)
Get the breakdown field in Townsend.
integer, public, protected n_reactions
Number of reactions present.
Module for the coupling between the gas dynamics and the fluid model.
subroutine, public coupling_add_fluid_source(tree, dt)
Add source terms form the fluid model to the Euler equations.
subroutine, public coupling_update_gas_density(tree)
Update gas number density in the fluid model.
Module with settings and routines to handle dielectrics.
integer, parameter, public i_surf_dens
logical, public, protected surface_output
subroutine, public dielectric_initialize(tree, cfg)
type(surfaces_t), public diel
To store dielectric surface.
integer, parameter, public i_photon_flux
Module to set the time step.
integer, public, protected time_integrator
Which time integrator is used.
real(dp), public, protected dt_max
real(dp), public, protected dt_safety_factor
subroutine, public dt_initialize(cfg)
Initialize the time step module.
real(dp), public, protected dt_min
real(dp), public, protected dt_max_growth_factor
Maximal relative increase dt for the next iteration.
Module to compute electric fields.
real(dp), public, protected field_pulse_period
Time of one complete voltage pulse.
subroutine, public field_compute_energy(tree, field_energy)
Compute total field energy in Joule, defined as the volume integral over 1/2 * epsilon * E^2.
real(dp), public, protected current_voltage
The current applied voltage.
subroutine, public field_initialize(tree, cfg, mg)
Initialize this module.
subroutine, public field_compute(tree, mg, s_in, time, have_guess)
Compute electric field on the tree. First perform multigrid to get electric potential,...
subroutine, public forward_euler(tree, dt, dt_stiff, dt_lim, time, s_deriv, n_prev, s_prev, w_prev, s_out, i_step, n_steps)
Advance fluid model using forward Euler step. If the equation is written as y' = f(y),...
Module that stores parameters related to the gas.
subroutine, public gas_forward_euler(tree, dt, dt_stiff, dt_lim, time, s_deriv, n_prev, s_prev, w_prev, s_out, i_step, n_steps)
A forward-Euler method for the Euler equations.
logical, public, protected gas_dynamics
Whether the gas dynamics are simulated.
integer, public, protected i_gas_dens
integer, dimension(n_vars_euler), public, protected gas_vars
logical, public, protected gas_constant_density
Whether the gas has a constant density.
subroutine, public gas_initialize(tree, cfg)
Initialize this module.
Module to help setting up initial conditions.
subroutine, public init_cond_set_box(box)
Sets the initial condition.
subroutine, public init_cond_initialize(tree, cfg)
Set the initial conditions from the configuration.
Module to set the type of model.
subroutine, public model_initialize(cfg)
Initialize the module.
Module with methods and parameters related to writing output.
integer, public, protected output_dt_factor_pulse_off
subroutine, public output_status(tree, time, wc_time, it, dt)
Print short status message.
real(dp), public, protected output_dt
subroutine, public output_initial_summary(tree)
Write a summary of the model and parameters used.
character(len=string_len), public, protected output_name
real(dp), public, protected output_status_delay
integer, public, protected output_dit
subroutine, public output_initialize(tree, cfg)
subroutine, public output_write(tree, output_cnt, wc_time, write_sim_data)
Top-module for photoionization, which can make use of different methods.
subroutine, public photoi_set_src(tree, dt)
Sets the photoionization.
logical, public, protected photoi_enabled
Whether photoionization is enabled.
subroutine, public photoi_initialize(tree, cfg)
Initialize photoionization parameters.
integer, public, protected photoi_per_steps
Update photoionization every N time step.
Module with grid refinement settings and routines.
integer, public, protected refine_buffer_width
subroutine, public refine_initialize(cfg)
Initialize the grid refinement options.
real(dp), public current_electrode_dx
real(dp), public refine_electrode_dx
real(dp), public, protected refine_prepulse_time
integer, public, protected refine_per_steps
procedure(af_subr_ref), pointer, public refine_routine
real(dp), public, protected electrode_derefine_factor
real(dp), public, protected refine_max_dx
This module contains several pre-defined variables like:
logical, public, protected st_cylindrical
Whether cylindrical coordinates are used.
real(dp), public wc_time_refine
procedure(af_subr_prolong), pointer, public, protected st_prolongation_method
Method used to prolong (interpolate) densities.
real(dp), public wc_time_output
integer, public, protected i_eps
Index can be set to include a dielectric.
real(dp), dimension(:), allocatable, public st_global_rates
Global sum of reaction rates.
integer, public, protected current_update_per_steps
Per how many iterations the electric current is computed.
integer, public, protected st_initial_streamer_pos_steps_wait
Wait n steps before initializing streamer begin position.
real(dp), public wc_time_field
real(dp), public st_global_jdote_current
Electric current through electrodes due to J.E.
real(dp), public global_dt
Global time step.
type(mg_t), public mg
Multigrid option structure.
integer, dimension(ndim), public, protected st_coarse_grid_size
Size of the coarse grid.
subroutine, public st_initialize(tree, cfg, ndim)
Create the configuration file with default values.
integer, public, protected st_box_size
The size of the boxes that we use to construct our mesh.
real(dp), dimension(:, :), allocatable, public st_current_jdote
Current sum of J.E per thread.
logical, dimension(ndim), public, protected st_periodic
Whether the domain is periodic (per dimension)
real(dp), public global_time
Global time.
procedure(af_subr_bc), pointer, public, protected bc_species
Boundary condition for the plasma species.
logical, public, protected st_use_electrode
Whether to include an electrode.
real(dp), public wc_time_source
integer, public, protected i_lsf
Index can be set to include an electrode.
logical, public, protected st_use_dielectric
Whether a dielectric is used.
integer, public, protected i_electron
Index of electron density.
integer, public, protected i_1pos_ion
Index of first positive ion species.
logical, public, protected st_abort_axisymmetric_if_branching
Abort axisymmetric simulations if there is branching.
integer, public, protected i_phi
Index of electrical potential.
real(dp), public wc_time_flux
real(dp), public st_global_jdote
Global sum of J.E.
integer, public, protected i_electric_fld
Index of electric field norm.
real(dp), public wc_time_copy_state
real(dp), public, protected st_end_streamer_length
Streamer length at which the simulation will stop.
real(dp), dimension(ndim), public, protected st_domain_len
Domain length per dimension.
real(dp), public, protected st_end_time
End time of the simulation.
real(dp), dimension(ndim), public, protected st_domain_origin
Origin of domain.
type(rng_t), public st_rng
Random number generator.
integer, dimension(:), allocatable, public, protected all_densities
Index of all densities that evolve in time.
real(dp), public wc_time_photoi
logical, public, protected st_use_end_streamer_length
Whether streamer length is used as a simulation stopping.
real(dp), dimension(:, :), allocatable, public st_current_rates
Current sum of reaction rates per thread.
real(dp), public st_global_displ_current
Electric current through electrodes due to displacement current.
Module with settings and routines for tabulated data.
subroutine, public table_data_initialize(cfg)
Initialize this module.
Module that provides routines for reading in arbritrary transport data.
subroutine, public transport_data_initialize(cfg)
Initialize the transport coefficients.
real(dp), parameter huge_real
Huge number.
character(len= *), parameter undefined_str
Undefined string.
Module that contains physical and numerical constants.
This module contains all the methods that users can customize.
procedure(generic_subr), pointer user_generic_method
Generic procedure that is called every time step.
procedure(af_subr), pointer user_initial_conditions
If defined, call this routine after setting initial conditions.
procedure(af_subr), pointer user_new_pulse_conditions
If defined, call this routine after a new voltage pulse starts.
Template for user code, this one simply uses the default routines.
subroutine, public user_initialize(cfg, tree)
program streamer
Program to perform streamer simulations with AMR.
subroutine initialize_modules(cfg, tree, mg, restart)