Afivo  0.3
dielectric_2d.f90
1 
6 program dielectric_test
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_eps = 4
20 
21  ! The dielectric constant used in this example
22  double precision, parameter :: epsilon_high = 10.0_dp
23 
24  type(a2_t) :: tree
25  type(ref_info_t) :: ref_info
26  integer :: mg_iter
27  integer :: ix_list(2, n_boxes_base)
28  real(dp) :: dr, residu(2)
29  character(len=100) :: fname
30  type(mg2_t) :: mg
31  integer :: count_rate, t_start, t_end
32 
33  print *, "****************************************"
34  print *, "Warning: functionality demonstrated here is not fully ready"
35  print *, "For large epsilon, convergence will probably be slow"
36  print *, "****************************************"
37  print *, "Number of threads", af_get_max_threads()
38 
39  ! The cell spacing at the coarsest grid level
40  dr = 1.0_dp / box_size
41 
42  ! Initialize tree
43  call a2_init(tree, & ! Tree to initialize
44  box_size, & ! A box contains box_size**DIM cells
45  n_var_cell, & ! Number of cell-centered variables
46  0, & ! Number of face-centered variables
47  dr, & ! Distance between cells on base level
48  coarsen_to=2, & ! Add coarsened levels for multigrid
49  cc_names=["phi", "rhs", "tmp", "eps"]) ! Variable names
50 
51  ! Set up geometry. These indices are used to define the coordinates of a box,
52  ! by default the box at [1,1] touches the origin (x,y) = (0,0)
53  ix_list(:, 1) = 1 ! Set index of box 1
54 
55  ! Create the base mesh, using the box indices and their neighbor information
56  call a2_set_base(tree, 1, ix_list)
57  call a2_print_info(tree)
58 
59  call system_clock(t_start, count_rate)
60  do
61  ! For each box, set the initial conditions
62  call a2_loop_box(tree, set_init_cond)
63 
64  ! This updates the refinement of the tree, by at most one level per call.
65  call a2_adjust_refinement(tree, ref_routine, ref_info)
66 
67  ! If no new boxes have been added, exit the loop
68  if (ref_info%n_add == 0) exit
69  end do
70  call system_clock(t_end, count_rate)
71 
72  ! Average epsilon on coarse grids. In the future, it could be better to define
73  ! epsilon on cell faces, and to perform this restriction in a matrix fashion:
74  ! A_coarse = M_restrict * A_fine * M_prolong (A = matrix operator, M = matrix)
75  call a2_restrict_tree(tree, i_eps)
76 
77  write(*,"(A,Es10.3,A)") " Wall-clock time generating AMR grid: ", &
78  (t_end-t_start) / real(count_rate,dp), " seconds"
79 
80  call a2_print_info(tree)
81 
82  ! Set the multigrid options.
83  mg%i_phi = i_phi ! Solution variable
84  mg%i_rhs = i_rhs ! Right-hand side variable
85  mg%i_tmp = i_tmp ! Variable for temporary space
86  mg%i_eps = i_eps ! Variable for epsilon coefficient
87  mg%sides_bc => sides_bc
88 
89  ! Automatically detect the right methods
90  mg%box_op => mg2_auto_op
91  mg%box_gsrb => mg2_auto_gsrb
92  mg%box_corr => mg2_auto_corr
93 
94  ! Initialize the multigrid options. This performs some basics checks and sets
95  ! default values where necessary.
96  ! This routine does not initialize the multigrid variables i_phi, i_rhs
97  ! and i_tmp. These variables will be initialized at the first call of mg2_fas_fmg
98  call mg2_init_mg(mg)
99 
100  print *, "Multigrid iteration | max residual | max error"
101  call system_clock(t_start, count_rate)
102 
103  do mg_iter = 1, n_iterations
104  ! Perform a FAS-FMG cycle (full approximation scheme, full multigrid). The
105  ! third argument controls whether the residual is stored in i_tmp. The
106  ! fourth argument controls whether to improve the current solution.
107  call mg2_fas_fmg(tree, mg, .true., mg_iter>1)
108 
109  ! Determine the minimum and maximum residual and error
110  call a2_tree_min_cc(tree, i_tmp, residu(1))
111  call a2_tree_max_cc(tree, i_tmp, residu(2))
112  write(*,"(I8,Es14.5)") mg_iter, maxval(abs(residu))
113 
114  write(fname, "(A,I0)") "dielectric_2d_", mg_iter
115  call a2_write_silo(tree, trim(fname), dir="output")
116  end do
117  call system_clock(t_end, count_rate)
118 
119  write(*, "(A,I0,A,E10.3,A)") &
120  " Wall-clock time after ", n_iterations, &
121  " iterations: ", (t_end-t_start) / real(count_rate, dp), &
122  " seconds"
123 
124  ! This call is not really necessary here, but cleaning up the data in a tree
125  ! is important if your program continues with other tasks.
126  call a2_destroy(tree)
127 
128 contains
129 
130  ! Set refinement flags for box
131  subroutine ref_routine(box, cell_flags)
132  type(box2_t), intent(in) :: box
133  integer, intent(out) :: cell_flags(box%n_cell, box%n_cell)
134  real(dp) :: eps_min, eps_max
135 
136  eps_min = minval(box%cc(:, :, i_eps))
137  eps_max = maxval(box%cc(:, :, i_eps))
138 
139  if ((box%lvl < 7 .and. eps_max > eps_min) .or. box%lvl < 3) then
140  cell_flags(:, :) = af_do_ref
141  else
142  cell_flags(:, :) = af_keep_ref
143  end if
144  end subroutine ref_routine
145 
146  ! This routine sets the initial conditions for each box
147  subroutine set_init_cond(box)
148  type(box2_t), intent(inout) :: box
149  integer :: i, j, nc
150  real(dp) :: rr(2)
151  real(dp) :: ellips_fac(2)
152 
153  nc = box%n_cell
154  dr = box%dr
155 
156  ! Create ellipsoidal shape
157  ellips_fac(2:) = 3.0_dp
158  ellips_fac(1) = 1.0_dp
159 
160  do j = 0, nc+1; do i = 0, nc+1
161  rr = a2_r_cc(box, [i, j])
162 
163  ! Change epsilon in part of the domain
164  if (norm2((rr - 0.5_dp) * ellips_fac) < 0.25_dp) then
165  box%cc(i, j, i_eps) = epsilon_high
166  else
167  box%cc(i, j, i_eps) = 1.0_dp
168  end if
169 
170  box%cc(i, j, i_rhs) = 0.0d0
171  box%cc(i, j, i_phi) = 0.0d0
172  end do; end do
173 
174  end subroutine set_init_cond
175 
176  subroutine sides_bc(box, nb, iv, bc_type)
177  use m_a2_ghostcell
178  type(box2_t), intent(inout) :: box
179  integer, intent(in) :: nb ! Direction for the boundary condition
180  integer, intent(in) :: iv ! Index of variable
181  integer, intent(out) :: bc_type ! Type of boundary condition
182  integer :: nc
183 
184  nc = box%n_cell
185 
186  select case (nb)
187  case (a2_neighb_lowx)
188  call a2_bc_dirichlet_zero(box, nb, iv, bc_type)
189  case (a2_neighb_highx)
190  bc_type = af_bc_dirichlet
191  box%cc(nc+1, 1:nc, iv) = 1.0_dp
192  case default
193  call a2_bc_neumann_zero(box, nb, iv, bc_type)
194  end select
195  end subroutine sides_bc
196 
197 end program dielectric_test
This module can be used to construct solutions consisting of one or more Gaussians.
Definition: m_gaussians.f90:3
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
Module which contains all Afivo modules, so that a user does not have to include them separately...
Definition: m_a2_all.f90:4