Afivo  0.3
test_refinement_buffer.f90

This example tests the refinement buffer option, which extends the refinement a given number of cells

1 #include "../src/cpp_macros.h"
2 !> \example test_refinement_buffer.f90
3 !>
4 !> This example tests the refinement buffer option, which extends the
5 !> refinement a given number of cells
6 program random_refinement
7  use m_af_all
8 
9  implicit none
10 
11  type(af_t) :: tree
12  integer :: grid_size(NDIM)
13  real(dp) :: domain_size(NDIM)
14  logical :: periodic(NDIM) = .true.
15  integer :: iter
16  integer, parameter :: coord_type = af_xyz ! af_xyz or af_cyl
17 
18  integer, parameter :: box_size = 8
19  integer :: i_phi
20  type(ref_info_t) :: ref_info
21  character(len=100) :: fname
22 
23  write(*,'(A,I0,A)') 'program test_refinement_buffer_', ndim, "d"
24 
25  print *, "Number of threads", af_get_max_threads()
26 
27  call af_add_cc_variable(tree, "phi", ix=i_phi)
28 
29  call af_set_cc_methods(tree, 1, af_bc_dirichlet_zero, &
30  prolong=af_prolong_linear)
31 
32  ! Initialize tree
33  grid_size(:) = box_size
34  domain_size(:) = 1.0_dp
35 
36  call af_init(tree, & ! Tree to initialize
37  box_size, & ! A box contains box_size**DIM cells
38  domain_size, &
39  grid_size, &
40  periodic=periodic, &
41  coord=coord_type)
42 
43  do iter = 1, 10
44  call af_adjust_refinement(tree, ref_routine, ref_info, ref_buffer=0)
45  if (ref_info%n_add == 0) exit
46  end do
47 
48  call af_loop_box(tree, set_init_cond)
49  call af_gc_tree(tree, [i_phi])
50 
51  do iter = 0, box_size
52  call af_adjust_refinement(tree, ref_routine, ref_info, ref_buffer=iter)
53 
54  call af_loop_box(tree, set_init_cond)
55  call af_gc_tree(tree, [i_phi])
56 
57  write(fname, "(A,I0)") "output/test_refinement_buffer_" // dimname // "_", iter
58  call af_write_silo(tree, trim(fname), n_cycle=iter)
59  end do
60 
61 contains
62 
63  ! Return the refinement flag for boxes(id)
64  subroutine ref_routine(box, cell_flags)
65  type(box_t), intent(in) :: box ! A list of all boxes in the tree
66  integer, intent(out) :: cell_flags(DTIMES(box%n_cell))
67  real(dp) :: rr(NDIM)
68  integer :: nc, IJK
69 
70  nc = box%n_cell
71 
72  do kji_do(1,nc)
73  rr = af_r_cc(box, [ijk])
74  if (norm2(rr - 0.5_dp) < 0.25_dp .and. box%lvl < 5) then
75  cell_flags(ijk) = af_do_ref ! Add refinement
76  else
77  cell_flags(ijk) = af_keep_ref
78  end if
79  end do; close_do
80 
81  end subroutine ref_routine
82 
83  ! This routine sets the initial conditions for each box
84  subroutine set_init_cond(box)
85  type(box_t), intent(inout) :: box
86  integer :: IJK, nc
87  real(dp) :: rr(NDIM)
88 
89  nc = box%n_cell
90  do kji_do(1,nc)
91  ! Get the coordinate of the cell center at i,j
92  rr = af_r_cc(box, [ijk])
93 
94  if (norm2(rr - 0.5_dp) < 0.25_dp) then
95  box%cc(ijk, i_phi) = 1
96  else
97  box%cc(ijk, i_phi) = 0
98  end if
99 
100  end do; close_do
101  end subroutine set_init_cond
102 
103 end program random_refinement
Module which contains all Afivo modules, so that a user does not have to include them separately.
Definition: m_af_all.f90:3