Example showing how to use multigrid in combination with Neumann boundary conditions. There is also a test case for a cylindrical geometry.
1 #include "../src/cpp_macros.h"
6 program poisson_neumann
11 integer,
parameter :: box_size = 8
12 integer,
parameter :: n_iterations = 10
16 type(ref_info_t) :: refine_info
18 character(len=100) :: fname, arg
20 integer :: count_rate, t_start, t_end, coord
21 logical :: cylindrical = .false.
23 print *,
"Running poisson_neumann_" // dimname //
""
24 print *,
"Number of threads", af_get_max_threads()
26 if (command_argument_count() > 0)
then
27 call get_command_argument(1, arg)
28 read (arg, *) cylindrical
37 call af_add_cc_variable(tree,
"phi", ix=mg%i_phi)
38 call af_add_cc_variable(tree,
"rhs", ix=mg%i_rhs)
39 call af_add_cc_variable(tree,
"tmp", ix=mg%i_tmp)
40 call af_add_cc_variable(tree,
"err", ix=i_err)
49 call system_clock(t_start, count_rate)
51 call af_loop_box(tree, set_initial_condition)
52 call af_adjust_refinement(tree, refine_routine, refine_info, 0)
53 if (refine_info%n_add == 0)
exit
55 call system_clock(t_end, count_rate)
57 write(*,
"(A,Es10.3,A)")
" Wall-clock time generating AMR grid: ", &
58 (t_end-t_start) / real(count_rate,dp),
" seconds"
60 call af_print_info(tree)
64 mg%sides_bc => sides_bc_cyl
66 mg%sides_bc => sides_bc
69 call mg_init(tree, mg)
71 print *,
"Multigrid iteration | max residual | max error"
72 call system_clock(t_start, count_rate)
74 do mg_iter = 1, n_iterations
75 call mg_fas_fmg(tree, mg, set_residual=.true., have_guess=(mg_iter>1))
77 call af_loop_box(tree, set_error)
79 write(fname,
"(A,I0)")
"output/poisson_neumann_" // dimname //
"_", mg_iter
80 call af_write_silo(tree, trim(fname))
82 call system_clock(t_end, count_rate)
84 write(*,
"(A,I0,A,E10.3,A)") &
85 " Wall-clock time after ", n_iterations, &
86 " iterations: ", (t_end-t_start) / real(count_rate, dp), &
92 subroutine refine_routine(box, cell_flags)
93 type(box_t),
intent(in) :: box
94 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
97 if (box%lvl <= 4 .and. all(box%r_min < 0.25_dp))
then
98 cell_flags(dtimes(:)) = af_do_ref
100 cell_flags(dtimes(:)) = af_keep_ref
102 end subroutine refine_routine
105 subroutine set_initial_condition(box)
106 type(box_t),
intent(inout) :: box
108 if (cylindrical)
then
111 box%cc(dtimes(:), mg%i_rhs) = 4.0_dp
113 box%cc(dtimes(:), mg%i_rhs) = 0.0_dp
115 end subroutine set_initial_condition
118 subroutine set_error(box)
119 type(box_t),
intent(inout) :: box
125 rr = af_r_cc(box, [ijk])
126 box%cc(ijk, i_err) = box%cc(ijk, mg%i_phi) - solution(rr(1))
128 end subroutine set_error
130 real(dp)
elemental function solution(x)
131 real(dp),
intent(in) :: x
133 if (cylindrical)
then
138 end function solution
141 subroutine sides_bc(box, nb, iv, coords, bc_val, bc_type)
142 type(box_t),
intent(in) :: box
143 integer,
intent(in) :: nb
144 integer,
intent(in) :: iv
145 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
146 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
147 integer,
intent(out) :: bc_type
154 case (af_neighb_lowx)
155 bc_type = af_bc_dirichlet
157 case (af_neighb_highx)
158 bc_type = af_bc_neumann
161 bc_type = af_bc_neumann
164 end subroutine sides_bc
167 subroutine sides_bc_cyl(box, nb, iv, coords, bc_val, bc_type)
168 type(box_t),
intent(in) :: box
169 integer,
intent(in) :: nb
170 integer,
intent(in) :: iv
171 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
172 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
173 integer,
intent(out) :: bc_type
180 case (af_neighb_lowx)
181 bc_type = af_bc_neumann
183 case (af_neighb_highx)
184 bc_type = af_bc_neumann
188 bc_type = af_bc_dirichlet
189 bc_val = solution(coords(1, :))
191 end subroutine sides_bc_cyl
Module which contains all Afivo modules, so that a user does not have to include them separately.