Example showing how to use multigrid and compare with an analytic solution, using the method of manufactured solutions. A standard 5-point Laplacian is used in cylindrical coordinates.
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 = 4
16 integer,
parameter :: i_phi = 1
17 integer,
parameter :: i_rhs = 2
18 integer,
parameter :: i_err = 3
19 integer,
parameter :: i_tmp = 4
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" 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.0_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"])
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)
81 mg%sides_bc => sides_bc
84 mg%box_op => mg2_auto_op
85 mg%box_gsrb => mg2_auto_gsrb
86 mg%box_corr => mg2_auto_corr
94 print *,
"Multigrid iteration | max residual | max error" 95 call system_clock(t_start, count_rate)
97 do mg_iter = 1, n_iterations
101 call mg2_fas_fmg(tree, mg, .true., mg_iter>1)
104 call a2_loop_box(tree, set_err)
107 call a2_tree_min_cc(tree, i_tmp, residu(1))
108 call a2_tree_max_cc(tree, i_tmp, residu(2))
109 call a2_tree_min_cc(tree, i_err, anal_err(1))
110 call a2_tree_max_cc(tree, i_err, anal_err(2))
111 write(*,
"(I8,2Es14.5)") mg_iter, maxval(abs(residu)), &
112 maxval(abs(anal_err))
114 write(fname,
"(A,I0)")
"poisson_cyl_", mg_iter
115 call a2_write_vtk(tree, trim(fname), dir=
"output")
117 call system_clock(t_end, count_rate)
119 write(*,
"(A,I0,A,E10.3,A)") &
120 " Wall-clock time after ", n_iterations, &
121 " iterations: ", (t_end-t_start) /
real(count_rate, dp), &
126 call a2_destroy(tree)
131 subroutine ref_routine(box, cell_flags)
132 type(box2_t),
intent(in) :: box
133 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
142 crv = box%dr**2 * abs(box%cc(i, j, i_rhs))
145 if (crv > 5.0e-4_dp)
then 146 cell_flags(i, j) = af_do_ref
148 cell_flags(i, j) = af_keep_ref
152 end subroutine ref_routine
155 subroutine set_init_cond(box)
156 type(box2_t),
intent(inout) :: box
164 rz = a2_r_cc(box, [i,j])
165 box%cc(i, j, i_rhs) = gauss_laplacian_cyl(gs, rz)
168 end subroutine set_init_cond
171 subroutine set_err(box)
172 type(box2_t),
intent(inout) :: box
179 rz = a2_r_cc(box, [i,j])
180 box%cc(i, j, i_err) = box%cc(i, j, i_phi) - gauss_value(gs, rz)
183 end subroutine set_err
189 subroutine sides_bc(box, nb, iv, bc_type)
190 type(box2_t),
intent(inout) :: box
191 integer,
intent(in) :: nb
192 integer,
intent(in) :: iv
193 integer,
intent(out) :: bc_type
200 case (a2_neighb_lowx)
201 bc_type = af_bc_neumann
202 box%cc(0, 1:nc, iv) = 0
203 case (a2_neighb_highx)
204 bc_type = af_bc_dirichlet
206 rz = a2_rr_cc(box, [nc+0.5_dp,
real(n, dp)])
207 box%cc(nc+1, n, iv) = gauss_value(gs, rz)
209 case (a2_neighb_lowy)
210 bc_type = af_bc_dirichlet
212 rz = a2_rr_cc(box, [
real(n, dp), 0.5_dp])
213 box%cc(n, 0, iv) = gauss_value(gs, rz)
215 case (a2_neighb_highy)
216 bc_type = af_bc_dirichlet
218 rz = a2_rr_cc(box, [
real(n, dp), nc+0.5_dp])
219 box%cc(n, nc+1, iv) = gauss_value(gs, rz)
222 end subroutine sides_bc
224 end program poisson_cyl