6 program poisson_basic_2d
12 integer,
parameter :: box_size = 8
13 integer,
parameter :: n_boxes_base = 1
14 integer,
parameter :: n_iterations = 10
15 integer,
parameter :: n_var_cell = 6
16 integer,
parameter :: i_phi = 1
17 integer,
parameter :: i_rhs = 2
18 integer,
parameter :: i_err = 3
19 integer,
parameter :: i_tmp = 4
20 integer,
parameter :: i_gradx = 5
21 integer,
parameter :: i_egradx = 6
24 type(ref_info_t) :: refine_info
26 integer :: ix_list(2, n_boxes_base)
27 real(dp) :: dr, residu(2), anal_err(2)
28 character(len=100) :: fname
31 integer :: count_rate,t_start,t_end
33 print *,
"Running poisson_basic_2d" 34 print *,
"Number of threads", af_get_max_threads()
41 call gauss_init(gs, [1.0_dp, 1.0_dp], [0.04_dp, 0.04_dp], &
42 reshape([0.25_dp, 0.25_dp, &
43 0.75_dp, 0.75_dp], [2,2]))
47 dr = 1.0_dp / box_size
57 cc_names=[
"phi",
"rhs",
"err",
"tmp",
"Dx ",
"eDx"])
63 ix_list(:, 1) = [1, 1]
66 call a2_set_base(tree, 1, ix_list)
69 call system_clock(t_start, count_rate)
73 call a2_loop_box(tree, set_initial_condition)
79 call a2_adjust_refinement(tree, refine_routine, refine_info)
82 if (refine_info%n_add == 0)
exit 85 call system_clock(t_end, count_rate)
87 write(*,
"(A,Es10.3,A)")
" Wall-clock time generating AMR grid: ", &
88 (t_end-t_start) /
real(count_rate,dp),
" seconds" 90 call a2_print_info(tree)
97 mg%sides_bc => sides_bc
104 print *,
"Multigrid iteration | max residual | max error" 105 call system_clock(t_start, count_rate)
107 do mg_iter = 1, n_iterations
111 call mg2_fas_fmg(tree, mg, set_residual=.true., have_guess=(mg_iter>1))
114 call a2_loop_box(tree, set_error)
117 call a2_tree_min_cc(tree, i_tmp, residu(1))
118 call a2_tree_max_cc(tree, i_tmp, residu(2))
119 call a2_tree_min_cc(tree, i_err, anal_err(1))
120 call a2_tree_max_cc(tree, i_err, anal_err(2))
121 write(*,
"(I8,13x,2(Es14.5))") mg_iter, maxval(abs(residu)), &
122 maxval(abs(anal_err))
127 write(fname,
"(A,I0)")
"poisson_basic_2d_", mg_iter
128 call a2_write_silo(tree, trim(fname), dir=
"output")
131 call system_clock(t_end, count_rate)
133 write(*,
"(A,I0,A,E10.3,A)") &
134 " Wall-clock time after ", n_iterations, &
135 " iterations: ", (t_end-t_start) /
real(count_rate, dp), &
140 call a2_destroy(tree)
146 subroutine refine_routine(box, cell_flags)
147 type(box2_t),
intent(in) :: box
148 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
150 real(dp) :: rr(2), dr2, drhs
155 do j = 1, nc;
do i = 1, nc
156 rr = a2_r_cc(box, [i, j])
160 drhs = dr2 * box%cc(i, j, i_rhs)
162 if (abs(drhs) > 1e-3_dp .and. box%lvl < 5)
then 163 cell_flags(i, j) = af_do_ref
165 cell_flags(i, j) = af_keep_ref
168 end subroutine refine_routine
173 subroutine set_initial_condition(box)
174 type(box2_t),
intent(inout) :: box
176 real(dp) :: rr(2), grad(2)
180 do j = 1, nc;
do i = 1, nc
182 rr = a2_r_cc(box, [i, j])
185 box%cc(i, j, i_rhs) = gauss_laplacian(gs, rr)
186 call gauss_gradient(gs, rr, grad)
187 box%cc(i, j, i_gradx) = grad(1)
189 end subroutine set_initial_condition
194 subroutine set_error(box)
195 type(box2_t),
intent(inout) :: box
197 real(dp) :: rr(2), dx, gradx
202 do j = 1, nc;
do i = 1, nc
203 rr = a2_r_cc(box, [i, j])
204 box%cc(i, j, i_err) = box%cc(i, j, i_phi) - gauss_value(gs, rr)
205 gradx = 0.5_dp * (box%cc(i+1, j, i_phi) - &
206 box%cc(i-1, j, i_phi)) / dx
207 box%cc(i, j, i_egradx) = gradx - box%cc(i, j, i_gradx)
210 end subroutine set_error
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 bc_type = af_bc_dirichlet
230 case (a2_neighb_lowx)
232 rr = a2_rr_cc(box, [0.5_dp,
real(n, dp)])
233 box%cc(0, n, iv) = gauss_value(gs, rr)
235 case (a2_neighb_highx)
237 rr = a2_rr_cc(box, [nc+0.5_dp,
real(n, dp)])
238 box%cc(nc+1, n, iv) = gauss_value(gs, rr)
240 case (a2_neighb_lowy)
242 rr = a2_rr_cc(box, [
real(n, dp), 0.5_dp])
243 box%cc(n, 0, iv) = gauss_value(gs, rr)
245 case (a2_neighb_highy)
247 rr = a2_rr_cc(box, [
real(n, dp), nc+0.5_dp])
248 box%cc(n, nc+1, iv) = gauss_value(gs, rr)
251 end subroutine sides_bc
253 end program poisson_basic_2d
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...