Afivo  0.3
random_refinement_2d.f90

This example shows how to create an AMR tree, perform random refinement, write output files, and how to fill ghost cells.

1 
6 program random_refinement_2d
7  use m_a2_all
8 
9  implicit none
10 
11  type(a2_t) :: tree
12  integer :: id, iter, boxes_used
13  integer, parameter :: n_boxes_base = 2
14  integer :: ix_list(2, n_boxes_base)
15  integer :: nb_list(a2_num_neighbors, n_boxes_base)
16  integer, parameter :: box_size = 8
17  integer, parameter :: i_phi = 1
18  type(ref_info_t) :: ref_info
19  real(dp) :: dr, sum_phi_t0, sum_phi
20  character(len=40) :: fname
21  integer :: count_rate,t_start,t_end
22 
23  write(*,'(A)') 'program random_refinement_2d'
24 
25  print *, "Number of threads", af_get_max_threads()
26 
27  ! The cell spacing at the coarsest grid level
28  dr = 2 * acos(-1.0_dp) / box_size ! 2 * pi / box_size
29 
30  ! Initialize tree
31  call a2_init(tree, & ! Tree to initialize
32  box_size, & ! A box contains box_size**DIM cells
33  1, & ! Number of cell-centered variables
34  0, & ! Number of face-centered variables
35  dr, & ! Distance between cells on base level
36  cc_names = ["phi"]) ! Optional: names of cell-centered variables
37 
38  call a2_set_cc_methods(tree, 1, a2_bc_dirichlet_zero, &
39  prolong=a2_prolong_linear_cons)
40 
41  ! Set up geometry.
42  ! Neighbors for the boxes are stored in nb_list
43  nb_list(:, :) = af_no_box ! Default value
44 
45  ! Spatial indices are used to define the coordinates of a box.
46  do id = 1, n_boxes_base
47  ix_list(:, id) = [id, 1] ! Boxes at (1,1), (2,1), etc.
48  nb_list(a2_neighb_highy, id) = id ! Periodic in y-direction
49  end do
50 
51  ! Periodic boundaries only have to be specified from one side. The lower-x
52  ! neighbor of box 1 is the last box
53  nb_list(a2_neighb_lowx, 1) = n_boxes_base
54 
55  ! Create the base mesh, using the box indices and their neighbor information
56  call a2_set_base(tree, n_boxes_base, ix_list, nb_list)
57  call a2_print_info(tree)
58 
59  ! Set variables on base by using the helper functions a2_loop_box(tree, sub)
60  ! and a2_loop_boxes(tree, sub). These functions call the subroutine sub for
61  ! each box in the tree, with a slightly different syntax.
62  call a2_loop_box(tree, set_init_cond)
63 
64  ! Fill ghost cells for phi. The third argument is a subroutine that fills
65  ! ghost cells near refinement boundaries, and the fourth argument fill ghost
66  ! cells near physical boundaries.
67  call a2_gc_tree(tree, i_phi, a2_gc_interp, a2_bc_dirichlet_zero)
68 
69  call a2_tree_sum_cc(tree, i_phi, sum_phi_t0)
70 
71  call system_clock(t_start, count_rate)
72  boxes_used = 1
73  do iter = 1, 15
74  ! This writes a VTK output file containing the cell-centered values of the
75  ! leaves of the tree (the boxes not covered by refinement).
76  ! Variables are the names given as the third argument.
77  write(fname, "(A,I0)") "random_refinement_2d_", iter
78  call a2_write_silo(tree, trim(fname), dir="output", n_cycle=iter)
79 
80  ! This updates the refinement of the tree, by at most one level per call.
81  ! The second argument is a subroutine that is called for each box that can
82  ! be refined or derefined, and it should set refinement flags. Information
83  ! about the changes in refinement are returned in the third argument.
84  !
85  ! Within each (de)refinement step the subroutine a2_tidy_up is called. This
86  ! means that every now and then whether too much holes appear in the tree
87  ! box list. Therefore reordering of tree is not necessary at the end of the
88  ! iteration process
89  call a2_adjust_refinement(tree, ref_routine, ref_info, ref_buffer=0)
90 
91  call a2_tree_sum_cc(tree, i_phi, sum_phi)
92  write(*, "(A,E10.2)") " conservation error: ", sum_phi - sum_phi_t0
93  boxes_used = boxes_used + ref_info%n_add - ref_info%n_rm
94  write(*,'(4(3x,A,1x,i6))') "# new boxes", ref_info%n_add, &
95  "# removed boxes", ref_info%n_rm, &
96  "# boxes used ", boxes_used, &
97  " highest level ", tree%highest_lvl
98  end do
99  call system_clock(t_end, count_rate)
100 
101  write(*, '(A,i3,1x,A,f8.2,1x,A,/)') &
102  ' Wall-clock time after ',iter, &
103  ' iterations: ', (t_end-t_start) / real(count_rate, dp), &
104  ' seconds'
105 
106  call a2_print_info(tree)
107 
108  ! This call is not really necessary here, but cleaning up the data in a tree
109  ! is important if your program continues with other tasks.
110  call a2_destroy(tree)
111 
112 contains
113 
114  ! Return the refinement flag for boxes(id)
115  subroutine ref_routine(box, cell_flags)
116  type(box2_t), intent(in) :: box ! A list of all boxes in the tree
117  integer, intent(out) :: cell_flags(box%n_cell, box%n_cell)
118  real(dp) :: rr
119 
120  ! Draw a [0, 1) random number
121  call random_number(rr)
122 
123  if (rr < 0.5_dp**2 .and. box%lvl < 5) then
124  cell_flags = af_do_ref ! Add refinement
125  else
126  cell_flags = af_rm_ref ! Ask to remove this box, which will not always
127  ! happen (see documentation)
128  end if
129  end subroutine ref_routine
130 
131  ! This routine sets the initial conditions for each box
132  subroutine set_init_cond(box)
133  type(box2_t), intent(inout) :: box
134  integer :: i, j, nc
135  real(dp) :: rr(2)
136 
137  nc = box%n_cell
138  do j = 1, nc; do i = 1, nc
139  ! Get the coordinate of the cell center at i,j
140  rr = a2_r_cc(box, [i, j])
141 
142  ! Set the values at each cell according to some function
143  box%cc(i, j, i_phi) = sin(0.5_dp * rr(1)) * cos(rr(2))
144  end do; end do
145  end subroutine set_init_cond
146 
147 end program random_refinement_2d