1 #include "../src/cpp_macros.h"
5 program dielectric_surface
10 integer,
parameter :: box_size = 8
11 integer,
parameter :: n_iterations = 10
16 integer :: i_fld_norm_fc
17 integer :: i_fld_cc(NDIM)
19 real(dp) :: fac = 1.0_dp
22 double precision,
parameter :: epsilon_high = 2.0_dp
25 real(dp),
parameter :: interface_location = 0.25_dp
27 integer,
parameter :: interface_dimension = 1
30 type(ref_info_t) :: ref_info
31 type(surfaces_t) :: dielectric
34 character(len=100) :: fname
36 integer,
allocatable :: ref_links(:, :)
38 call af_add_cc_variable(tree,
"phi", ix=i_phi)
39 call af_add_cc_variable(tree,
"rhs", ix=i_rhs)
40 call af_add_cc_variable(tree,
"tmp", ix=i_tmp)
41 call af_add_cc_variable(tree,
"eps", ix=i_eps)
44 call af_add_cc_variable(tree,
"fld_cc_x", ix=i_fld_cc(1))
46 call af_add_cc_variable(tree,
"fld_cc_y", ix=i_fld_cc(2))
49 call af_add_cc_variable(tree,
"fld_cc_z", ix=i_fld_cc(3))
51 call af_add_cc_variable(tree,
"fld_norm_fc", ix=i_fld_norm_fc)
52 call af_add_fc_variable(tree,
"fld_fc", i_fld_fc)
54 call af_set_cc_methods(tree, i_eps, af_bc_neumann_zero, &
55 af_gc_prolong_copy, prolong=af_prolong_zeroth)
56 call af_set_cc_methods(tree, i_rhs, af_bc_neumann_zero)
65 call af_adjust_refinement(tree, ref_routine, ref_info)
66 if (ref_info%n_add == 0)
exit
69 call af_loop_box(tree, set_init_cond)
71 call surface_initialize(tree, i_eps, dielectric, 1)
73 call surface_set_values(tree, dielectric, 1, sigma_function)
76 call surface_get_refinement_links(dielectric, ref_links)
77 call af_adjust_refinement(tree, ref_random, ref_info, ref_links=ref_links)
78 call surface_update_after_refinement(tree, dielectric, ref_info)
81 call surface_surface_charge_to_rhs(tree, dielectric, 1, i_rhs, fac)
89 call mg_init(tree, mg)
91 do mg_iter = 1, n_iterations
92 call mg_fas_fmg(tree, mg, .true., mg_iter>1)
93 call af_loop_box(tree, compute_fields)
94 call surface_correct_field_fc(tree, dielectric, 1, i_fld_fc, i_phi, -fac)
95 call af_loop_box(tree, compute_field_norms)
98 call af_tree_maxabs_cc(tree, i_tmp, residu)
99 write(*,
"(I8,Es14.5)") mg_iter, residu
101 write(fname,
"(A,I0)")
"output/dielectric_surface_" // &
102 dimname //
"_", mg_iter
103 call af_write_silo(tree, trim(fname))
109 subroutine ref_routine(box, cell_flags)
110 type(box_t),
intent(in) :: box
111 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
113 if (box%lvl < 3)
then
114 cell_flags(dtimes(:)) = af_do_ref
116 cell_flags(dtimes(:)) = af_keep_ref
118 end subroutine ref_routine
120 subroutine ref_random(box, cell_flags)
121 type(box_t),
intent(in) :: box
122 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
126 call random_number(rr)
128 if (rr < 0.5_dp**ndim .and. box%lvl < 6)
then
129 cell_flags = af_do_ref
130 else if (box%lvl > 3)
then
131 cell_flags = af_rm_ref
133 cell_flags = af_keep_ref
135 end subroutine ref_random
138 subroutine set_init_cond(box)
139 type(box_t),
intent(inout) :: box
146 rr = af_r_cc(box, [ijk])
149 if (rr(interface_dimension) < interface_location)
then
150 box%cc(ijk, i_eps) = epsilon_high
152 box%cc(ijk, i_eps) = 1.0_dp
159 end subroutine set_init_cond
162 subroutine bc_phi(box, nb, iv, coords, bc_val, bc_type)
163 type(box_t),
intent(in) :: box
164 integer,
intent(in) :: nb
165 integer,
intent(in) :: iv
166 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
167 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
168 integer,
intent(out) :: bc_type
175 case (2 * interface_dimension - 1)
176 bc_type = af_bc_dirichlet
178 case (2 * interface_dimension)
179 bc_type = af_bc_dirichlet
182 bc_type = af_bc_neumann
185 end subroutine bc_phi
187 real(dp) function sigma_function(r)
188 real(dp),
intent(in) :: r(NDIM)
200 sigma_function = 0.0_dp
201 end function sigma_function
203 subroutine compute_fields(box)
204 type(box_t),
intent(inout) :: box
206 real(dp) :: inv_dr(NDIM)
212 box%fc(1:nc+1, 1:nc, 1, i_fld_fc) = inv_dr(1) * &
213 (box%cc(0:nc, 1:nc, i_phi) - box%cc(1:nc+1, 1:nc, i_phi))
214 box%fc(1:nc, 1:nc+1, 2, i_fld_fc) = inv_dr(2) * &
215 (box%cc(1:nc, 0:nc, i_phi) - box%cc(1:nc, 1:nc+1, i_phi))
216 box%cc(1:nc, 1:nc, i_fld_cc(1)) = 0.5_dp * inv_dr(1) * &
217 (box%cc(0:nc-1, 1:nc, i_phi) - box%cc(2:nc+1, 1:nc, i_phi))
218 box%cc(1:nc, 1:nc, i_fld_cc(2)) = 0.5_dp * inv_dr(2) * &
219 (box%cc(1:nc, 0:nc-1, i_phi) - box%cc(1:nc, 2:nc+1, i_phi))
221 box%fc(1:nc+1, 1:nc, 1:nc, 1, i_fld_fc) = inv_dr(1) * &
222 (box%cc(0:nc, 1:nc, 1:nc, i_phi) - &
223 box%cc(1:nc+1, 1:nc, 1:nc, i_phi))
224 box%fc(1:nc, 1:nc+1, 1:nc, 2, i_fld_fc) = inv_dr(2) * &
225 (box%cc(1:nc, 0:nc, 1:nc, i_phi) - &
226 box%cc(1:nc, 1:nc+1, 1:nc, i_phi))
227 box%fc(1:nc, 1:nc, 1:nc+1, 3, i_fld_fc) = inv_dr(3) * &
228 (box%cc(1:nc, 1:nc, 0:nc, i_phi) - &
229 box%cc(1:nc, 1:nc, 1:nc+1, i_phi))
232 end subroutine compute_fields
234 subroutine compute_field_norms(box)
235 type(box_t),
intent(inout) :: box
241 box%cc(1:nc, 1:nc, i_fld_norm_fc) = 0.5_dp * sqrt(&
242 (box%fc(1:nc, 1:nc, 1, i_fld_fc) + &
243 box%fc(2:nc+1, 1:nc, 1, i_fld_fc))**2 + &
244 (box%fc(1:nc, 1:nc, 2, i_fld_fc) + &
245 box%fc(1:nc, 2:nc+1, 2, i_fld_fc))**2)
247 box%cc(1:nc, 1:nc, 1:nc, i_fld_norm_fc) = 0.5_dp * sqrt(&
248 (box%fc(1:nc, 1:nc, 1:nc, 1, i_fld_fc) + &
249 box%fc(2:nc+1, 1:nc, 1:nc, 1, i_fld_fc))**2 + &
250 (box%fc(1:nc, 1:nc, 1:nc, 2, i_fld_fc) + &
251 box%fc(1:nc, 2:nc+1, 1:nc, 2, i_fld_fc))**2 + &
252 (box%fc(1:nc, 1:nc, 1:nc, 3, i_fld_fc) + &
253 box%fc(1:nc, 1:nc, 2:nc+1, 3, i_fld_fc))**2)
255 end subroutine compute_field_norms
257 end program dielectric_surface
Module which contains all Afivo modules, so that a user does not have to include them separately.