Afivo  0.3
boundary_conditions_2d.f90

Example showing how to use different types of boundary conditions

1 
5 program boundary_conditions_2d
6  use m_a2_all
7 
8  implicit none
9 
10  integer, parameter :: box_size = 8
11  integer, parameter :: n_boxes_base = 1
12  integer, parameter :: n_iterations = 20
13  integer, parameter :: i_phi = 1
14 
15  type(a2_t) :: tree
16  integer :: iter
17  integer :: ix_list(2, n_boxes_base)
18  real(dp) :: dr
19  character(len=100) :: fname
20 
21  print *, "Running boundary_conditions_2d"
22  print *, "Number of threads", af_get_max_threads()
23 
24  ! The cell spacing at the coarsest grid level
25  dr = 1.0_dp / box_size
26 
27  ! Initialize tree
28  call a2_init(tree, & ! Tree to initialize
29  box_size, & ! A box contains box_size**DIM cells
30  1, & ! Number of cell-centered variables
31  0, & ! Number of face-centered variables
32  dr, & ! Distance between cells on base level
33  cc_names=["phi"]) ! Variable names
34 
35  ! Set up geometry. This creates a single box, for which we need boundary
36  ! conditions on all sides
37  ix_list(:, 1) = [1, 1] ! Set index of box 1
38 
39  ! Create the base mesh, using the box indices and their neighbor information
40  call a2_set_base(tree, 1, ix_list)
41  call a2_print_info(tree)
42 
43  do iter = 1, n_iterations
44  if (iter == 1) then
45  ! Set initial conditions on first iteration
46  call a2_loop_box(tree, set_phi_zero)
47  else
48  ! Replace phi by average over neighbors
49  call a2_loop_box(tree, average_phi)
50  end if
51 
52  call a2_gc_tree(tree, i_phi, a2_gc_interp, boundary_method)
53 
54  write(fname, "(A,I0)") "boundary_conditions_2d_", iter
55  call a2_write_vtk(tree, trim(fname), dir="output")
56  end do
57 
58 contains
59 
60  ! This routine sets the initial conditions for each box
61  subroutine set_phi_zero(box)
62  type(box2_t), intent(inout) :: box
63 
64  box%cc(:, :, i_phi) = 0.0_dp
65  end subroutine set_phi_zero
66 
67  subroutine average_phi(box)
68  type(box2_t), intent(inout) :: box
69  integer :: i, j, nc
70  real(dp) :: tmp(box%n_cell, box%n_cell)
71 
72  nc = box%n_cell
73 
74  do j = 1, nc; do i = 1, nc
75  tmp(i, j) = 0.25_dp * ( &
76  box%cc(i+1, j, i_phi) + &
77  box%cc(i-1, j, i_phi) + &
78  box%cc(i, j+1, i_phi) + &
79  box%cc(i, j-1, i_phi))
80  end do; end do
81 
82  box%cc(1:nc, 1:nc, i_phi) = tmp(:, :)
83  end subroutine average_phi
84 
86  subroutine boundary_method(box, nb, iv, bc_type)
87  type(box2_t), intent(inout) :: box ! Box to operate on
88  integer, intent(in) :: nb ! Direction for the boundary condition
89  integer, intent(in) :: iv ! Index of variable
90  integer, intent(out) :: bc_type ! Type of boundary condition
91  integer :: nc
92 
93  nc = box%n_cell
94 
95  ! Below the solution is specified in the approriate ghost cells
96  select case (nb)
97  case (a2_neighb_lowx) ! Lower-x direction
98  bc_type = af_bc_dirichlet
99  box%cc(0, 1:nc, iv) = 1.0_dp
100  case (a2_neighb_highx) ! Higher-x direction
101  bc_type = af_bc_neumann
102  box%cc(nc+1, 1:nc, iv) = 0.0_dp
103  case (a2_neighb_lowy) ! Lower-y direction
104  bc_type = af_bc_dirichlet
105  box%cc(1:nc, 0, iv) = 1.0_dp
106  case (a2_neighb_highy) ! Higher-y direction
107  bc_type = af_bc_neumann
108  box%cc(1:nc, nc+1, iv) = 0.0_dp
109  end select
110 
111  end subroutine boundary_method
113 
114 end program boundary_conditions_2d