6 program dielectric_test
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
22 double precision,
parameter :: epsilon_high = 10.0_dp
25 type(ref_info_t) :: ref_info
27 integer :: ix_list(3, n_boxes_base)
28 real(dp) :: dr, residu(2)
29 character(len=100) :: fname
31 integer :: count_rate, t_start, t_end
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()
40 dr = 1.0_dp / box_size
49 cc_names=[
"phi",
"rhs",
"tmp",
"eps"])
56 call a3_set_base(tree, 1, ix_list)
57 call a3_print_info(tree)
59 call system_clock(t_start, count_rate)
62 call a3_loop_box(tree, set_init_cond)
65 call a3_adjust_refinement(tree, ref_routine, ref_info)
68 if (ref_info%n_add == 0)
exit 70 call system_clock(t_end, count_rate)
75 call a3_restrict_tree(tree, i_eps)
77 write(*,
"(A,Es10.3,A)")
" Wall-clock time generating AMR grid: ", &
78 (t_end-t_start) /
real(count_rate,dp),
" seconds" 80 call a3_print_info(tree)
87 mg%sides_bc => sides_bc
90 mg%box_op => mg3_auto_op
91 mg%box_gsrb => mg3_auto_gsrb
92 mg%box_corr => mg3_auto_corr
100 print *,
"Multigrid iteration | max residual | max error" 101 call system_clock(t_start, count_rate)
103 do mg_iter = 1, n_iterations
107 call mg3_fas_fmg(tree, mg, .true., mg_iter>1)
110 call a3_tree_min_cc(tree, i_tmp, residu(1))
111 call a3_tree_max_cc(tree, i_tmp, residu(2))
112 write(*,
"(I8,Es14.5)") mg_iter, maxval(abs(residu))
114 write(fname,
"(A,I0)")
"dielectric_3d_", mg_iter
115 call a3_write_silo(tree, trim(fname), dir=
"output")
117 call system_clock(t_end, count_rate)
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), &
126 call a3_destroy(tree)
131 subroutine ref_routine(box, cell_flags)
132 type(box3_t),
intent(in) :: box
133 integer,
intent(out) :: cell_flags(box%n_cell, box%n_cell, box%n_cell)
134 real(dp) :: eps_min, eps_max
136 eps_min = minval(box%cc(:, :, :, i_eps))
137 eps_max = maxval(box%cc(:, :, :, i_eps))
139 if ((box%lvl < 7 .and. eps_max > eps_min) .or. box%lvl < 3)
then 140 cell_flags(:, :, :) = af_do_ref
142 cell_flags(:, :, :) = af_keep_ref
144 end subroutine ref_routine
147 subroutine set_init_cond(box)
148 type(box3_t),
intent(inout) :: box
149 integer :: i, j, k, nc
151 real(dp) :: ellips_fac(3)
157 ellips_fac(2:) = 3.0_dp
158 ellips_fac(1) = 1.0_dp
160 do k = 0, nc+1;
do j = 0, nc+1;
do i = 0, nc+1
161 rr = a3_r_cc(box, [i, j, k])
164 if (norm2((rr - 0.5_dp) * ellips_fac) < 0.25_dp)
then 165 box%cc(i, j, k, i_eps) = epsilon_high
167 box%cc(i, j, k, i_eps) = 1.0_dp
170 box%cc(i, j, k, i_rhs) = 0.0d0
171 box%cc(i, j, k, i_phi) = 0.0d0
172 end do; end do; end do
174 end subroutine set_init_cond
176 subroutine sides_bc(box, nb, iv, bc_type)
178 type(box3_t),
intent(inout) :: box
179 integer,
intent(in) :: nb
180 integer,
intent(in) :: iv
181 integer,
intent(out) :: bc_type
187 case (a3_neighb_lowx)
188 call a3_bc_dirichlet_zero(box, nb, iv, bc_type)
189 case (a3_neighb_highx)
190 bc_type = af_bc_dirichlet
191 box%cc(nc+1, 1:nc, 1:nc, iv) = 1.0_dp
193 call a3_bc_neumann_zero(box, nb, iv, bc_type)
195 end subroutine sides_bc
197 end program dielectric_test
Module which contains all Afivo modules, so that a user does not have to include them separately...
This module can be used to construct solutions consisting of one or more Gaussians.
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...