1 #include "../src/cpp_macros.h"
5 program electrode_dielectric
10 integer,
parameter :: box_size = 8
11 integer,
parameter :: n_iterations = 10
18 integer :: i_field_norm
21 type(ref_info_t) :: ref_info
22 integer :: n, mg_iter, coord, n_args
24 character(len=100) :: fname
27 if (ndim /= 2) error stop
"Example only set up for 2D"
29 call af_add_cc_variable(tree,
"phi", ix=i_phi)
30 call af_add_cc_variable(tree,
"rhs", ix=i_rhs)
31 call af_add_cc_variable(tree,
"tmp", ix=i_tmp)
32 call af_add_cc_variable(tree,
"lsf", ix=i_lsf)
33 call af_add_cc_variable(tree,
"eps", ix=i_eps)
34 call af_add_cc_variable(tree,
"field_norm", ix=i_field_norm)
35 call af_add_fc_variable(tree,
"field", ix=i_field)
37 call af_set_cc_methods(tree, i_lsf, funcval=set_lsf)
38 call af_set_cc_methods(tree, i_rhs, af_bc_neumann_zero)
41 n_args = command_argument_count()
42 if (ndim == 2 .and. n_args == 1)
then
52 4 * [dtimes(box_size)], &
56 call af_adjust_refinement(tree, ref_routine, ref_info)
57 if (ref_info%n_add == 0)
exit
65 mg%sides_bc => af_bc_dirichlet_zero
67 mg%lsf_boundary_value = 1.0_dp
69 call mg_init(tree, mg)
71 do mg_iter = 1, n_iterations
72 call mg_fas_fmg(tree, mg, .true., mg_iter>1)
73 call mg_compute_phi_gradient(tree, mg, i_field, 1.0_dp, i_field_norm)
76 call af_tree_maxabs_cc(tree, i_tmp, residu)
77 write(*,
"(I8,Es14.5)") mg_iter, residu
79 write(fname,
"(A,I0)")
"output/electrode_dielectric_" // &
80 dimname //
"_", mg_iter
81 call af_write_silo(tree, trim(fname))
87 subroutine ref_routine(box, cell_flags)
88 type(box_t),
intent(in) :: box
89 integer,
intent(out) :: cell_flags(DTIMES(box%n_cell))
91 if (box%lvl < 9 - 2 * ndim .and. box%r_min(2) < 0.5_dp)
then
92 cell_flags(dtimes(:)) = af_do_ref
94 cell_flags(dtimes(:)) = af_keep_ref
96 end subroutine ref_routine
99 subroutine set_lsf(box, iv)
100 type(box_t),
intent(inout) :: box
101 integer,
intent(in) :: iv
108 rr = af_r_cc(box, [ijk])
109 box%cc(ijk, iv) = get_lsf(rr)
110 if (rr(1) < 0.5_dp)
then
111 box%cc(ijk, i_eps) = 20.0_dp
113 box%cc(ijk, i_eps) = 1.0_dp
116 end subroutine set_lsf
118 real(dp) function get_lsf(r)
119 real(dp),
intent(in) :: r(NDIM)
120 real(dp) :: r0(NDIM), r1(NDIM)
127 get_lsf = dist_vec_line(r, r0, r1, ndim) - radius
132 function dist_vec_line(r, r0, r1, n_dim)
result(dist)
133 integer,
intent(in) :: n_dim
134 real(dp),
intent(in) :: r(n_dim), r0(n_dim), r1(n_dim)
135 real(dp) :: dist_vec(n_dim), frac
136 real(dp) :: dist, line_len2
138 line_len2 = sum((r1 - r0)**2)
139 frac = sum((r - r0) * (r1 - r0))
141 if (frac <= 0.0_dp)
then
143 else if (frac >= line_len2)
then
146 dist_vec = r - (r0 + frac/line_len2 * (r1 - r0))
149 dist = norm2(dist_vec)
150 end function dist_vec_line
152 end program electrode_dielectric
Module which contains all Afivo modules, so that a user does not have to include them separately.