1 #include "../src/cpp_macros.h"
5 program dielectric_surface
10 integer,
parameter :: box_size = 8
11 integer,
parameter :: n_iterations = 10
17 integer :: i_field_norm
20 type(ref_info_t) :: ref_info
21 integer :: n, mg_iter, coord, n_args
23 character(len=100) :: fname
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)
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)
37 n_args = command_argument_count()
38 if (ndim == 2 .and. n_args == 1)
then
48 4 * [dtimes(box_size)], &
52 call af_adjust_refinement(tree, ref_routine, ref_info)
53 if (ref_info%n_add == 0)
exit
60 mg%sides_bc => af_bc_dirichlet_zero
61 mg%lsf_boundary_value = 1.0_dp
64 call mg_init(tree, mg)
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)
71 call af_tree_maxabs_cc(tree, i_tmp, residu)
72 write(*,
"(I8,Es14.5)") mg_iter, residu
74 write(fname,
"(A,I0)")
"output/electrode_example_" // dimname //
"_", mg_iter
75 call af_write_silo(tree, trim(fname))
81 subroutine ref_routine(box, cell_flags)
82 type(box_t),
intent(in) :: box
83 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
85 if (box%lvl < 9 - 2 * ndim .and. box%r_min(1) < 0.5_dp)
then
86 cell_flags(dtimes(:)) = af_do_ref
88 cell_flags(dtimes(:)) = af_keep_ref
90 end subroutine ref_routine
93 subroutine set_lsf(box, iv)
94 type(box_t),
intent(inout) :: box
95 integer,
intent(in) :: iv
102 rr = af_r_cc(box, [ijk])
103 box%cc(ijk, iv) = get_lsf(rr)
106 end subroutine set_lsf
108 real(dp) function get_lsf(r)
109 real(dp),
intent(in) :: r(NDIM)
110 real(dp) :: r0(NDIM), r1(NDIM)
118 get_lsf = dist_vec_line(r, r0, r1, ndim) - radius
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
129 line_len2 = sum((r1 - r0)**2)
130 frac = sum((r - r0) * (r1 - r0))
132 if (frac <= 0.0_dp)
then
134 else if (frac >= line_len2)
then
137 dist_vec = r - (r0 + frac/line_len2 * (r1 - r0))
140 dist = norm2(dist_vec)
141 end function dist_vec_line
143 end program dielectric_surface
Module which contains all Afivo modules, so that a user does not have to include them separately.