Afivo  0.3
poisson_neumann_2d.f90
1 
6 program poisson_neumann_2d
7  use m_a2_all
8  use m_gaussians
9 
10  implicit none
11 
12  integer, parameter :: box_size = 8
13  integer, parameter :: n_boxes_base = 1
14  integer, parameter :: n_iterations = 10
15  integer, parameter :: n_var_cell = 4
16  integer, parameter :: i_phi = 1
17  integer, parameter :: i_rhs = 2
18  integer, parameter :: i_tmp = 3
19  integer, parameter :: i_err = 4
20 
21  type(a2_t) :: tree
22  type(ref_info_t) :: refine_info
23  integer :: mg_iter
24  integer :: ix_list(2, n_boxes_base)
25  real(dp) :: dr
26  character(len=100) :: fname
27  type(mg2_t) :: mg
28  integer :: count_rate,t_start,t_end
29 
30  print *, "Running poisson_neumann_2d"
31  print *, "Number of threads", af_get_max_threads()
32 
33  ! The cell spacing at the coarsest grid level
34  dr = 1.0_dp / box_size
35 
36  ! Initialize tree
37  call a2_init(tree, & ! Tree to initialize
38  box_size, & ! A box contains box_size**DIM cells
39  n_var_cell, & ! Number of cell-centered variables
40  0, & ! Number of face-centered variables
41  dr, & ! Distance between cells on base level
42  coarsen_to=2, & ! Add coarsened levels for multigrid
43  cc_names=["phi", "rhs", "tmp", "err"]) ! Variable names
44 
45  ! Set up geometry. These indices are used to define the coordinates of a box,
46  ! by default the box at [1,1] touches the origin (x,y) = (0,0)
47  ix_list(:, 1) = [1, 1] ! Set index of box 1
48 
49  ! Create the base mesh, using the box indices and their neighbor information
50  call a2_set_base(tree, 1, ix_list)
51 
52  call system_clock(t_start, count_rate)
53  do
54  ! For each box, set the initial conditions
55  call a2_loop_box(tree, set_initial_condition)
56 
57  ! This updates the refinement of the tree, by at most one level per call.
58  ! The second argument is a subroutine that is called for each box that can
59  ! be refined or derefined, and it should set refinement flags. Information
60  ! about the changes in refinement are returned in the third argument.
61  call a2_adjust_refinement(tree, refine_routine, refine_info, 0)
62 
63  ! If no new boxes have been added, exit the loop
64  if (refine_info%n_add == 0) exit
65  end do
66  call system_clock(t_end, count_rate)
67 
68  write(*,"(A,Es10.3,A)") " Wall-clock time generating AMR grid: ", &
69  (t_end-t_start) / real(count_rate,dp), " seconds"
70 
71  call a2_print_info(tree)
72 
73  ! Initialize the multigrid options. This performs some basics checks and sets
74  ! default values where necessary.
75  mg%i_phi = i_phi ! Solution variable
76  mg%i_rhs = i_rhs ! Right-hand side variable
77  mg%i_tmp = i_tmp ! Variable for temporary space
78  mg%sides_bc => sides_bc ! Method for boundary conditions
79 
80  ! This routine does not initialize the multigrid fields boxes%i_phi,
81  ! boxes%i_rhs and boxes%i_tmp. These fileds will be initialized at the
82  ! first call of mg2_fas_fmg
83  call mg2_init_mg(mg)
84 
85  print *, "Multigrid iteration | max residual | max error"
86  call system_clock(t_start, count_rate)
87 
88  do mg_iter = 1, n_iterations
89  call mg2_fas_fmg(tree, mg, set_residual=.true., have_guess=(mg_iter>1))
90  ! call mg2_fas_vcycle(tree, mg, tree%highest_lvl, set_residual=.true.)
91 
92  call a2_loop_box(tree, set_error)
93 
94  write(fname, "(A,I0)") "poisson_neumann_2d_", mg_iter
95  call a2_write_silo(tree, trim(fname), dir="output")
96  end do
97  call system_clock(t_end, count_rate)
98 
99  write(*, "(A,I0,A,E10.3,A)") &
100  " Wall-clock time after ", n_iterations, &
101  " iterations: ", (t_end-t_start) / real(count_rate, dp), &
102  " seconds"
103 
104 contains
105 
106  ! Return the refinement flags for box
107  subroutine refine_routine(box, cell_flags)
108  type(box2_t), intent(in) :: box
109  integer, intent(out) :: cell_flags(box%n_cell, box%n_cell)
110 
111  ! Refine around one corner
112  if (box%lvl <= 4 .and. all(box%r_min < 0.25_dp)) then
113  cell_flags(:, :) = af_do_ref
114  else
115  cell_flags(:, :) = af_keep_ref
116  end if
117  end subroutine refine_routine
118 
119  ! This routine sets the initial conditions for each box
120  subroutine set_initial_condition(box)
121  type(box2_t), intent(inout) :: box
122 
123  box%cc(:, :, i_rhs) = 0.0_dp
124  end subroutine set_initial_condition
125 
126  ! Set the error compared to the analytic solution
127  subroutine set_error(box)
128  type(box2_t), intent(inout) :: box
129  integer :: i, j, nc
130  real(dp) :: rr(2)
131 
132  nc = box%n_cell
133  do j = 1, nc; do i = 1, nc
134  rr = a2_r_cc(box, [i, j])
135  box%cc(i, j, i_err) = box%cc(i, j, i_phi) - rr(1)
136  end do; end do
137  end subroutine set_error
138 
139  ! This routine sets boundary conditions for a box, by filling its ghost cells
140  ! with approriate values.
141  subroutine sides_bc(box, nb, iv, bc_type)
142  type(box2_t), intent(inout) :: box
143  integer, intent(in) :: nb ! Direction for the boundary condition
144  integer, intent(in) :: iv ! Index of variable
145  integer, intent(out) :: bc_type ! Type of boundary condition
146  integer :: nc
147 
148  nc = box%n_cell
149 
150  ! Below the solution is specified in the approriate ghost cells
151  select case (nb)
152  case (a2_neighb_lowx) ! Lower-x direction
153  call a2_bc_dirichlet_zero(box, nb, iv, bc_type)
154  case (a2_neighb_highx) ! Higher-x direction
155  bc_type = af_bc_dirichlet
156  box%cc(nc+1, 1:nc, iv) = 1
157  case default
158  call a2_bc_neumann_zero(box, nb, iv, bc_type)
159  end select
160  end subroutine sides_bc
161 
162 end program
This module can be used to construct solutions consisting of one or more Gaussians.
Definition: m_gaussians.f90:3
Module which contains all Afivo modules, so that a user does not have to include them separately...
Definition: m_a2_all.f90:4