Afivo  0.3
electrode_example.f90
1 #include "../src/cpp_macros.h"
2 !> \example dielectric_surface.f90
3 !>
4 !> Example showing how to include a dielectric surface
5 program dielectric_surface
6  use m_af_all
7 
8  implicit none
9 
10  integer, parameter :: box_size = 8
11  integer, parameter :: n_iterations = 10
12  integer :: i_phi
13  integer :: i_rhs
14  integer :: i_tmp
15  integer :: i_lsf
16  integer :: i_field
17  integer :: i_field_norm
18 
19  type(af_t) :: tree
20  type(ref_info_t) :: ref_info
21  integer :: n, mg_iter, coord, n_args
22  real(dp) :: residu
23  character(len=100) :: fname
24  type(mg_t) :: mg
25 
26  call af_add_cc_variable(tree, "phi", ix=i_phi)
27  call af_add_cc_variable(tree, "rhs", ix=i_rhs)
28  call af_add_cc_variable(tree, "tmp", ix=i_tmp)
29  call af_add_cc_variable(tree, "lsf", ix=i_lsf)
30  call af_add_cc_variable(tree, "field_norm", ix=i_field_norm)
31  call af_add_fc_variable(tree, "field", ix=i_field)
32 
33  call af_set_cc_methods(tree, i_lsf, funcval=set_lsf)
34  call af_set_cc_methods(tree, i_rhs, af_bc_neumann_zero)
35 
36  ! If an argument is given, switch to cylindrical coordinates in 2D
37  n_args = command_argument_count()
38  if (ndim == 2 .and. n_args == 1) then
39  coord = af_cyl
40  else
41  coord = af_xyz
42  end if
43 
44  ! Initialize tree
45  call af_init(tree, & ! Tree to initialize
46  box_size, & ! A box contains box_size**DIM cells
47  [dtimes(1.0_dp)], &
48  4 * [dtimes(box_size)], &
49  coord=coord)
50 
51  do n = 1, 4
52  call af_adjust_refinement(tree, ref_routine, ref_info)
53  if (ref_info%n_add == 0) exit
54  end do
55 
56  tree%mg_i_lsf = i_lsf
57  mg%i_phi = i_phi
58  mg%i_rhs = i_rhs
59  mg%i_tmp = i_tmp
60  mg%sides_bc => af_bc_dirichlet_zero
61  mg%lsf_boundary_value = 1.0_dp
62  mg%lsf => get_lsf
63 
64  call mg_init(tree, mg)
65 
66  do mg_iter = 1, n_iterations
67  call mg_fas_fmg(tree, mg, .true., mg_iter>1)
68  call mg_compute_phi_gradient(tree, mg, i_field, 1.0_dp, i_field_norm)
69 
70  ! Determine the minimum and maximum residual and error
71  call af_tree_maxabs_cc(tree, i_tmp, residu)
72  write(*, "(I8,Es14.5)") mg_iter, residu
73 
74  write(fname, "(A,I0)") "output/electrode_example_" // dimname // "_", mg_iter
75  call af_write_silo(tree, trim(fname))
76  end do
77 
78 contains
79 
80  ! Set refinement flags for box
81  subroutine ref_routine(box, cell_flags)
82  type(box_t), intent(in) :: box
83  integer, intent(out) :: cell_flags(DTIMES(box%n_cell))
84 
85  if (box%lvl < 9 - 2 * ndim .and. box%r_min(1) < 0.5_dp) then
86  cell_flags(dtimes(:)) = af_do_ref
87  else
88  cell_flags(dtimes(:)) = af_keep_ref
89  end if
90  end subroutine ref_routine
91 
92  ! This routine sets the level set function
93  subroutine set_lsf(box, iv)
94  type(box_t), intent(inout) :: box
95  integer, intent(in) :: iv
96  integer :: IJK, nc
97  real(dp) :: rr(NDIM)
98 
99  nc = box%n_cell
100 
101  do kji_do(0,nc+1)
102  rr = af_r_cc(box, [ijk])
103  box%cc(ijk, iv) = get_lsf(rr)
104  end do; close_do
105 
106  end subroutine set_lsf
107 
108  real(dp) function get_lsf(r)
109  real(dp), intent(in) :: r(NDIM)
110  real(dp) :: r0(NDIM), r1(NDIM)
111  real(dp) :: radius
112 
113  ! Start and end point of line
114  r0(:) = 0.4_dp
115  r1(:) = 0.6_dp
116  radius = 0.02_dp
117 
118  get_lsf = dist_vec_line(r, r0, r1, ndim) - radius
119  end function get_lsf
120 
121  !> Compute distance vector between point and its projection onto a line
122  !> between r0 and r1
123  function dist_vec_line(r, r0, r1, n_dim) result(dist)
124  integer, intent(in) :: n_dim
125  real(dp), intent(in) :: r(n_dim), r0(n_dim), r1(n_dim)
126  real(dp) :: dist_vec(n_dim), frac
127  real(dp) :: dist, line_len2
128 
129  line_len2 = sum((r1 - r0)**2)
130  frac = sum((r - r0) * (r1 - r0))
131 
132  if (frac <= 0.0_dp) then
133  dist_vec = r - r0
134  else if (frac >= line_len2) then
135  dist_vec = r - r1
136  else
137  dist_vec = r - (r0 + frac/line_len2 * (r1 - r0))
138  end if
139 
140  dist = norm2(dist_vec)
141  end function dist_vec_line
142 
143 end program dielectric_surface
Module which contains all Afivo modules, so that a user does not have to include them separately.
Definition: m_af_all.f90:3