Afivo  0.3
m_af_surface.f90
1 #include "cpp_macros.h"
2 
6  use m_af_types
7 
8  implicit none
9  private
10 
12  type surface_t
13  logical :: in_use
14  integer :: id_in
15  integer :: id_out
17  integer :: direction
18  integer :: ix_parent
19  real(dp) :: eps
20  integer :: offset_parent(ndim-1)
21  real(dp) :: dr(ndim-1)
23  real(dp), allocatable :: sd(dtimes(:))
24  end type surface_t
25 
27  integer, public, parameter :: surface_none = -1
28 
32  logical :: initialized = .false.
34  integer :: n_variables = 0
36  integer :: n_cell = -1
38  integer :: i_eps = -1
40  integer :: i_phi = -1
42  integer :: max_ix = 0
44  integer :: surface_limit = -1
46  logical :: cylindrical = .false.
47 
49  type(surface_t), allocatable :: surfaces(:)
51  integer :: n_removed = 0
53  integer, allocatable :: removed_surfaces(:)
55  integer, allocatable :: box_id_out_to_surface_ix(:)
57  integer, allocatable :: box_id_in_to_surface_ix(:)
58  end type surfaces_t
59 
60  interface
61  function value_func(x) result(my_val)
62  import
63  real(dp), intent(in) :: x(ndim)
64  real(dp) :: my_val
65  end function value_func
66  end interface
67 
68  public :: surface_t
69  public :: surfaces_t
70  public :: surface_initialize
71  public :: surface_set_values
72  public :: surface_set_weighted_sum
73  public :: surface_copy_variable
74  public :: surface_get_integral
75  public :: surface_update_after_refinement
76  public :: surface_inside_layer_to_surface
77  public :: surface_surface_charge_to_rhs
78  public :: surface_correct_field_fc
79  public :: surface_get_refinement_links
80  public :: surface_get_surface_cell
81  public :: surface_write_output
82 #if NDIM == 2
83  public :: surface_correct_field_cc
84 #endif
85 
86 contains
87 
89  subroutine surface_initialize(tree, i_eps, diel, n_variables)
90  use m_af_ghostcell, only: af_gc_tree
91  type(af_t), intent(inout) :: tree
92  integer, intent(in) :: i_eps
93  type(surfaces_t), intent(inout) :: diel
94  integer, intent(in) :: n_variables
95  integer :: nc, id, i, j, ix
96  integer :: nb, nb_id
97  real(dp) :: curr_eps, nb_eps
98 
99  if (.not. tree%ready) error stop "Tree not ready"
100 
101  if (tree%highest_lvl > 1) then
102  ! Check that there are no partially refined grid levels
103  if (size(tree%lvls(tree%highest_lvl-1)%leaves) > 0) &
104  error stop "Tree not uniformly refined"
105  end if
106 
107  diel%initialized = .true.
108  diel%cylindrical = (tree%coord_t == af_cyl)
109  diel%i_eps = i_eps
110  nc = tree%n_cell
111  diel%max_ix = 0
112  diel%n_removed = 0
113  diel%n_variables = n_variables
114  diel%n_cell = tree%n_cell
115  ! Assume that at most 1/10th of boxes has a surface
116  diel%surface_limit = tree%box_limit / 10
117 
118  allocate(diel%surfaces(diel%surface_limit))
119  allocate(diel%removed_surfaces(diel%surface_limit))
120  allocate(diel%box_id_out_to_surface_ix(tree%box_limit))
121  allocate(diel%box_id_in_to_surface_ix(tree%box_limit))
122  diel%box_id_out_to_surface_ix = surface_none
123  diel%box_id_in_to_surface_ix = surface_none
124 
125  ! Construct a list of surfaces
126  do i = 1, size(tree%lvls(tree%highest_lvl)%ids)
127  id = tree%lvls(tree%highest_lvl)%ids(i)
128 
129  ! Check whether epsilon is set consistently on this box
130  if (maxval(tree%boxes(id)%cc(dtimes(1:nc), i_eps)) > &
131  minval(tree%boxes(id)%cc(dtimes(1:nc), i_eps))) &
132  error stop "epsilon not uniform on a box"
133 
134  ! Get eps on the interior of the current box
135  curr_eps = tree%boxes(id)%cc(dtimes(1), i_eps)
136 
137  do nb = 1, af_num_neighbors
138  nb_id = tree%boxes(id)%neighbors(nb)
139 
140  ! Skip physical boundaries
141  if (nb_id <= af_no_box) cycle
142 
143  ! Get eps on the interior of the neighbour box
144  nb_eps = tree%boxes(nb_id)%cc(dtimes(1), i_eps)
145 
146  if (nb_eps > curr_eps) then
147  ! There is a surface
148  ix = get_new_surface_ix(diel)
149 
150  diel%surfaces(ix)%in_use = .true.
151  diel%surfaces(ix)%ix_parent = surface_none
152  diel%surfaces(ix)%offset_parent(:) = 0
153  diel%surfaces(ix)%id_in = nb_id
154  diel%surfaces(ix)%id_out = id
155  diel%surfaces(ix)%direction = nb
156  diel%surfaces(ix)%eps = nb_eps
157 
158  ! Extract the grid spacing parallel to the surface
159  diel%surfaces(ix)%dr = &
160  pack(tree%boxes(id)%dr, [(j, j=1,ndim)] /= af_neighb_dim(nb))
161 
162  if (diel%box_id_out_to_surface_ix(id) /= surface_none) &
163  error stop "box with multiple adjacent surfaces"
164  if (diel%box_id_in_to_surface_ix(nb_id) /= surface_none) &
165  error stop "box with multiple adjacent surfaces"
166  diel%box_id_out_to_surface_ix(id) = ix
167  diel%box_id_in_to_surface_ix(nb_id) = ix
168  end if
169  end do
170  end do
171 
172  end subroutine surface_initialize
173 
175  function get_new_surface_ix(diel) result(ix)
176  type(surfaces_t), intent(inout) :: diel
177  logical :: use_removed
178  integer :: ix, nc
179 
180  !$omp critical
181  use_removed = (diel%n_removed > 0)
182  if (use_removed) then
183  ix = diel%removed_surfaces(diel%n_removed)
184  diel%n_removed = diel%n_removed - 1
185  else
186  diel%max_ix = diel%max_ix + 1
187  ix = diel%max_ix
188  end if
189  !$omp end critical
190 
191  if (.not. use_removed) then
192  nc = diel%n_cell
193 #if NDIM == 1
194  allocate(diel%surfaces(ix)%sd(diel%n_variables))
195 #elif NDIM == 2
196  allocate(diel%surfaces(ix)%sd(nc, diel%n_variables))
197 #elif NDIM == 3
198  allocate(diel%surfaces(ix)%sd(nc, nc, diel%n_variables))
199 #endif
200  diel%surfaces(ix)%sd = 0.0_dp
201  end if
202  end function get_new_surface_ix
203 
205  subroutine surface_set_values(tree, diel, iv, user_func)
206  type(af_t), intent(in) :: tree
207  type(surfaces_t), intent(inout) :: diel
208  integer, intent(in) :: iv
209  procedure(value_func) :: user_func
210  integer :: ix, id_out
211 #if NDIM > 1
212  integer :: i
213 #endif
214 #if NDIM == 3
215  integer :: j, n
216 #endif
217  real(dp) :: coords(ndim, tree%n_cell**(ndim-1))
218 
219  if (.not. diel%initialized) error stop "dielectric not initialized"
220 
221  ! Loop over the surfaces and call the user function to set values
222  do ix = 1, diel%max_ix
223  if (diel%surfaces(ix)%in_use) then
224  id_out = diel%surfaces(ix)%id_out
225  call af_get_face_coords(tree%boxes(id_out), &
226  diel%surfaces(ix)%direction, coords)
227 #if NDIM == 1
228  diel%surfaces(ix)%sd(iv) = user_func(coords(1, 1))
229 #elif NDIM == 2
230  do i = 1, tree%n_cell
231  diel%surfaces(ix)%sd(i, iv) = user_func(coords(:, i))
232  end do
233 #elif NDIM == 3
234  n = 0
235  do j = 1, tree%n_cell
236  do i = 1, tree%n_cell
237  n = n + 1
238  diel%surfaces(ix)%sd(i, j, iv) = user_func(coords(:, n))
239  end do
240  end do
241 #endif
242  end if
243  end do
244  end subroutine surface_set_values
245 
247  subroutine surface_set_weighted_sum(diel, i_out, i_in, w_in)
248  type(surfaces_t), intent(inout) :: diel
249  integer, intent(in) :: i_out
250  integer, intent(in) :: i_in(:)
251  real(dp), intent(in) :: w_in(:)
252  integer :: ix
253 #if NDIM == 3
254  integer :: ii
255 #endif
256 
257  do ix = 1, diel%max_ix
258  if (diel%surfaces(ix)%in_use) then
259 #if NDIM == 2
260  diel%surfaces(ix)%sd(:, i_out) = &
261  matmul(diel%surfaces(ix)%sd(:, i_in), w_in)
262 #elif NDIM == 3
263  ! TODO improve by removing do-loop??
264  diel%surfaces(ix)%sd(:, :, i_out) = 0.0_dp
265  do ii=1, size(i_in)
266  diel%surfaces(ix)%sd(:, :, i_out) = diel%surfaces(ix)%sd(:, :, i_out) &
267  + w_in(ii) * diel%surfaces(ix)%sd(:, :, i_in(ii))
268  end do
269 #endif
270  end if
271  end do
272  end subroutine surface_set_weighted_sum
273 
275  subroutine surface_copy_variable(diel, i_in, i_out)
276  type(surfaces_t), intent(inout) :: diel
277  integer, intent(in) :: i_in
278  integer, intent(in) :: i_out
279  integer :: ix
280 
281  do ix = 1, diel%max_ix
282  if (diel%surfaces(ix)%in_use) then
283 #if NDIM == 2
284  diel%surfaces(ix)%sd(:, i_out) = diel%surfaces(ix)%sd(:, i_in)
285 #elif NDIM == 3
286  diel%surfaces(ix)%sd(:, :, i_out) = diel%surfaces(ix)%sd(:, :, i_in)
287 #endif
288  end if
289  end do
290  end subroutine surface_copy_variable
291 
293  subroutine surface_get_integral(tree, diel, i_surf, surf_int)
294  type(af_t), intent(in) :: tree
295  type(surfaces_t), intent(in) :: diel
296  integer, intent(in) :: i_surf
297  real(dp), intent(out) :: surf_int
298  integer :: ix
299 #if NDIM == 2
300  real(dp), parameter :: two_pi = 2 * acos(-1.0_dp)
301  real(dp) :: coords(ndim, diel%n_cell**(ndim-1))
302 #endif
303 
304  surf_int = 0
305  do ix = 1, diel%max_ix
306  if (diel%surfaces(ix)%in_use) then
307 #if NDIM == 2
308  if (diel%cylindrical) then
309  ! Multiply surface elements with 2 * pi * r
310  call af_get_face_coords(tree%boxes(diel%surfaces(ix)%id_out), &
311  diel%surfaces(ix)%direction, coords)
312  surf_int = surf_int + two_pi * diel%surfaces(ix)%dr(1) * &
313  sum(diel%surfaces(ix)%sd(:, i_surf) * coords(1, :))
314  else
315  surf_int = surf_int + product(diel%surfaces(ix)%dr) * &
316  sum(diel%surfaces(ix)%sd(:, i_surf))
317  end if
318 #elif NDIM == 3
319  surf_int = surf_int + product(diel%surfaces(ix)%dr) * &
320  sum(diel%surfaces(ix)%sd(:, :, i_surf))
321 #endif
322  end if
323  end do
324  end subroutine surface_get_integral
325 
327  subroutine surface_update_after_refinement(tree, diel, ref_info)
328  type(af_t), intent(in) :: tree
329  type(surfaces_t), intent(inout) :: diel
330  type(ref_info_t), intent(in) :: ref_info
331  integer :: lvl, i, id, p_id, ix, p_ix, nc
332 
333  if (.not. diel%initialized) error stop "dielectric not initialized"
334  nc = tree%n_cell
335 
336  ! Handle removed surfaces
337  do i = 1, size(ref_info%rm)
338  id = ref_info%rm(i)
339  ix = diel%box_id_out_to_surface_ix(id)
340  if (ix > surface_none) then
341  call restrict_surface_to_parent(tree, diel, ix)
342  end if
343  end do
344 
345  ! Add new surfaces
346  do lvl = 1, tree%highest_lvl
347  do i = 1, size(ref_info%lvls(lvl)%add)
348  id = ref_info%lvls(lvl)%add(i)
349 
350  ! Check if parent has a surface
351  p_id = tree%boxes(id)%parent
352  p_ix = diel%box_id_out_to_surface_ix(p_id)
353 
354  if (p_ix > surface_none) then
355  ! We only need to prolong once per parent surface
356  if (af_ix_to_ichild(tree%boxes(id)%ix) == 1) then
357  call prolong_surface_from_parent(tree, diel, p_ix, p_id)
358  end if
359  end if
360  end do
361  end do
362 
363  end subroutine surface_update_after_refinement
364 
366  subroutine prolong_surface_from_parent(tree, diel, p_ix, p_id)
367  type(af_t), intent(in) :: tree
368  type(surfaces_t), intent(inout) :: diel
369  integer, intent(in) :: p_ix
370  integer, intent(in) :: p_id
371  integer :: i, n, ix, ix_offset(NDIM)
372  integer :: nc, dix(NDIM-1), direction
373  integer :: i_child, id_in, id_out
374 
375  nc = tree%n_cell
376 
377  do n = 1, size(af_child_adj_nb, 1)
378  ix = get_new_surface_ix(diel)
379  direction = diel%surfaces(p_ix)%direction
380 
381  ! Select a child adjacent to the neighbor containing the dielectric
382  i_child = af_child_adj_nb(n, direction)
383  id_out = tree%boxes(p_id)%children(i_child)
384  ! The neighbor inside the dielectric should exist
385  id_in = tree%boxes(id_out)%neighbors(direction)
386  if (id_out <= af_no_box .or. id_in <= af_no_box) &
387  error stop "Error in prolongation: children do not exist"
388 
389  diel%surfaces(ix)%in_use = .true.
390  diel%surfaces(ix)%ix_parent = p_ix
391  diel%surfaces(ix)%id_in = id_in
392  diel%surfaces(ix)%id_out = id_out
393  diel%surfaces(ix)%direction = direction
394  diel%surfaces(ix)%eps = diel%surfaces(p_ix)%eps
395  diel%box_id_out_to_surface_ix(id_out) = ix
396  diel%box_id_in_to_surface_ix(id_in) = ix
397 
398  ix_offset = af_get_child_offset(tree%boxes(id_out))
399 
400  ! Compute index offset of child surface relative to parent. This is nc/2
401  ! times the child offset, but only in the dimensions parallel to the
402  ! surface.
403  dix = nc/2 * pack(af_child_dix(:, i_child), &
404  [(i, i=1,ndim)] /= af_neighb_dim(direction))
405  diel%surfaces(ix)%offset_parent(:) = dix
406 
407  ! Extract the grid spacing parallel to the surface
408  diel%surfaces(ix)%dr = pack(tree%boxes(id_out)%dr, &
409  [(i, i=1,ndim)] /= af_neighb_dim(direction))
410 
411  associate(sd => diel%surfaces(ix)%sd, &
412  sd_p => diel%surfaces(p_ix)%sd)
413 #if NDIM == 2
414  ! Copy the values from the parent
415  sd(1:nc:2, :) = sd_p(dix(1)+1:dix(1)+nc/2, :)
416  sd(2:nc:2, :) = sd(1:nc:2, :)
417 #elif NDIM == 3
418  ! Copy the values from the parent
419 !FLAG
420  sd(1:nc:2, 1:nc:2, :) = sd_p(dix(1)+1:dix(1)+nc/2, dix(2)+1:dix(2)+nc/2, :)
421  sd(2:nc:2, 1:nc:2, :) = sd(1:nc:2, 1:nc:2, :)
422  sd(:, 2:nc:2, :) = sd(:, 1:nc:2, :)
423 #endif
424  end associate
425  end do
426 
427  diel%surfaces(p_ix)%in_use = .false.
428 
429  end subroutine prolong_surface_from_parent
430 
432  subroutine restrict_surface_to_parent(tree, diel, ix)
433  type(af_t), intent(in) :: tree
434  type(surfaces_t), intent(inout) :: diel
435  integer, intent(in) :: ix
436  integer :: p_ix, nc, dix(NDIM-1), id_out, id_in
437 
438  p_ix = diel%surfaces(ix)%ix_parent
439  if (p_ix == surface_none) error stop "Too much derefinement on surface"
440  dix = diel%surfaces(ix)%offset_parent
441  nc = tree%n_cell
442 
443  associate(sd => diel%surfaces(ix)%sd, &
444  sd_p => diel%surfaces(p_ix)%sd)
445 #if NDIM == 2
446  ! Average the value on the children
447  sd_p(dix(1)+1:dix(1)+nc/2, :) = 0.5_dp * (sd(1:nc:2, :) + sd(2:nc:2, :))
448 #elif NDIM == 3
449 !FLAG
450  sd_p(dix(1)+1:dix(1)+nc/2, dix(2)+1:dix(2)+nc/2, :) = 0.25_dp * &
451  (sd(1:nc:2, 1:nc:2, :) + sd(2:nc:2, 1:nc:2, :) + sd(1:nc:2, 2:nc:2, :) + sd(2:nc:2, 1:nc:2, :) + sd(2:nc:2, 2:nc:2, :))
452 #endif
453  end associate
454 
455  !$omp critical
456  ! Remove child surface
457  diel%n_removed = diel%n_removed + 1
458  diel%removed_surfaces(diel%n_removed) = ix
459  diel%surfaces(p_ix)%in_use = .true.
460  !$omp end critical
461  id_out = diel%surfaces(ix)%id_out
462  id_in = diel%surfaces(ix)%id_in
463  diel%box_id_out_to_surface_ix(id_out) = surface_none
464  diel%box_id_in_to_surface_ix(id_in) = surface_none
465  diel%surfaces(ix)%in_use = .false.
466 
467  end subroutine restrict_surface_to_parent
468 
472  subroutine surface_get_refinement_links(diel, refinement_links)
473  type(surfaces_t), intent(in) :: diel
475  integer, allocatable, intent(inout) :: refinement_links(:, :)
476  integer :: max_ix, n, ix
477 
478  if (allocated(refinement_links)) deallocate(refinement_links)
479  max_ix = diel%max_ix
480  n = count(diel%surfaces(1:max_ix)%in_use)
481  allocate(refinement_links(2, n))
482 
483  n = 0
484  do ix = 1, max_ix
485  if (diel%surfaces(ix)%in_use) then
486  n = n + 1
487  refinement_links(:, n) = [diel%surfaces(ix)%id_out, &
488  diel%surfaces(ix)%id_in]
489  end if
490  end do
491  end subroutine surface_get_refinement_links
492 
494  subroutine surface_surface_charge_to_rhs(tree, diel, i_sigma, i_rhs, fac)
495  type(af_t), intent(inout) :: tree
496  type(surfaces_t), intent(in) :: diel
497  integer, intent(in) :: i_sigma
498  integer, intent(in) :: i_rhs
499  real(dp), intent(in) :: fac
500  integer :: n
501 
502  if (.not. diel%initialized) error stop "dielectric not initialized"
503 
504  do n = 1, diel%max_ix
505  if (.not. diel%surfaces(n)%in_use) cycle
506 
507  call surface_charge_to_rhs(tree%boxes, diel%surfaces(n), &
508  i_sigma, i_rhs, fac)
509  end do
510  end subroutine surface_surface_charge_to_rhs
511 
514  subroutine surface_charge_to_rhs(boxes, surface, i_sigma, i_rhs, fac)
515  type(box_t), intent(inout) :: boxes(:)
516  type(surface_t), intent(in) :: surface
517  integer, intent(in) :: i_sigma
518  integer, intent(in) :: i_rhs
519  real(dp), intent(in) :: fac
520  integer :: nb, nc, id_in, id_out
521  integer :: glo(NDIM), ghi(NDIM)
522  integer :: dlo(NDIM), dhi(NDIM)
523  real(dp) :: frac_gas, dr, fac_to_volume
524 
525  nb = surface%direction
526  id_in = surface%id_in
527  id_out = surface%id_out
528  nc = boxes(id_out)%n_cell
529  dr = boxes(id_out)%dr(af_neighb_dim(nb))
530 
531  ! Factor to convert surface charge density to volume density
532  fac_to_volume = 1 / dr
533 
534  ! Get index range for cells adjacent to the boundary of the dielectric box
535  call af_get_index_bc_inside(af_neighb_rev(nb), nc, 1, dlo, dhi)
536 
537  ! Get similar index range for gas-phase box
538  call af_get_index_bc_inside(nb, nc, 1, glo, ghi)
539 
540  ! How much of the rhs to put on the gas and dielectric side
541  frac_gas = 1.0_dp / (1.0_dp + surface%eps)
542 
543 #if NDIM == 2
544  boxes(id_out)%cc(glo(1):ghi(1), glo(2):ghi(2), i_rhs) = &
545  boxes(id_out)%cc(glo(1):ghi(1), glo(2):ghi(2), i_rhs) &
546  + frac_gas * fac * fac_to_volume * &
547  reshape(surface%sd(:, i_sigma), [ghi(1)-glo(1)+1, ghi(2)-glo(2)+1])
548 
549  boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), i_rhs) = &
550  boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), i_rhs) &
551  + (1-frac_gas) * fac * fac_to_volume * &
552  reshape(surface%sd(:, i_sigma), [dhi(1)-dlo(1)+1, dhi(2)-dlo(2)+1])
553 #elif NDIM == 3
554 !FLAG
555  boxes(id_out)%cc(glo(1):ghi(1), glo(2):ghi(2), glo(3):ghi(3), i_rhs) = &
556  boxes(id_out)%cc(glo(1):ghi(1), glo(2):ghi(2), glo(3):ghi(3), i_rhs) &
557  + frac_gas * fac * fac_to_volume * &
558  reshape(surface%sd(:, :, i_sigma), [ghi(1)-glo(1)+1, ghi(2)-glo(2)+1, ghi(3)-glo(3)+1])
559 
560  boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), i_rhs) = &
561  boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), i_rhs) &
562  + (1-frac_gas) * fac * fac_to_volume * &
563  reshape(surface%sd(:, :, i_sigma), [dhi(1)-dlo(1)+1, dhi(2)-dlo(2)+1, dhi(3)-dlo(3)+1])
564 #endif
565 
566  end subroutine surface_charge_to_rhs
567 
569  subroutine surface_inside_layer_to_surface(tree, diel, i_cc, i_sigma, &
570  fac, clear_cc, clear_surf)
571  type(af_t), intent(inout) :: tree
572  type(surfaces_t), intent(inout) :: diel
573  integer, intent(in) :: i_cc
574  integer, intent(in) :: i_sigma
575  real(dp), intent(in) :: fac
576  logical, intent(in) :: clear_cc
577  logical, intent(in) :: clear_surf
578  integer :: n
579 
580  if (.not. diel%initialized) error stop "dielectric not initialized"
581 
582  do n = 1, diel%max_ix
583  if (.not. diel%surfaces(n)%in_use) cycle
584 
585  call inside_layer_to_surface(tree%boxes, diel%surfaces(n), &
586  i_cc, i_sigma, fac, clear_cc, clear_surf)
587  end do
588  end subroutine surface_inside_layer_to_surface
589 
591  subroutine inside_layer_to_surface(boxes, surface, i_cc, i_sigma, fac, &
592  clear_cc, clear_surf)
593  type(box_t), intent(inout) :: boxes(:)
594  type(surface_t), intent(inout) :: surface
595  integer, intent(in) :: i_cc
596  integer, intent(in) :: i_sigma
597  real(dp), intent(in) :: fac
598  logical, intent(in) :: clear_cc
599  logical, intent(in) :: clear_surf
600  integer :: nb, nc, id_in
601  integer :: dlo(NDIM), dhi(NDIM)
602  real(dp) :: dr
603 
604  nb = surface%direction
605  id_in = surface%id_in
606  nc = boxes(id_in)%n_cell
607  dr = boxes(id_in)%dr(af_neighb_dim(nb))
608 
609  ! Get index range for cells adjacent to the boundary of the dielectric box
610  call af_get_index_bc_inside(af_neighb_rev(nb), nc, 1, dlo, dhi)
611 
612 #if NDIM == 2
613  if (clear_surf) surface%sd(:, i_sigma) = 0.0_dp
614 
615  surface%sd(:, i_sigma) = surface%sd(:, i_sigma) + fac * dr * &
616  pack(boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), i_cc), .true.)
617  if (clear_cc) boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), i_cc) = 0.0_dp
618 #elif NDIM == 3
619  if (clear_surf) surface%sd(:, :, i_sigma) = 0.0_dp
620 !FLAG
621  surface%sd(:,:, i_sigma) = surface%sd(:, :, i_sigma) + fac * dr * &
622  reshape(boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), i_cc), [nc, nc])
623  if (clear_cc) boxes(id_in)%cc(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), i_cc) = 0.0_dp
624 #endif
625 
626  end subroutine inside_layer_to_surface
627 
629  subroutine surface_correct_field_fc(tree, diel, i_sigma, i_fld, i_phi, fac)
630  type(af_t), intent(inout) :: tree
631  type(surfaces_t), intent(in) :: diel
632  integer, intent(in) :: i_sigma
633  integer, intent(in) :: i_fld
634  integer, intent(in) :: i_phi
635  real(dp), intent(in) :: fac
636  integer :: id_in, id_out, nc, nb, ix
637  real(dp) :: eps, fac_fld(2), fac_charge, inv_dr(ndim)
638  integer :: dlo(ndim), dhi(ndim)
639  integer :: glo(ndim), ghi(ndim)
640  integer :: dlo_shft(ndim), dhi_shft(ndim)
641  integer :: glo_shft(ndim), ghi_shft(ndim)
642  integer :: shft(ndim, 2*ndim)
643 
644  if (.not. diel%initialized) error stop "dielectric not initialized"
645 
646  nc = tree%n_cell
647 
648  do ix = 1, diel%max_ix
649  if (diel%surfaces(ix)%in_use) then
650  nb = diel%surfaces(ix)%direction
651  eps = diel%surfaces(ix)%eps
652  id_out = diel%surfaces(ix)%id_out
653  id_in = diel%surfaces(ix)%id_in
654  inv_dr = 1/tree%boxes(id_out)%dr
655 
656  fac_fld = [2 * eps, 2.0_dp] / (1 + eps)
657  fac_charge = fac / (1 + eps)
658 
659  ! Get index range for cells adjacent to the boundary of the dielectric box
660  call af_get_index_bface_inside(af_neighb_rev(nb), nc, 1, dlo, dhi)
661 
662  ! Get similar index range for gas-phase box
663  call af_get_index_bface_inside(nb, nc, 1, glo, ghi)
664 
665  ! Compute the offset for finite difference
666  shft = - abs(af_neighb_dix(:, :))
667 
668  dlo_shft(:) = dlo(:) + shft(:, nb)
669  dhi_shft(:) = dhi(:) + shft(:, nb)
670  glo_shft(:) = glo(:) + shft(:, nb)
671  ghi_shft(:) = ghi(:) + shft(:, nb)
672 
673  associate(fc_out => tree%boxes(id_out)%fc, &
674  fc_in => tree%boxes(id_in)%fc, &
675  cc_out => tree%boxes(id_out)%cc, &
676  cc_in => tree%boxes(id_in)%cc, &
677  sd => diel%surfaces(ix)%sd)
678  select case (nb)
679 #if NDIM == 2
680  case (af_neighb_lowx)
681  fc_out(1, 1:nc, 1, i_fld) = fac_fld(1) * inv_dr(1) * &
682  (cc_out(0, 1:nc, i_phi) - cc_out(1, 1:nc, i_phi)) &
683  + fac_charge * sd(:, i_sigma)
684  fc_in(nc+1, 1:nc, 1, i_fld) = fac_fld(2) * inv_dr(1) * &
685  (cc_in(nc, 1:nc, i_phi) - cc_in(nc+1, 1:nc, i_phi)) &
686  - fac_charge * sd(:, i_sigma)
687  case (af_neighb_highx)
688  fc_out(nc+1, 1:nc, 1, i_fld) = fac_fld(1) * inv_dr(1) * &
689  (cc_out(nc, 1:nc, i_phi) - cc_out(nc+1, 1:nc, i_phi)) &
690  - fac_charge * sd(:, i_sigma)
691  fc_in(1, 1:nc, 1, i_fld) = fac_fld(2) * inv_dr(1) * &
692  (cc_in(0, 1:nc, i_phi) - cc_in(1, 1:nc, i_phi)) &
693  + fac_charge * sd(:, i_sigma)
694  case (af_neighb_lowy)
695  fc_out(1:nc, 1, 2, i_fld) = fac_fld(1) * inv_dr(2) * &
696  (cc_out(1:nc, 0, i_phi) - cc_out(1:nc, 1, i_phi)) &
697  + fac_charge * sd(:, i_sigma)
698  fc_in(1:nc, nc+1, 2, i_fld) = fac_fld(2) * inv_dr(2) * &
699  (cc_in(1:nc, nc, i_phi) - cc_in(1:nc, nc+1, i_phi)) &
700  - fac_charge * sd(:, i_sigma)
701  case (af_neighb_highy)
702  fc_out(1:nc, nc+1, 2, i_fld) = fac_fld(1) * inv_dr(2) * &
703  (cc_out(1:nc, nc, i_phi) - cc_out(1:nc, nc+1, i_phi)) &
704  - fac_charge * sd(:, i_sigma)
705  fc_in(1:nc, 1, 2, i_fld) = fac_fld(2) * inv_dr(2) * &
706  (cc_in(1:nc, 0, i_phi) - cc_in(1:nc, 1, i_phi)) &
707  + fac_charge * sd(:, i_sigma)
708 #elif NDIM == 3
709  case default
710  fc_out(glo(1):ghi(1), glo(2):ghi(2), glo(3):ghi(3), af_neighb_dim(nb), i_fld) = &
711  fac_fld(1) * inv_dr(af_neighb_dim(nb)) * &
712  (cc_out(glo_shft(1):ghi_shft(1), glo_shft(2):ghi_shft(2), glo_shft(3):ghi_shft(3), i_phi) - &
713  cc_out(glo(1):ghi(1), glo(2):ghi(2), glo(3):ghi(3), i_phi)) &
714  + fac_charge * reshape(sd(:, :, i_sigma), [ghi(1)-glo(1)+1, ghi(2)-glo(2)+1, ghi(3)-glo(3)+1])
715 
716  fc_in(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), af_neighb_dim(nb), i_fld) = &
717  fac_fld(2) * inv_dr(af_neighb_dim(nb)) * &
718  (cc_in(dlo_shft(1):dhi_shft(1), dlo_shft(2):dhi_shft(2), dlo_shft(3):dhi_shft(3), i_phi) - &
719  cc_in(dlo(1):dhi(1), dlo(2):dhi(2), dlo(3):dhi(3), i_phi)) &
720  - fac_charge * reshape(sd(:, :, i_sigma), [dhi(1)-dlo(1)+1, dhi(2)-dlo(2)+1, dhi(3)-dlo(3)+1])
721 #endif
722  end select
723  end associate
724  end if
725  end do
726 
727  end subroutine surface_correct_field_fc
728 
729 #if NDIM == 2
730 
731  subroutine surface_correct_field_cc(tree, diel, i_sigma, i_fld, i_phi, fac)
732  type(af_t), intent(inout) :: tree
733  type(surfaces_t), intent(in) :: diel
734  integer, intent(in) :: i_sigma
735  integer, intent(in) :: i_fld(ndim)
736  integer, intent(in) :: i_phi
737  real(dp), intent(in) :: fac
738  integer :: id_in, id_out, nc, nb, ix
739  real(dp) :: eps, fac_fld(2), fac_charge, inv_dr(ndim)
740  real(dp) :: f_out(tree%n_cell), f_in(tree%n_cell)
741 
742  if (.not. diel%initialized) error stop "dielectric not initialized"
743 
744  nc = tree%n_cell
745 
746  do ix = 1, diel%max_ix
747  if (diel%surfaces(ix)%in_use) then
748  nb = diel%surfaces(ix)%direction
749  eps = diel%surfaces(ix)%eps
750  id_out = diel%surfaces(ix)%id_out
751  id_in = diel%surfaces(ix)%id_in
752  inv_dr = 1/tree%boxes(id_out)%dr
753 
754  fac_fld = [2 * eps, 2.0_dp] / (1 + eps)
755  fac_charge = fac / (1 + eps)
756 
757  associate(cc_out => tree%boxes(id_out)%cc, &
758  cc_in => tree%boxes(id_in)%cc, &
759  sd => diel%surfaces(ix)%sd)
760 
761  ! Compute field at two cell faces and average them
762  select case (nb)
763  case (af_neighb_lowx)
764  f_out = fac_fld(1) * inv_dr(1) * &
765  (cc_out(0, 1:nc, i_phi) - cc_out(1, 1:nc, i_phi)) &
766  + fac_charge * sd(:, i_sigma)
767  f_in = fac_fld(2) * inv_dr(1) * &
768  (cc_in(nc, 1:nc, i_phi) - cc_in(nc+1, 1:nc, i_phi)) &
769  - fac_charge * sd(:, i_sigma)
770  cc_out(1, 1:nc, i_fld(1)) = 0.5_dp * (f_out + inv_dr(1) * &
771  (cc_out(1, 1:nc, i_phi) - cc_out(2, 1:nc, i_phi)))
772  cc_out(0, 1:nc, i_fld(1)) = f_out
773  cc_in(nc, 1:nc, i_fld(1)) = 0.5_dp * (f_in + inv_dr(1) * &
774  (cc_in(nc-1, 1:nc, i_phi) - cc_in(nc, 1:nc, i_phi)))
775  cc_in(nc+1, 1:nc, i_fld(1)) = f_in
776  case (af_neighb_highx)
777  f_out = fac_fld(1) * inv_dr(1) * &
778  (cc_out(nc, 1:nc, i_phi) - cc_out(nc+1, 1:nc, i_phi)) &
779  - fac_charge * sd(:, i_sigma)
780  f_in = fac_fld(2) * inv_dr(1) * &
781  (cc_in(0, 1:nc, i_phi) - cc_in(1, 1:nc, i_phi)) &
782  + fac_charge * sd(:, i_sigma)
783  cc_out(nc, 1:nc, i_fld(1)) = 0.5_dp * (f_out + inv_dr(1) * &
784  (cc_out(nc-1, 1:nc, i_phi) - cc_out(nc, 1:nc, i_phi)))
785  cc_out(nc+1, 1:nc, i_fld(1)) = f_out
786  cc_in(1, 1:nc, i_fld(1)) = 0.5_dp * (f_in + inv_dr(1) * &
787  (cc_in(1, 1:nc, i_phi) - cc_in(2, 1:nc, i_phi)))
788  cc_in(0, 1:nc, i_fld(1)) = f_in
789  case (af_neighb_lowy)
790  f_out = fac_fld(1) * inv_dr(2) * &
791  (cc_out(1:nc, 0, i_phi) - cc_out(1:nc, 1, i_phi)) &
792  + fac_charge * sd(:, i_sigma)
793  f_in = fac_fld(2) * inv_dr(2) * &
794  (cc_in(1:nc, nc, i_phi) - cc_in(1:nc, nc+1, i_phi)) &
795  - fac_charge * sd(:, i_sigma)
796  cc_out(1:nc, 1, i_fld(2)) = 0.5_dp * (f_out + inv_dr(2) * &
797  (cc_out(1:nc, 1, i_phi) - cc_out(1:nc, 2, i_phi)))
798  cc_out(1:nc, 0, i_fld(2)) = f_out
799  cc_in(1:nc, nc, i_fld(2)) = 0.5_dp * (f_in + inv_dr(2) * &
800  (cc_in(1:nc, nc-1, i_phi) - cc_in(1:nc, nc, i_phi)))
801  cc_in(1:nc, nc+1, i_fld(2)) = f_in
802  case (af_neighb_highy)
803  f_out = fac_fld(1) * inv_dr(2) * &
804  (cc_out(1:nc, nc, i_phi) - cc_out(1:nc, nc+1, i_phi)) &
805  - fac_charge * sd(:, i_sigma)
806  f_in = fac_fld(2) * inv_dr(2) * &
807  (cc_in(1:nc, 0, i_phi) - cc_in(1:nc, 1, i_phi)) &
808  + fac_charge * sd(:, i_sigma)
809  cc_out(1:nc, nc, i_fld(2)) = 0.5_dp * (f_out + inv_dr(2) * &
810  (cc_out(1:nc, nc-1, i_phi) - cc_out(1:nc, nc, i_phi)))
811  cc_out(1:nc, nc+1, i_fld(2)) = f_out
812  cc_in(1:nc, 1, i_fld(2)) = 0.5_dp * (f_in + inv_dr(2) * &
813  (cc_in(1:nc, 1, i_phi) - cc_in(1:nc, 2, i_phi)))
814  cc_in(1:nc, 0, i_fld(2)) = f_in
815  end select
816  end associate
817  end if
818  end do
819  end subroutine surface_correct_field_cc
820 #endif
821 
822  subroutine surface_get_surface_cell(tree, diel, x, ix_surf, ix_cell)
823  use m_af_utils
824  type(af_t), intent(in) :: tree
825  type(surfaces_t), intent(in) :: diel
826  real(dp), intent(in) :: x(ndim)
827  integer, intent(out) :: ix_surf
828  integer, intent(out) :: ix_cell(ndim-1)
829  type(af_loc_t) :: loc
830  integer :: i, id, dim, direction
831 
832  ! Find location in box
833  loc = af_get_loc(tree, x)
834 
835  ! Check whether id is valid
836  id = loc%id
837  if (id == -1) error stop "Coordinate out of domain"
838 
839  ! Check for a surface from both sides
840  ix_surf = diel%box_id_out_to_surface_ix(id)
841  if (ix_surf == surface_none) ix_surf = diel%box_id_in_to_surface_ix(id)
842  if (ix_surf == surface_none) error stop "No surface found"
843 
844  ! Extract the cell index but only parallel to the surface
845  direction = diel%surfaces(ix_surf)%direction
846  dim = af_neighb_dim(direction)
847  ix_cell = pack(loc%ix, [(i, i=1,ndim)] /= dim)
848  end subroutine surface_get_surface_cell
849 
851  subroutine surface_write_output(tree, diel, i_vars, var_names, output_name, &
852  output_cnt)
853  use m_npy
854  type(af_t), intent(inout) :: tree
855  type(surfaces_t), intent(in) :: diel
856  integer, intent(in) :: i_vars(:)
857  character(len=*), intent(in) :: var_names(:)
858  character(len=*), intent(in) :: output_name
859  integer, intent(in) :: output_cnt
860  integer :: n, i, ix, nc, n_vars
861  integer :: lo(ndim-1), hi(ndim-1)
862  character(len=200) :: tmpname, filename
863 
864  real(dp) :: coords(ndim, tree%n_cell**(ndim-1))
865  real(dp), allocatable :: r(:, :), dr(:, :)
866  real(dp), allocatable :: surface_vars(:, :)
867  integer, allocatable :: surf_dim(:)
868 
869  nc = diel%n_cell
870  n = count(diel%surfaces(1:diel%max_ix)%in_use)
871  n_vars = size(i_vars)
872  allocate(r(ndim, n*nc), surface_vars(n*nc, n_vars))
873  allocate(dr(ndim-1, n), surf_dim(n))
874 
875  write(filename, "(A,I6.6,A)") trim(output_name) // "_", &
876  output_cnt, "_surface.npz"
877  write(tmpname, "(A,I6.6,A)") trim(output_name) // "_", &
878  output_cnt, "_tmp.npy"
879 
880  i = 0
881  do ix = 1, diel%max_ix
882  if (diel%surfaces(ix)%in_use) then
883  i = i + 1
884  lo = (i-1) * nc + 1
885  hi = i * nc
886 
887  associate(box => tree%boxes(diel%surfaces(ix)%id_out), &
888  surf => diel%surfaces(ix))
889  call af_get_face_coords(box, surf%direction, coords)
890 #if NDIM == 2
891  r(:, lo(1):hi(1)) = coords
892  dr(:, i) = surf%dr
893  surf_dim(i) = af_neighb_dim(surf%direction)
894  surface_vars(lo(1):hi(1), :) = surf%sd(:, i_vars)
895 #elif NDIM == 3
896  error stop "not implemented"
897 #endif
898  end associate
899  end if
900  end do
901 
902  call save_npy(tmpname, [n])
903  call add_to_zip(filename, tmpname, .false., "n_surfaces")
904  call save_npy(tmpname, [nc])
905  call add_to_zip(filename, tmpname, .false., "n_cell")
906  call save_npy(tmpname, r)
907  call add_to_zip(filename, tmpname, .false., "r")
908  call save_npy(tmpname, dr)
909  call add_to_zip(filename, tmpname, .false., "dr")
910  call save_npy(tmpname, surf_dim)
911  call add_to_zip(filename, tmpname, .false., "normal_dim")
912 
913  do n = 1, n_vars
914  call save_npy(tmpname, surface_vars(:, n))
915  call add_to_zip(filename, tmpname, .false., trim(var_names(n)))
916  end do
917  print *, "surface_write_output: written " // trim(filename)
918 
919  end subroutine surface_write_output
920 
921 end module m_af_surface
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
This module contains routines for including flat surfaces between changes in epsilon (some material p...
Definition: m_af_surface.f90:5
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Definition: m_af_types.f90:4
This module contains all kinds of different 'helper' routines for Afivo. If the number of routines fo...
Definition: m_af_utils.f90:5
Definition: m_npy.f90:1
Type for a single surface.
Type for storing all the surfaces on a mesh.
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
Definition: m_af_types.f90:326
Type that contains the refinement changes in a tree.
Definition: m_af_types.f90:88