5 program poisson_cyl_analytic
10 integer,
parameter :: box_size = 8
11 integer,
parameter :: n_boxes_base = 1
12 integer,
parameter :: n_iterations = 10
13 integer,
parameter :: n_var_cell = 5
14 integer,
parameter :: i_phi = 1
15 integer,
parameter :: i_rhs = 2
16 integer,
parameter :: i_err = 3
17 integer,
parameter :: i_sol = 4
18 integer,
parameter :: i_tmp = 5
20 real(dp),
parameter :: domain_len = 1.25e-2_dp
21 real(dp),
parameter :: pi = acos(-1.0_dp)
22 real(dp),
parameter :: sigma = 4e-4_dp * sqrt(0.5_dp)
23 real(dp),
parameter :: rz_source(2) = [0.0_dp, 0.5_dp] * domain_len
24 real(dp),
parameter :: epsilon_source = 8.85e-12_dp
25 real(dp),
parameter :: q_source = 3e18_dp * 1.6022e-19_dp * &
26 sigma**3 * sqrt(2 * pi)**3
29 type(ref_info_t) :: ref_info
31 integer :: ix_list(2, n_boxes_base)
32 real(dp) :: dr, residu(2), anal_err(2), max_val
33 character(len=100) :: fname
35 integer :: count_rate,t_start, t_end
37 print *,
"Running poisson_cyl_analytic" 38 print *,
"Number of threads", af_get_max_threads()
41 dr = 1.25e-2_dp / box_size
51 cc_names=[
"phi",
"rhs",
"err",
"sol",
"tmp"])
58 call a2_set_base(tree, 1, ix_list)
59 call a2_print_info(tree)
61 call system_clock(t_start, count_rate)
64 call a2_loop_box(tree, set_init_cond)
67 call a2_adjust_refinement(tree, ref_routine, ref_info)
70 if (ref_info%n_add == 0)
exit 72 call system_clock(t_end, count_rate)
74 write(*,
"(A,Es10.3,A)")
" Wall-clock time generating AMR grid: ", &
75 (t_end-t_start) /
real(count_rate,dp),
" seconds" 77 call a2_print_info(tree)
83 mg%sides_bc => sides_bc
86 mg%box_op => mg2_auto_op
87 mg%box_gsrb => mg2_auto_gsrb
88 mg%box_corr => mg2_auto_corr
96 print *,
"Multigrid iteration | max residual | max rel. error" 97 call system_clock(t_start, count_rate)
99 do mg_iter = 1, n_iterations
103 call mg2_fas_fmg(tree, mg, .true., mg_iter>1)
106 call a2_loop_box(tree, set_err)
109 call a2_tree_min_cc(tree, i_tmp, residu(1))
110 call a2_tree_max_cc(tree, i_tmp, residu(2))
111 call a2_tree_min_cc(tree, i_err, anal_err(1))
112 call a2_tree_max_cc(tree, i_err, anal_err(2))
113 call a2_tree_max_cc(tree, i_phi, max_val)
114 anal_err = anal_err / max_val
116 write(*,
"(I8,2Es14.5)") mg_iter, maxval(abs(residu)), &
117 maxval(abs(anal_err))
119 write(fname,
"(A,I0)")
"poisson_cyl_analytic_", mg_iter
120 call a2_write_vtk(tree, trim(fname), dir=
"output")
122 call system_clock(t_end, count_rate)
124 write(*,
"(A,I0,A,E10.3,A)") &
125 " Wall-clock time after ", n_iterations, &
126 " iterations: ", (t_end-t_start) /
real(count_rate, dp), &
131 call a2_destroy(tree)
136 subroutine ref_routine(box, cell_flags)
137 type(box2_t),
intent(in) :: box
138 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
147 crv = box%dr**2 * abs(box%cc(i, j, i_rhs))
150 if (crv > 1e-1_dp .and. box%lvl < 10)
then 151 cell_flags(i, j) = af_do_ref
153 cell_flags(i, j) = af_keep_ref
157 end subroutine ref_routine
160 subroutine set_init_cond(box)
161 type(box2_t),
intent(inout) :: box
169 rz = a2_r_cc(box, [i,j]) - rz_source
172 box%cc(i, j, i_rhs) = - q_source * &
173 exp(-sum(rz**2) / (2 * sigma**2)) / &
174 (sigma**3 * sqrt(2 * pi)**3 * epsilon_source)
177 end subroutine set_init_cond
179 real(dp) function solution(rz_in)
180 real(dp),
intent(in) :: rz_in(2)
183 d = norm2(rz_in - rz_source)
186 if (d < sqrt(epsilon(1.0d0)))
then 187 erf_d = sqrt(2.0_dp/pi) / sigma
189 erf_d = erf(d * sqrt(0.5_dp) / sigma) / d
192 solution = erf_d * q_source / (4 * pi * epsilon_source)
193 end function solution
196 subroutine set_err(box)
197 type(box2_t),
intent(inout) :: box
204 rz = a2_r_cc(box, [i,j])
205 box%cc(i, j, i_sol) = solution(rz)
206 box%cc(i, j, i_err) = box%cc(i, j, i_phi) - box%cc(i, j, i_sol)
209 end subroutine set_err
215 subroutine sides_bc(box, nb, iv, bc_type)
216 type(box2_t),
intent(inout) :: box
217 integer,
intent(in) :: nb
218 integer,
intent(in) :: iv
219 integer,
intent(out) :: bc_type
226 case (a2_neighb_lowx)
227 bc_type = af_bc_neumann
228 box%cc(0, 1:nc, iv) = 0
229 case (a2_neighb_highx)
230 bc_type = af_bc_dirichlet
232 rz = a2_rr_cc(box, [nc+0.5_dp,
real(n, dp)])
233 box%cc(nc+1, n, iv) = solution(rz)
235 case (a2_neighb_lowy)
236 bc_type = af_bc_dirichlet
238 rz = a2_rr_cc(box, [
real(n, dp), 0.5_dp])
239 box%cc(n, 0, iv) = solution(rz)
241 case (a2_neighb_highy)
242 bc_type = af_bc_dirichlet
244 rz = a2_rr_cc(box, [
real(n, dp), nc+0.5_dp])
245 box%cc(n, nc+1, iv) = solution(rz)
248 end subroutine sides_bc
250 end program poisson_cyl_analytic
Module which contains all Afivo modules, so that a user does not have to include them separately...