Afivo  0.3
space_plasma_2d.f90
1 
2 program space_plasma_2d
3  use m_a2_all
4 
5  implicit none
6 
7  ! Cell-centered variables
8  integer, parameter :: i_cs_ion = 1
9  integer, parameter :: i_cs = 2
10  integer, parameter :: i_o = 3
11  integer, parameter :: i_o_min = 4
12  integer, parameter :: i_cs_ion_old = 5
13  integer, parameter :: i_cs_old = 6
14  integer, parameter :: i_o_old = 7
15  integer, parameter :: i_o_min_old = 8
16 
17  ! Domain
18  integer, parameter :: box_size = 8
19  real(dp), parameter :: domain_len = 3000.0_dp
20  real(dp), parameter :: dr = domain_len / box_size
21  integer, parameter :: min_refinement_lvl = 4
22  integer, parameter :: max_refinement_lvl = 7
23 
24  ! Time steps
25  real(dp) :: dt_adapt = 60.0_dp
26  real(dp) :: dt_output = 60.0_dp
27  real(dp) :: end_time = 3600.0_dp
28 
29  ! Physical parameters
30  real(dp), parameter :: dc_cs_ion_parallel = 10.0_dp
31  real(dp), parameter :: dc_cs_ion_perpendicular = 2.0_dp
32  real(dp), parameter :: dc_cs = 20.0_dp
33  real(dp), parameter :: k_attachment = 1.0e-3_dp
34  real(dp), parameter :: k_recombination = 1.0e-3_dp
35  real(dp), parameter :: k_photoionization = 5.0e-4_dp
36  real(dp), parameter :: initial_cs_density = 1.0_dp
37  real(dp), parameter :: initial_ionization_fraction = 0.5_dp
38  real(dp), parameter :: initial_o_density = 1e-7
39  real(dp), parameter :: initial_sigma = 50.0_dp
40 
41  real(dp) :: max_diffusion_coeff
42  type(a2_t) :: tree
43  type(ref_info_t) :: refine_info
44  integer :: ix_list(2, 1)
45  integer :: time_steps, output_cnt
46  integer :: i, n, n_steps
47  real(dp) :: dt, time
48  real(dp) :: dr_min(2)
49  character(len=100) :: fname
50 
51  print *, "Running space_plasma_2d"
52  print *, "Number of threads", af_get_max_threads()
53 
54  ! Initialize tree
55  call a2_init(tree, & ! Tree to initialize
56  box_size, & ! A box contains box_size**DIM cells
57  n_var_cell=8, & ! Number of cell-centered variables
58  n_var_face=2, & ! Number of face-centered variables
59  dr=dr, & ! Distance between cells on base level
60  r_min=[-0.5_dp, -0.5_dp] * domain_len, &
61  cc_names=[character(len=10) :: "Cs_ion", "Cs", "O", "O_min", &
62  "Cs_ion_old", "Cs_old", "O_old", "O_min_old"]) ! Variable names
63 
64  ! Set default methods
65  call a2_set_cc_methods(tree, i_cs_ion, a2_bc_neumann_zero, a2_gc_interp_lim)
66  call a2_set_cc_methods(tree, i_cs, a2_bc_neumann_zero, a2_gc_interp_lim)
67  call a2_set_cc_methods(tree, i_o, a2_bc_neumann_zero, a2_gc_interp_lim)
68  call a2_set_cc_methods(tree, i_o_min, a2_bc_neumann_zero, a2_gc_interp_lim)
69 
70  ! Set index of box
71  ix_list(:, 1) = [1, 1]
72 
73  ! Create the base mesh
74  call a2_set_base(tree, 1, ix_list)
75 
76  output_cnt = 0
77  time = 0
78 
79  do
80  call a2_loop_box(tree, set_initial_condition)
81  call a2_adjust_refinement(tree, & ! tree
82  refine_routine, & ! Refinement function
83  refine_info, & ! Information about refinement
84  1) ! Buffer width (in cells)
85  if (refine_info%n_add == 0) exit
86  end do
87 
88  call a2_gc_tree(tree, i_cs_ion, a2_gc_interp_lim, a2_bc_neumann_zero)
89  call a2_gc_tree(tree, i_cs, a2_gc_interp_lim, a2_bc_neumann_zero)
90 
91  call a2_print_info(tree)
92 
93  time_steps = 0
94 
95  max_diffusion_coeff = max(dc_cs_ion_parallel, dc_cs_ion_perpendicular, &
96  dc_cs)
97 
98  ! Starting simulation
99  do
100  time_steps = time_steps + 1
101  dr_min = a2_min_dr(tree)
102 
103  ! Limit time step for diffusion
104  dt = 0.25_dp * minval(dr_min)**2 / max_diffusion_coeff
105 
106  ! Limit for reactions (using the maximum density, this can be improved)
107  dt = min(dt, 0.5_dp/(k_attachment*initial_cs_density), &
108  0.5_dp/(k_recombination*initial_cs_density), &
109  0.5_dp/k_photoionization)
110 
111  n_steps = ceiling(dt_adapt/dt)
112  dt = dt_adapt / n_steps
113 
114  if (output_cnt * dt_output <= time) then
115  output_cnt = output_cnt + 1
116  write(fname, "(A,I0)") "space_plasma_2d_", output_cnt
117  call a2_gc_tree(tree, i_cs_ion, a2_gc_interp_lim, a2_bc_neumann_zero)
118  call a2_gc_tree(tree, i_cs, a2_gc_interp_lim, a2_bc_neumann_zero)
119  call a2_gc_tree(tree, i_o, a2_gc_interp_lim, a2_bc_neumann_zero)
120  call a2_gc_tree(tree, i_o_min, a2_gc_interp_lim, a2_bc_neumann_zero)
121 
122  call a2_write_silo(tree, trim(fname), output_cnt, time, &
123  ixs_cc=[i_cs_ion, i_cs, i_o, i_o_min], dir="output")
124  end if
125 
126  if (time > end_time) exit
127 
128  do n = 1, n_steps
129  ! Copy previous solution
130  call a2_tree_copy_cc(tree, i_cs_ion, i_cs_ion_old)
131  call a2_tree_copy_cc(tree, i_cs, i_cs_old)
132  call a2_tree_copy_cc(tree, i_o, i_o_old)
133  call a2_tree_copy_cc(tree, i_o_min, i_o_min_old)
134 
135  ! Explicit trapezoidal rule. First take two forward Euler steps over dt
136  do i = 1, 2
137  ! Call procedure get_fluxes for each id in tree, giving the list of boxes
138  call a2_loop_boxes(tree, get_fluxes)
139 
140  ! Restrict fluxes from children to parents on refinement boundaries.
141  call a2_consistent_fluxes(tree, [i_cs_ion, i_cs])
142 
143  ! Call procedure update_solution (see below) for each box in tree, with argument dt
144  call a2_loop_box_arg(tree, update_solution, [dt])
145 
146  ! Restrict variables i_Cs_ion to all parent boxes
147  call a2_restrict_tree(tree, i_cs_ion)
148  call a2_restrict_tree(tree, i_cs)
149  end do
150 
151  ! Now take average of phi_old and phi (this is the explicit trapezoidal rule)
152  call a2_loop_box(tree, average_solution)
153  time = time + dt
154  end do
155 
156  ! Fill ghost cells before refinement
157  call a2_gc_tree(tree, i_cs_ion, a2_gc_interp_lim, a2_bc_neumann_zero)
158  call a2_gc_tree(tree, i_cs, a2_gc_interp_lim, a2_bc_neumann_zero)
159 
160  call a2_adjust_refinement(tree, refine_routine, refine_info, 1)
161  end do
162 
163 contains
164 
166  subroutine refine_routine(box, cell_flags)
167  type(box2_t), intent(in) :: box
168  integer, intent(out) :: cell_flags(box%n_cell, box%n_cell)
169  real(dp) :: diff
170  integer :: i, j, nc
171 
172  nc = box%n_cell
173 
174  do j = 1, nc; do i = 1, nc
175  diff = abs(box%cc(i+1, j, i_cs_ion) + &
176  box%cc(i-1, j, i_cs_ion) + &
177  box%cc(i, j+1, i_cs_ion) + &
178  box%cc(i, j-1, i_cs_ion) - &
179  4 * box%cc(i, j, i_cs_ion)) * box%dr / &
180  (initial_cs_density * initial_ionization_fraction)
181 
182  if (box%lvl < min_refinement_lvl .or. diff > 0.1_dp &
183  .and. box%lvl < max_refinement_lvl) then
184  cell_flags(i, j) = af_do_ref
185  else if (diff < 0.01_dp) then
186  cell_flags(i, j) = af_rm_ref
187  else
188  cell_flags(i, j) = af_keep_ref
189  end if
190  end do; end do
191  end subroutine refine_routine
192 
194  subroutine set_initial_condition(box)
195  type(box2_t), intent(inout) :: box
196  integer :: i, j, nc
197  real(dp) :: distance, density
198 
199  nc = box%n_cell
200  do j = 0, nc+1; do i = 0, nc+1
201  distance = norm2(a2_r_cc(box, [i, j]))
202  density = exp(-0.5_dp * distance**2/initial_sigma**2)
203 
204  box%cc(i, j, i_cs_ion) = initial_ionization_fraction * density
205  box%cc(i, j, i_cs) = (1 - initial_ionization_fraction) * density
206  box%cc(i, j, i_o) = initial_o_density
207  end do; end do
208  end subroutine set_initial_condition
209 
210  subroutine get_fluxes(boxes, id)
211  use m_flux_schemes
212  type(box2_t), intent(inout) :: boxes(:)
213  integer, intent(in) :: id
214  integer :: nc
215  real(dp) :: inv_dx
216  real(dp), allocatable :: diff_coeff(:, :, :)
217 
218  nc = boxes(id)%n_cell
219  inv_dx = 1/boxes(id)%dr
220 
221  allocate(diff_coeff(1:nc+1, 1:nc+1, 2))
222 
223  call a2_gc_box(boxes, id, i_cs_ion, a2_gc_interp_lim, a2_bc_neumann_zero)
224  call a2_gc_box(boxes, id, i_cs, a2_gc_interp_lim, a2_bc_neumann_zero)
225 
226  diff_coeff(:, :, 1) = dc_cs_ion_parallel
227  diff_coeff(:, :, 2) = dc_cs_ion_perpendicular
228  call flux_diff_2d(boxes(id)%cc(:, :, i_cs_ion), diff_coeff, inv_dx, nc, 1)
229  boxes(id)%fc(:, :, :, i_cs_ion) = diff_coeff
230 
231  diff_coeff(:, :, :) = dc_cs
232  call flux_diff_2d(boxes(id)%cc(:, :, i_cs), diff_coeff, inv_dx, nc, 1)
233  boxes(id)%fc(:, :, :, i_cs) = diff_coeff
234 
235  end subroutine get_fluxes
236 
238  subroutine update_solution(box, dt)
239  type(box2_t), intent(inout) :: box
240  real(dp), intent(in) :: dt(:)
241  real(dp) :: inv_dr, recombination, photoionization
242  real(dp) :: attachment, electron_density
243  integer :: i, j, nc
244 
245  nc = box%n_cell
246  inv_dr = 1/box%dr
247 
248  do j = 1, nc
249  do i = 1, nc
250  electron_density = box%cc(i, j, i_cs_ion) - box%cc(i, j, i_o_min)
251  recombination = dt(1) * k_recombination * box%cc(i, j, i_cs_ion)**2
252  photoionization = dt(1) * k_photoionization * box%cc(i, j, i_cs)
253  box%cc(i, j, i_cs_ion) = box%cc(i, j, i_cs_ion) - recombination + photoionization
254  box%cc(i, j, i_cs) = box%cc(i, j, i_cs) + recombination - photoionization
255 
256  attachment = dt(1) * k_attachment * box%cc(i, j, i_o) * electron_density
257  box%cc(i, j, i_o) = box%cc(i, j, i_o) - attachment
258  box%cc(i, j, i_o_min) = box%cc(i, j, i_o_min) + attachment
259 
260  ! Fluxes
261  box%cc(i, j, i_cs_ion) = box%cc(i, j, i_cs_ion) + dt(1) * inv_dr * ( &
262  box%fc(i, j, 1, i_cs_ion) - box%fc(i+1, j, 1, i_cs_ion) + &
263  box%fc(i, j, 2, i_cs_ion) - box%fc(i, j+1, 2, i_cs_ion))
264  box%cc(i, j, i_cs) = box%cc(i, j, i_cs) + dt(1) * inv_dr * ( &
265  box%fc(i, j, 1, i_cs) - box%fc(i+1, j, 1, i_cs) + &
266  box%fc(i, j, 2, i_cs) - box%fc(i, j+1, 2, i_cs))
267  end do
268  end do
269  end subroutine update_solution
270 
272  subroutine average_solution(box)
273  type(box2_t), intent(inout) :: box
274 
275  box%cc(:, :, :) = 0.5_dp * (box%cc(:, :, :) + &
276  box%cc(:, :, :))
277  end subroutine average_solution
278 
279 end program space_plasma_2d
Module containing a couple simple flux schemes for scalar advection/diffusion problems.
Module which contains all Afivo modules, so that a user does not have to include them separately...
Definition: m_a2_all.f90:4