Afivo  0.3
ghostcell_benchmark_3d.f90
1 
5 program ghostcell_benchmark_3d
6  use m_a3_all
7 
8  implicit none
9 
10  integer, parameter :: n_boxes_base = 1
11  integer, parameter :: n_var_cell = 1
12  integer, parameter :: i_phi = 1
13 
14  type(a3_t) :: tree
15  type(ref_info_t) :: ref_info
16  integer :: n_args
17  integer :: n_cell, it, n_iterations, max_ref_lvl
18  integer :: ix_list(3, n_boxes_base)
19  real(dp) :: dr, time
20  character(len=100) :: arg_string
21  integer :: count_rate, t_start, t_end
22 
23  print *, "Running ghostcell_benchmark_3d"
24  print *, "Number of threads", af_get_max_threads()
25 
26  ! Get box size and mesh size from command line argument
27  n_args = command_argument_count()
28 
29  if (n_args >= 1) then
30  call get_command_argument(1, arg_string)
31  read(arg_string, *) n_cell
32  else
33  print *, "No arguments specified, using default values"
34  print *, "Usage: ./ghostcell_benchmark_3d n_cell max_ref_lvl n_iterations"
35  print *, ""
36  n_cell = 16
37  end if
38 
39  if (n_args >= 2) then
40  call get_command_argument(2, arg_string)
41  read(arg_string, *) max_ref_lvl
42  else
43  max_ref_lvl = 4
44  end if
45 
46  if (n_args >= 3) then
47  call get_command_argument(3, arg_string)
48  read(arg_string, *) n_iterations
49  else
50  n_iterations = 100
51  end if
52 
53  print *, "Box size: ", n_cell
54  print *, "Max refinement lvl: ", max_ref_lvl
55  print *, "Num iterations: ", n_iterations
56 
57  dr = 1.0_dp / n_cell
58 
59  ! Initialize tree
60  call a3_init(tree, & ! Tree to initialize
61  n_cell, & ! A box contains box_size**DIM cells
62  n_var_cell, & ! Number of cell-centered variables
63  0, & ! Number of face-centered variables
64  dr, & ! Distance between cells on base level
65  cc_names=["phi"]) ! Variable names
66 
67  ! Set up geometry. These indices are used to define the coordinates of a box,
68  ! by default the box at [1,1] touches the origin (x,y) = (0,0)
69  ix_list(:, 1) = [1, 1, 1] ! Set index of box 1
70 
71  ! Create the base mesh, using the box indices and their neighbor information
72  call a3_set_base(tree, 1, ix_list)
73 
74  call system_clock(t_start, count_rate)
75  do
76  ! For each box, set the initial conditions
77  call a3_loop_box(tree, set_init_cond)
78 
79  ! This updates the refinement of the tree, by at most one level per call.
80  ! The second argument is a subroutine that is called for each box that can
81  ! be refined or derefined, and it should set refinement flags. Information
82  ! about the changes in refinement are returned in the third argument.
83  call a3_adjust_refinement(tree, ref_routine, ref_info)
84 
85  ! If no new boxes have been added, exit the loop
86  if (ref_info%n_add == 0) exit
87  end do
88  call system_clock(t_end, count_rate)
89 
90  write(*,"(A,Es10.3,A)") " Wall-clock time generating AMR grid: ", &
91  (t_end-t_start) / real(count_rate,dp), " seconds"
92 
93  call a3_print_info(tree)
94 
95  ! Do the actual benchmarking
96  call system_clock(t_start, count_rate)
97  do it = 1, n_iterations
98  call a3_gc_tree(tree, i_phi, a3_gc_interp, a3_bc_dirichlet_zero, .false.)
99  end do
100  call system_clock(t_end, count_rate)
101 
102  time = (t_end-t_start) / real(count_rate, dp)
103  write(*, "(A,I0,A,E10.3,A)") &
104  " Wall-clock time after ", n_iterations, &
105  " iterations: ", time, " seconds"
106  write(*, "(A,E10.3,A)") " Per iteration: ", time/n_iterations, " seconds"
107 
108 contains
109 
110  ! Set refinement flags for box
111  subroutine ref_routine(box, cell_flags)
112  type(box3_t), intent(in) :: box ! A list of all boxes in the tree
113  integer, intent(out) :: cell_flags(box%n_cell, box%n_cell, box%n_cell)
114 
115  ! Fully refine up to max_ref_lvl
116  if (box%lvl < max_ref_lvl) then
117  cell_flags = af_do_ref
118  else
119  cell_flags = af_keep_ref
120  end if
121  end subroutine ref_routine
122 
123  ! This routine sets the initial conditions for each box
124  subroutine set_init_cond(box)
125  type(box3_t), intent(inout) :: box
126  integer :: nc
127 
128  nc = box%n_cell
129  box%cc(1:nc, 1:nc, 1:nc, i_phi) = 1.0_dp
130  end subroutine set_init_cond
131 
132 end program ghostcell_benchmark_3d
Module which contains all Afivo modules, so that a user does not have to include them separately...
Definition: m_a3_all.f90:4