An implicit diffusion example, showing how the multigrid methods can be used to solve the diffusion equation with a backward Euler scheme.
1 #include "../src/cpp_macros.h"
6 program implicit_diffusion
11 integer,
parameter :: box_size = 8
17 real(dp),
parameter :: domain_len = 2 * acos(-1.0_dp)
18 real(dp),
parameter :: diffusion_coeff = 1.0_dp
20 real(dp),
parameter :: solution_modes(NDIM) = [1]
22 real(dp),
parameter :: solution_modes(NDIM) = [1, 1]
24 real(dp),
parameter :: solution_modes(NDIM) = [1, 1, 0]
29 type(ref_info_t) :: refine_info
30 integer :: time_steps, output_cnt
31 real(dp) :: time, end_time
32 real(dp),
parameter :: dt = 0.1_dp
33 character(len=100) :: fname
35 print *,
"Running implicit_diffusion_" // dimname //
""
36 print *,
"Number of threads", af_get_max_threads()
38 call af_add_cc_variable(tree,
"phi", ix=i_phi)
39 call af_add_cc_variable(tree,
"rhs", ix=i_rhs)
40 call af_add_cc_variable(tree,
"tmp", ix=i_tmp)
41 call af_add_cc_variable(tree,
"err", ix=i_err)
46 [dtimes(domain_len)], &
48 periodic=[dtimes(.true.)])
53 mg%sides_bc => af_bc_dirichlet_zero
54 mg%helmholtz_lambda = 1/(diffusion_coeff * dt)
59 call mg_init(tree, mg)
61 call af_print_info(tree)
67 call af_loop_box(tree, set_initial_condition)
75 call af_gc_tree(tree, [i_phi])
83 call af_adjust_refinement(tree, &
88 if (refine_info%n_add == 0)
exit
99 time_steps = time_steps + 1
101 output_cnt = output_cnt + 1
102 write(fname,
"(A,I0)")
"output/implicit_diffusion_" // dimname //
"_", output_cnt
105 call af_loop_box_arg(tree, set_error, [time])
109 call af_write_silo(tree, trim(fname), output_cnt, time)
111 if (time > end_time)
exit
114 call af_loop_box(tree, set_rhs)
117 call mg_fas_fmg(tree, &
119 set_residual = .true., &
128 subroutine refine_routine(box, cell_flags)
129 type(box_t),
intent(in) :: box
130 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
132 if (maxval(box%dr) > 2e-2_dp * domain_len)
then
133 cell_flags = af_do_ref
135 cell_flags = af_keep_ref
137 end subroutine refine_routine
140 subroutine set_initial_condition(box)
141 type(box_t),
intent(inout) :: box
148 rr = af_r_cc(box, [ijk])
149 box%cc(ijk, i_phi) = solution(rr, 0.0_dp)
151 end subroutine set_initial_condition
154 subroutine set_error(box, time)
155 type(box_t),
intent(inout) :: box
156 real(dp),
intent(in) :: time(:)
162 rr = af_r_cc(box, [ijk])
163 box%cc(ijk, i_err) = &
164 box%cc(ijk, i_phi) - solution(rr, time(1))
166 end subroutine set_error
168 function solution(rr, t)
result(sol)
169 real(dp),
intent(in) :: rr(NDIM), t
170 real(dp) :: sol, tmp(NDIM)
172 tmp = solution_modes * rr
173 sol = 1 + product(cos(tmp)) * &
174 exp(-sum(solution_modes**2) * diffusion_coeff * t)
175 end function solution
178 subroutine set_rhs(box)
179 type(box_t),
intent(inout) :: box
183 box%cc(dtimes(1:nc), i_rhs) = -1/(dt * diffusion_coeff) * &
184 box%cc(dtimes(1:nc), i_phi)
185 end subroutine set_rhs
187 end program implicit_diffusion
Module which contains all Afivo modules, so that a user does not have to include them separately.