5 program poisson_cyl_dielectric
11 integer,
parameter :: box_size = 8
12 integer,
parameter :: n_boxes_base = 1
13 integer,
parameter :: n_iterations = 10
14 integer,
parameter :: n_var_cell = 5
15 integer,
parameter :: i_phi = 1
16 integer,
parameter :: i_rhs = 2
17 integer,
parameter :: i_err = 3
18 integer,
parameter :: i_tmp = 4
19 integer,
parameter :: i_eps = 5
22 type(ref_info_t) :: ref_info
24 integer :: ix_list(2, n_boxes_base)
25 real(dp) :: dr, residu(2), anal_err(2)
26 character(len=100) :: fname
29 integer :: count_rate,t_start, t_end
31 print *,
"Running poisson_cyl_dielectric" 32 print *,
"Number of threads", af_get_max_threads()
35 call gauss_init(gs, [1.0_dp, 1.0_dp], [0.04_dp, 0.04_dp], &
36 reshape([0.25_dp, 0.25_dp, 0.75_dp, 0.75_dp], [2,2]))
39 dr = 1.0_dp / box_size
49 cc_names=[
"phi",
"rhs",
"err",
"tmp",
"eps"])
56 call a2_set_base(tree, 1, ix_list)
57 call a2_print_info(tree)
59 call system_clock(t_start, count_rate)
62 call a2_loop_box(tree, set_init_cond)
65 call a2_adjust_refinement(tree, ref_routine, ref_info)
68 if (ref_info%n_add == 0)
exit 70 call system_clock(t_end, count_rate)
72 write(*,
"(A,Es10.3,A)")
" Wall-clock time generating AMR grid: ", &
73 (t_end-t_start) /
real(count_rate,dp),
" seconds" 75 call a2_print_info(tree)
82 mg%sides_bc => sides_bc
85 mg%box_op => mg2_auto_op
86 mg%box_gsrb => mg2_auto_gsrb
87 mg%box_corr => mg2_auto_corr
95 print *,
"Multigrid iteration | max residual | max error" 96 call system_clock(t_start, count_rate)
98 do mg_iter = 1, n_iterations
102 call mg2_fas_fmg(tree, mg, .true., mg_iter>1)
105 call a2_loop_box(tree, set_err)
108 call a2_tree_min_cc(tree, i_tmp, residu(1))
109 call a2_tree_max_cc(tree, i_tmp, residu(2))
110 call a2_tree_min_cc(tree, i_err, anal_err(1))
111 call a2_tree_max_cc(tree, i_err, anal_err(2))
112 write(*,
"(I8,2Es14.5)") mg_iter, maxval(abs(residu)), &
113 maxval(abs(anal_err))
115 write(fname,
"(A,I0)")
"poisson_cyl_dielectric_", mg_iter
116 call a2_write_vtk(tree, trim(fname), dir=
"output")
118 call system_clock(t_end, count_rate)
120 write(*,
"(A,I0,A,E10.3,A)") &
121 " Wall-clock time after ", n_iterations, &
122 " iterations: ", (t_end-t_start) /
real(count_rate, dp), &
127 call a2_destroy(tree)
132 subroutine ref_routine(box, cell_flags)
133 type(box2_t),
intent(in) :: box
134 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
143 crv = box%dr**2 * abs(box%cc(i, j, i_rhs)) / box%cc(i, j, i_eps)
146 if (crv > 1.0e-3_dp)
then 147 cell_flags(i, j) = af_do_ref
149 cell_flags(i, j) = af_keep_ref
153 end subroutine ref_routine
156 subroutine set_init_cond(box)
157 type(box2_t),
intent(inout) :: box
159 real(dp) :: rz(2), grad(2), dr, qbnd, tmp
162 box%cc(:, :, i_phi) = 0
167 rz = a2_r_cc(box, [i,j])
170 if (rz(1) < 0.5_dp .and. rz(2) < 0.5_dp)
then 171 box%cc(i, j, i_eps) = 100.0_dp
173 box%cc(i, j, i_eps) = 1.0_dp
177 box%cc(i, j, i_rhs) = gauss_laplacian_cyl(gs, rz) * box%cc(i, j, i_eps)
185 rz = a2_rr_cc(box, [i + 0.5_dp,
real(j, dp)])
188 call gauss_gradient(gs, rz, grad)
189 qbnd = (box%cc(i+1, j, i_eps) - box%cc(i, j, i_eps)) * &
193 tmp = box%cc(i+1, j, i_eps) / &
194 (box%cc(i, j, i_eps) + box%cc(i+1, j, i_eps))
195 box%cc(i+1, j, i_rhs) = box%cc(i+1, j, i_rhs) + tmp * qbnd
196 box%cc(i, j, i_rhs) = box%cc(i, j, i_rhs) + (1-tmp) * qbnd
203 rz = a2_rr_cc(box, [
real(i, dp), j + 0.5_dp])
206 call gauss_gradient(gs, rz, grad)
207 qbnd = (box%cc(i, j+1, i_eps) - box%cc(i, j, i_eps)) * &
211 tmp = box%cc(i, j+1, i_eps) / &
212 (box%cc(i, j, i_eps) + box%cc(i, j+1, i_eps))
213 box%cc(i, j+1, i_rhs) = box%cc(i, j+1, i_rhs) + tmp * qbnd
214 box%cc(i, j, i_rhs) = box%cc(i, j, i_rhs) + (1-tmp) * qbnd
217 end subroutine set_init_cond
220 subroutine set_err(box)
221 type(box2_t),
intent(inout) :: box
228 rz = a2_r_cc(box, [i,j])
229 box%cc(i, j, i_err) = box%cc(i, j, i_phi) - gauss_value(gs, rz)
232 end subroutine set_err
238 subroutine sides_bc(box, nb, iv, bc_type)
239 type(box2_t),
intent(inout) :: box
240 integer,
intent(in) :: nb
241 integer,
intent(in) :: iv
242 integer,
intent(out) :: bc_type
249 case (a2_neighb_lowx)
250 bc_type = af_bc_neumann
251 box%cc(0, 1:nc, iv) = 0
252 case (a2_neighb_highx)
253 bc_type = af_bc_dirichlet
255 rz = a2_rr_cc(box, [nc+0.5_dp,
real(n, dp)])
256 box%cc(nc+1, n, iv) = gauss_value(gs, rz)
258 case (a2_neighb_lowy)
259 bc_type = af_bc_dirichlet
261 rz = a2_rr_cc(box, [
real(n, dp), 0.5_dp])
262 box%cc(n, 0, iv) = gauss_value(gs, rz)
264 case (a2_neighb_highy)
265 bc_type = af_bc_dirichlet
267 rz = a2_rr_cc(box, [
real(n, dp), nc+0.5_dp])
268 box%cc(n, nc+1, iv) = gauss_value(gs, rz)
271 end subroutine sides_bc
273 end program poisson_cyl_dielectric
This module can be used to construct solutions consisting of one or more Gaussians.
Module which contains all Afivo modules, so that a user does not have to include them separately...