Afivo  0.3
poisson_benchmark.f90

This program can be used to benchmark the multigrid routines. For simplicity, it does not compare results with a known solution.

1 #include "../src/cpp_macros.h"
2 
6 program poisson_benchmark
7  use m_af_all
8 
9  implicit none
10 
11  integer :: i_phi
12  integer :: i_rhs
13  integer :: i_tmp
14 
15  type(af_t) :: tree
16  type(ref_info_t) :: ref_info
17  integer :: mg_iter, n_args, coarse_grid_size
18  integer :: n_cell, n_iterations, max_ref_lvl
19  real(dp) :: time, runtime
20  character(len=100) :: arg_string !, fname
21  type(mg_t) :: mg
22  integer :: count_rate,t_start, t_end
23 
24  print *, "Running poisson_benchmark_" // dimname // ""
25  print *, "Number of threads", af_get_max_threads()
26 
27  ! Get box size and mesh size from command line argument
28  n_args = command_argument_count()
29 
30  if (n_args >= 1) then
31  call get_command_argument(1, arg_string)
32  read(arg_string, *) n_cell
33  else
34  print *, "No arguments specified, using default values"
35  print *, "Usage: ./poisson_benchmark_" // dimname // &
36  " n_cell coarse_grid_size max_ref_lvl runtime(s)"
37  n_cell = 16
38  end if
39 
40  if (n_args >= 2) then
41  call get_command_argument(2, arg_string)
42  read(arg_string, *) coarse_grid_size
43  else
44  coarse_grid_size = 16
45  end if
46 
47  if (n_args >= 3) then
48  call get_command_argument(3, arg_string)
49  read(arg_string, *) max_ref_lvl
50  else
51  max_ref_lvl = 2
52  end if
53 
54  if (n_args >= 4) then
55  call get_command_argument(4, arg_string)
56  read(arg_string, *) runtime
57  if (runtime < 1e-3_dp) stop "Run time should be > 1e-3 seconds"
58  else
59  runtime = 0.2_dp
60  end if
61 
62  print *, "Box size: ", n_cell
63  print *, "Coarse grid size: ", coarse_grid_size
64  print *, "Max refinement lvl: ", max_ref_lvl
65  print *, "Run time (s): ", runtime
66 
67  call af_add_cc_variable(tree, "phi", ix=i_phi)
68  call af_add_cc_variable(tree, "rhs", ix=i_rhs)
69  call af_add_cc_variable(tree, "tmp", ix=i_tmp)
70 
71  ! Initialize tree
72  call af_init(tree, & ! Tree to initialize
73  n_cell, & ! A box contains box_size**DIM cells
74  [dtimes(1.0_dp)], &
75  [dtimes(coarse_grid_size)])
76 
77  call system_clock(t_start, count_rate)
78  do
79  ! For each box, set the initial conditions
80  call af_loop_box(tree, set_init_cond)
81 
82  ! This updates the refinement of the tree, by at most one level per call.
83  ! The second argument is a subroutine that is called for each box that can
84  ! be refined or derefined, and it should set refinement flags. Information
85  ! about the changes in refinement are returned in the third argument.
86  call af_adjust_refinement(tree, ref_routine, ref_info)
87 
88  ! If no new boxes have been added, exit the loop
89  if (ref_info%n_add == 0) exit
90  end do
91  call system_clock(t_end, count_rate)
92 
93  write(*,"(A,Es10.3,A)") " Wall-clock time generating AMR grid: ", &
94  (t_end-t_start) / real(count_rate,dp), " seconds"
95 
96  call af_print_info(tree)
97 
98  ! Set the multigrid options.
99  mg%i_phi = i_phi ! Solution variable
100  mg%i_rhs = i_rhs ! Right-hand side variable
101  mg%i_tmp = i_tmp ! Variable for temporary space
102  mg%sides_bc => af_bc_dirichlet_zero ! Method for boundary conditions
103 
104  ! Initialize the multigrid options. This performs some basics checks and sets
105  ! default values where necessary.
106  ! This routine does not initialize the multigrid variables i_phi, i_rhs
107  ! and i_tmp. These variables will be initialized at the first call of mg_fas_fmg
108  call mg_init(tree, mg)
109 
110  ! Warm-up call
111  call mg_fas_fmg(tree, mg, .false., .false.)
112 
113  ! Test how long cycles take to determine the number of cycles
114  n_iterations = 1000
115  call system_clock(t_start, count_rate)
116  do mg_iter = 1, n_iterations
117  call mg_fas_fmg(tree, mg, .false., mg_iter>1)
118  call system_clock(t_end, count_rate)
119  time = (t_end-t_start) / real(count_rate, dp)
120  if (time > 0.2_dp * runtime) exit
121  end do
122 
123  n_iterations = ceiling((runtime/time) * mg_iter)
124 
125  ! Do the actual benchmarking
126  call system_clock(t_start, count_rate)
127  do mg_iter = 1, n_iterations
128  ! Perform a FAS-FMG cycle (full approximation scheme, full multigrid). The
129  ! third argument controls whether the residual is stored in i_tmp. The
130  ! fourth argument controls whether to improve the current solution.
131  call mg_fas_fmg(tree, mg, .false., mg_iter>1)
132 
133  ! This uses a V-cycle instead of an FMG-cycle
134  ! call mg_fas_vcycle(tree, mg, .false.)
135 
136  ! If uncommented, this writes Silo output files containing the
137  ! cell-centered values of the leaves of the tree
138  ! write(fname, "(A,I0)") "poisson_benchmark_" // DIMNAME // "_", mg_iter
139  ! call af_write_silo(tree, trim(fname), dir="output")
140  end do
141  call system_clock(t_end, count_rate)
142 
143  time = (t_end-t_start) / real(count_rate, dp)
144  write(*, "(A,I0,A,E10.3,A)") &
145  " Wall-clock time after ", n_iterations, &
146  " iterations: ", time, " seconds"
147  write(*, "(A,E10.3,A)") " Per iteration: ", time/n_iterations, " seconds"
148 
149  ! This writes a Silo output file containing the cell-centered values of the
150  ! leaves of the tree (the boxes not covered by refinement).
151  ! fname = "poisson_benchmark_" // DIMNAME // ""
152  ! call af_write_silo(tree, trim(fname), dir="output")
153 
154  ! This call is not really necessary here, but cleaning up the data in a tree
155  ! is important if your program continues with other tasks.
156  call af_destroy(tree)
157 
158 contains
159 
160  ! Set refinement flags for box
161  subroutine ref_routine(box, cell_flags)
162  type(box_t), intent(in) :: box ! A list of all boxes in the tree
163  integer, intent(out) :: cell_flags(DTIMES(box%n_cell))
164 
165  ! Fully refine up to max_ref_lvl
166  if (box%lvl < max_ref_lvl) then
167  cell_flags = af_do_ref
168  else
169  cell_flags = af_keep_ref
170  end if
171  end subroutine ref_routine
172 
173  ! This routine sets the initial conditions for each box
174  subroutine set_init_cond(box)
175  type(box_t), intent(inout) :: box
176  integer :: nc
177 
178  nc = box%n_cell
179  box%cc(dtimes(1:nc), i_rhs) = 1.0_dp
180  end subroutine set_init_cond
181 
182 end program poisson_benchmark
Module which contains all Afivo modules, so that a user does not have to include them separately.
Definition: m_af_all.f90:3