This example shows how to create an AMR tree, perform random refinement, write output files, and how to fill ghost cells.
6 program random_refinement_2d
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
23 write(*,
'(A)')
'program random_refinement_2d' 25 print *,
"Number of threads", af_get_max_threads()
28 dr = 2 * acos(-1.0_dp) / box_size
38 call a2_set_cc_methods(tree, 1, a2_bc_dirichlet_zero, &
39 prolong=a2_prolong_linear_cons)
43 nb_list(:, :) = af_no_box
46 do id = 1, n_boxes_base
47 ix_list(:, id) = [id, 1]
48 nb_list(a2_neighb_highy, id) = id
53 nb_list(a2_neighb_lowx, 1) = n_boxes_base
56 call a2_set_base(tree, n_boxes_base, ix_list, nb_list)
57 call a2_print_info(tree)
62 call a2_loop_box(tree, set_init_cond)
67 call a2_gc_tree(tree, i_phi, a2_gc_interp, a2_bc_dirichlet_zero)
69 call a2_tree_sum_cc(tree, i_phi, sum_phi_t0)
71 call system_clock(t_start, count_rate)
77 write(fname,
"(A,I0)")
"random_refinement_2d_", iter
78 call a2_write_silo(tree, trim(fname), dir=
"output", n_cycle=iter)
89 call a2_adjust_refinement(tree, ref_routine, ref_info, ref_buffer=0)
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
99 call system_clock(t_end, count_rate)
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), &
106 call a2_print_info(tree)
110 call a2_destroy(tree)
115 subroutine ref_routine(box, cell_flags)
116 type(box2_t),
intent(in) :: box
117 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell)
121 call random_number(rr)
123 if (rr < 0.5_dp**2 .and. box%lvl < 5)
then 124 cell_flags = af_do_ref
126 cell_flags = af_rm_ref
129 end subroutine ref_routine
132 subroutine set_init_cond(box)
133 type(box2_t),
intent(inout) :: box
138 do j = 1, nc;
do i = 1, nc
140 rr = a2_r_cc(box, [i, j])
143 box%cc(i, j, i_phi) = sin(0.5_dp * rr(1)) * cos(rr(2))
145 end subroutine set_init_cond
147 end program random_refinement_2d