Afivo 0.3
All Classes Namespaces Functions Variables Pages
poisson_coarse_solver.f90
1#include "../src/cpp_macros.h"
2program poisson_coarse_solver
3 use m_af_all
5
6 implicit none
7
8 integer, parameter :: box_size = 8
9 integer, parameter :: domain_size(NDIM) = 16
10 integer, parameter :: n_iterations = 10
11 integer :: i_phi
12 integer :: i_rhs
13 integer :: i_err
14 integer :: i_tmp
15 integer :: i_gradx
16 integer :: i_egradx
17
18 type(af_t) :: tree
19 type(ref_info_t) :: refine_info
20 integer :: mg_iter
21 real(dp) :: residu(2), anal_err(2)
22 character(len=100) :: fname
23 type(mg_t) :: mg
24 type(gauss_t) :: gs
25 integer :: count_rate,t_start,t_end
26
27 print *, "Running poisson_coarse_solver_" // dimname
28 print *, "Number of threads", af_get_max_threads()
29
30 !> [Gauss_init]
31 ! The manufactured solution exists of Gaussians
32 ! Amplitudes: [1.0_dp, 1.0_dp]
33 ! Sigmas : [0.04_dp, 0.04_dp]
34 ! Locations : x, y, z = 0.25 or x, y, z = 0.75
35 call gauss_init(gs, [1.0_dp, 1.0_dp], [0.1_dp, 0.1_dp], &
36 reshape([dtimes(0.25_dp), &
37 dtimes(0.6_dp)], [ndim,2]))
38 !> [Gauss_init]
39
40 !> [af_init]
41 call af_add_cc_variable(tree, "phi", ix=i_phi)
42 call af_add_cc_variable(tree, "rhs", ix=i_rhs)
43 call af_add_cc_variable(tree, "err", ix=i_err)
44 call af_add_cc_variable(tree, "tmp", ix=i_tmp)
45 call af_add_cc_variable(tree, "Dx", ix=i_gradx)
46 call af_add_cc_variable(tree, "eDx", ix=i_egradx)
47
48 ! Initialize tree
49 call af_init(tree, & ! Tree to initialize
50 box_size, & ! A box contains box_size**DIM cells
51 [dtimes(1.0_dp)], &
52 domain_size)
53 !> [af_init]
54
55 call system_clock(t_start, count_rate)
56 !> [set_refinement]
57 do
58 ! For each box, set the initial conditions
59 call af_loop_box(tree, set_initial_condition)
60
61 ! This updates the refinement of the tree, by at most one level per call.
62 ! The second argument is a subroutine that is called for each box that can
63 ! be refined or derefined, and it should set refinement flags. Information
64 ! about the changes in refinement are returned in the third argument.
65 call af_adjust_refinement(tree, refine_routine, refine_info)
66
67 ! If no new boxes have been added, exit the loop
68 if (refine_info%n_add == 0) exit
69 end do
70 !> [set_refinement]
71 call system_clock(t_end, count_rate)
72
73 write(*,"(A,Es10.3,A)") " Wall-clock time generating AMR grid: ", &
74 (t_end-t_start) / real(count_rate,dp), " seconds"
75
76 call af_print_info(tree)
77
78 ! Initialize the multigrid options. This performs some basics checks and sets
79 ! default values where necessary.
80 mg%i_phi = i_phi ! Solution variable
81 mg%i_rhs = i_rhs ! Right-hand side variable
82 mg%i_tmp = i_tmp ! Variable for temporary space
83 mg%sides_bc => sides_bc ! Method for boundary conditions
84
85 ! This routine does not initialize the multigrid fields boxes%i_phi,
86 ! boxes%i_rhs and boxes%i_tmp. These fileds will be initialized at the
87 ! first call of mg_fas_fmg
88 call mg_init(tree, mg)
89
90 print *, "Multigrid iteration | max residual | max error"
91 call system_clock(t_start, count_rate)
92
93 do mg_iter = 1, n_iterations
94 ! Perform a FAS-FMG cycle (full approximation scheme, full multigrid). The
95 ! third argument controls whether the residual is stored in i_tmp. The
96 ! fourth argument controls whether to improve the current solution.
97 call mg_fas_fmg(tree, mg, set_residual=.true., have_guess=(mg_iter>1))
98
99 ! Compute the error compared to the analytic solution
100 call af_loop_box(tree, set_error)
101
102 ! Determine the minimum and maximum residual and error
103 call af_tree_min_cc(tree, i_tmp, residu(1))
104 call af_tree_max_cc(tree, i_tmp, residu(2))
105 call af_tree_min_cc(tree, i_err, anal_err(1))
106 call af_tree_max_cc(tree, i_err, anal_err(2))
107 write(*,"(I8,13x,2(Es14.5))") mg_iter, maxval(abs(residu)), &
108 maxval(abs(anal_err))
109
110 !> [write_output]
111 ! This writes a Silo output file containing the cell-centered values of the
112 ! leaves of the tree (the boxes not covered by refinement).
113 write(fname, "(A,I0)") "output/poisson_coarse_solver_" // dimname // "_", mg_iter
114 call af_write_silo(tree, trim(fname))
115 !> [write_output]
116 end do
117 call system_clock(t_end, count_rate)
118
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), &
122 " seconds"
123
124 ! This call is not really necessary here, but cleaning up the data in a tree
125 ! is important if your program continues with other tasks.
126 call af_destroy(tree)
127
128contains
129
130 !> [refine_routine]
131 ! Return the refinement flags for box
132 subroutine refine_routine(box, cell_flags)
133 type(box_t), intent(in) :: box
134 integer, intent(out) :: cell_flags(DTIMES(box%n_cell))
135 integer :: IJK, nc
136 real(dp) :: rr(NDIM), dr2, drhs
137
138 nc = box%n_cell
139 dr2 = maxval(box%dr)**2
140
141 do kji_do(1,nc)
142 rr = af_r_cc(box, [ijk])
143
144 ! This is an estimate of the truncation error in the right-hand side,
145 ! which is related to the fourth derivative of the solution.
146 drhs = dr2 * box%cc(ijk, i_rhs)
147
148 if (abs(drhs) > 1e-3_dp .and. box%lvl < 5) then
149 cell_flags(ijk) = af_do_ref
150 else
151 cell_flags(ijk) = af_keep_ref
152 end if
153 end do; close_do
154 end subroutine refine_routine
155 !> [refine_routine]
156
157 !> [set_initial_condition]
158 ! This routine sets the initial conditions for each box
159 subroutine set_initial_condition(box)
160 type(box_t), intent(inout) :: box
161 integer :: IJK, nc
162 real(dp) :: rr(NDIM), grad(NDIM)
163
164 nc = box%n_cell
165
166 do kji_do(1,nc)
167 ! Get the coordinate of the cell center at IJK
168 rr = af_r_cc(box, [ijk])
169
170 ! And set the rhs values
171 box%cc(ijk, i_rhs) = gauss_laplacian(gs, rr)
172 call gauss_gradient(gs, rr, grad)
173 box%cc(ijk, i_gradx) = grad(1)
174 end do; close_do
175 end subroutine set_initial_condition
176 !> [set_initial_condition]
177
178 !> [set_error]
179 ! Set the error compared to the analytic solution
180 subroutine set_error(box)
181 type(box_t), intent(inout) :: box
182 integer :: IJK, nc
183 real(dp) :: rr(NDIM), gradx
184
185 nc = box%n_cell
186
187 do kji_do(1,nc)
188 rr = af_r_cc(box, [ijk])
189 box%cc(ijk, i_err) = box%cc(ijk, i_phi) - gauss_value(gs, rr)
190#if NDIM == 1
191 gradx = 0.5_dp * (box%cc(i+1, i_phi) - &
192 box%cc(i-1, i_phi)) / box%dr(1)
193#elif NDIM == 2
194 gradx = 0.5_dp * (box%cc(i+1, j, i_phi) - &
195 box%cc(i-1, j, i_phi)) / box%dr(1)
196#elif NDIM == 3
197 gradx = 0.5_dp * (box%cc(i+1, j, k, i_phi) - &
198 box%cc(i-1, j, k, i_phi)) / box%dr(1)
199#endif
200 box%cc(ijk, i_egradx) = gradx - box%cc(ijk, i_gradx)
201
202 end do; close_do
203 end subroutine set_error
204 !> [set_error]
205
206 ! This routine sets boundary conditions for a box
207 subroutine sides_bc(box, nb, iv, coords, bc_val, bc_type)
208 type(box_t), intent(in) :: box
209 integer, intent(in) :: nb
210 integer, intent(in) :: iv
211 real(dp), intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
212 real(dp), intent(out) :: bc_val(box%n_cell**(NDIM-1))
213 integer, intent(out) :: bc_type
214 integer :: n
215
216 ! We use dirichlet boundary conditions
217 bc_type = af_bc_dirichlet
218
219 ! Below the solution is specified in the approriate ghost cells
220 do n = 1, box%n_cell**(ndim-1)
221 bc_val(n) = gauss_value(gs, coords(:, n))
222 end do
223 end subroutine sides_bc
224
225end program poisson_coarse_solver
Module which contains all Afivo modules, so that a user does not have to include them separately.
Definition: m_af_all.f90:3
This module can be used to construct solutions consisting of one or more Gaussians.
Definition: m_gaussians.f90:3