Afivo  0.3
m_a3_core.f90
1 
4 module m_a3_core
5  use m_a3_types
6 
7  implicit none
8  private
9 
10  public :: a3_init
11  public :: a3_set_cc_methods
12  public :: a3_init_box
13  public :: a3_destroy
14  public :: a3_set_base
15  public :: a3_tidy_up
16  public :: a3_resize_box_storage
17  public :: a3_adjust_refinement
18  public :: a3_consistent_fluxes
19 
20 contains
21 
23  subroutine a3_init(tree, n_cell, n_var_cell, n_var_face, dr, r_min, &
24  lvl_limit, n_boxes, coarsen_to, coord, cc_names, fc_names, mem_limit_gb)
25  type(a3_t), intent(inout) :: tree
26  integer, intent(in) :: n_cell
27  integer, intent(in) :: n_var_cell
28  integer, intent(in) :: n_var_face
29  real(dp), intent(in) :: dr
31  real(dp), intent(in), optional :: r_min(3)
34  integer, intent(in), optional :: coarsen_to
36  integer, intent(in), optional :: lvl_limit
38  integer, intent(in), optional :: n_boxes
39  integer, intent(in), optional :: coord
40  real(dp), intent(in), optional :: mem_limit_gb
41 
43  character(len=*), intent(in), optional :: cc_names(:)
45  character(len=*), intent(in), optional :: fc_names(:)
46 
47  integer :: lvl_limit_a, n_boxes_a, coarsen_to_a
48  real(dp) :: r_min_a(3), gb_limit
49  integer :: n, lvl, min_lvl, coord_a, box_bytes
50 
51  ! Set default arguments if not present
52  lvl_limit_a = 30; if (present(lvl_limit)) lvl_limit_a = lvl_limit
53  n_boxes_a = 1000; if (present(n_boxes)) n_boxes_a = n_boxes
54  coarsen_to_a = -1; if (present(coarsen_to)) coarsen_to_a = coarsen_to
55  r_min_a = 0.0_dp; if (present(r_min)) r_min_a = r_min
56  coord_a = af_xyz; if (present(coord)) coord_a = coord
57  gb_limit = 16; if (present(mem_limit_gb)) gb_limit = mem_limit_gb
58 
59  if (tree%ready) stop "a3_init: tree was already initialized"
60  if (n_cell < 2) stop "a3_init: n_cell should be >= 2"
61  if (btest(n_cell, 0)) stop "a3_init: n_cell should be even"
62  if (n_var_cell <= 0) stop "a3_init: n_var_cell should be > 0"
63  if (n_var_face < 0) stop "a3_init: n_var_face should be >= 0"
64  if (n_boxes_a <= 0) stop "a3_init: n_boxes should be > 0"
65  if (lvl_limit_a <= 0) stop "a3_init: lvl_limit should be > 0"
66  if (gb_limit <= 0) stop "a3_init: mem_limit_gb should be > 0"
67  if (coord_a == af_cyl) stop "a3_init: cannot have 3d cyl coords"
68 
69  allocate(tree%boxes(n_boxes_a))
70 
71  if (coarsen_to_a > 0) then
72  ! Determine number of lvls for subtree
74  min_lvl = 1 - nint(log(real(n_cell, dp)/coarsen_to_a)/log(2.0_dp))
75 
76  if (2**(1-min_lvl) * coarsen_to_a /= n_cell) &
77  stop "a3_init: cannot coarsen to given value"
78  else
79  min_lvl = 1
80  end if
81 
82  ! up to lvl_limit_a+1 to add dummies that are always of size zero
83  allocate(tree%lvls(min_lvl:lvl_limit_a+1))
84 
85  do lvl = min_lvl, lvl_limit_a+1
86  allocate(tree%lvls(lvl)%ids(0))
87  allocate(tree%lvls(lvl)%leaves(0))
88  allocate(tree%lvls(lvl)%parents(0))
89  end do
90 
91  tree%n_cell = n_cell
92  tree%n_var_cell = n_var_cell
93  tree%n_var_face = n_var_face
94  tree%r_base = r_min_a
95  tree%dr_base = dr
96  tree%lvl_limit = lvl_limit_a
97  tree%highest_id = 0
98  tree%highest_lvl = 0
99  tree%coord_t = coord_a
100 
101  ! Calculate size of a box
102  box_bytes = a3_box_bytes(n_cell, n_var_cell, n_var_face)
103  tree%box_limit = nint(gb_limit * 2.0_dp**30 / box_bytes)
104 
105  ! Set variable names
106  allocate(tree%cc_names(n_var_cell))
107  allocate(tree%fc_names(n_var_face))
108 
109  if (present(cc_names)) then
110  if (size(cc_names) /= n_var_cell) &
111  stop "a3_init: size(cc_names) /= n_var_cell"
112  tree%cc_names = cc_names
113  else
114  do n = 1, n_var_cell
115  write(tree%cc_names(n), "(A,I0)") "cc_var", n
116  end do
117  end if
118 
119  if (present(fc_names)) then
120  if (size(fc_names) /= n_var_face) &
121  stop "a3_init: size(fc_names) /= n_var_face"
122  tree%fc_names = fc_names
123  else
124  do n = 1, n_var_face
125  write(tree%fc_names(n), "(A,I0)") "flux", n
126  end do
127  end if
128 
129  ! Initialize cell-centered methods (default is the null pointer)
130  allocate(tree%cc_methods(n_var_cell))
131  allocate(tree%has_cc_method(n_var_cell))
132  tree%has_cc_method(:) = .false.
133  allocate(tree%cc_method_vars(0))
134  end subroutine a3_init
135 
137  subroutine a3_set_cc_methods(tree, iv, bc, rb, prolong, restrict)
138  use m_a3_ghostcell, only: a3_gc_interp
139  use m_a3_prolong, only: a3_prolong_linear
140  use m_a3_restrict, only: a3_restrict_box
141  type(a3_t), intent(inout) :: tree
142  integer, intent(in) :: iv
143  procedure(a3_subr_bc) :: bc
144  procedure(a3_subr_rb), optional :: rb
145  procedure(a3_subr_prolong), optional :: prolong
146  procedure(a3_subr_restrict), optional :: restrict
147  integer :: i, n
148  tree%cc_methods(iv)%bc => bc
149 
150  if (present(rb)) then
151  tree%cc_methods(iv)%rb => rb
152  else
153  tree%cc_methods(iv)%rb => a3_gc_interp
154  end if
155 
156  if (present(prolong)) then
157  tree%cc_methods(iv)%prolong => prolong
158  else
159  tree%cc_methods(iv)%prolong => a3_prolong_linear
160  end if
161 
162  if (present(restrict)) then
163  tree%cc_methods(iv)%restrict => restrict
164  else
165  tree%cc_methods(iv)%restrict => a3_restrict_box
166  end if
167 
168  tree%has_cc_method(iv) = .true.
169  n = size(tree%cc_method_vars)
170  tree%cc_method_vars = [(tree%cc_method_vars(i), i=1,n), iv]
171  end subroutine a3_set_cc_methods
172 
175  subroutine a3_destroy(tree)
176  type(a3_t), intent(inout) :: tree
177 
178  if (.not. tree%ready) stop "a3_destroy: Tree not fully initialized"
179  deallocate(tree%boxes)
180  deallocate(tree%lvls)
181  deallocate(tree%cc_names)
182  deallocate(tree%fc_names)
183  deallocate(tree%cc_methods)
184  deallocate(tree%has_cc_method)
185  deallocate(tree%cc_method_vars)
186 
187  tree%highest_id = 0
188  tree%ready = .false.
189  end subroutine a3_destroy
190 
197  subroutine a3_set_base(tree, n_boxes, ix_list, nb_list)
199  type(a3_t), intent(inout) :: tree
201  integer, intent(in) :: n_boxes
203  integer, intent(in) :: ix_list(3, n_boxes)
205  integer, intent(inout), optional :: nb_list(a3_num_neighbors, n_boxes)
206  integer :: n, id, nb, nb_id
207  integer :: ix(3), lvl, offset
208  integer :: ix_min(3), ix_max(3), nb_ix(3)
209  integer :: nb_used(a3_num_neighbors, size(ix_list, 2))
210  integer, allocatable :: id_array(:, :, :)
211 
212  if (n_boxes < 1) stop "a3_set_base: need at least one box"
213  if (any(ix_list < 1)) stop "a3_set_base: need all ix_list > 0"
214  if (tree%highest_id > 0) stop "a3_set_base: this tree already has boxes"
215  if (.not. allocated(tree%lvls)) stop "a3_set_base: tree not initialized"
216 
217  if (present(nb_list)) then
218  nb_used = nb_list
219  else
220  nb_used = af_no_box ! Default value
221  end if
222 
223  ! Create an array covering the coarse grid, in which boxes are indicated by
224  ! their id.
225  do n = 1, 3
226  ! -1 and +1 to also include 'neighbors' of the coarse grid
227  ix_min(n) = minval(ix_list(n, :)) - 1
228  ix_max(n) = maxval(ix_list(n, :)) + 1
229  end do
230 
231  allocate(id_array(ix_min(1):ix_max(1), &
232  ix_min(2):ix_max(2), ix_min(3):ix_max(3)))
233 
234  ! Store box ids in the index array covering the coarse grid
235  id_array = af_no_box
236 
237  do id = 1, n_boxes
238  ix = ix_list(:, id)
239  id_array(ix(1), ix(2), ix(3)) = id
240  end do
241 
242  ! Loop over the boxes and set their neighbors
243  do id = 1, n_boxes
244  ix = ix_list(:, id)
245 
246  do nb = 1, a3_num_neighbors
247  ! Compute ix of neighbor
248  nb_ix = ix + a3_neighb_dix(:, nb)
249  nb_id = id_array(nb_ix(1), nb_ix(2), nb_ix(3))
250 
251  if (nb_id /= af_no_box) then
252  ! Neighbor present, so store id
253  nb_used(nb, id) = nb_id
254  else
255  ! A periodic or boundary condition
256  nb_id = nb_used(nb, id)
257 
258  if (nb_id > af_no_box) then
259  ! If periodic, copy connectivity information to other box
260  nb_used(a3_neighb_rev(nb), nb_id) = id
261  else if (nb_id == af_no_box) then
262  ! The value af_no_box is converted to af_phys_boundary,
263  ! indicating the default boundary condition
264  nb_used(nb, id) = af_phys_boundary
265  end if
266  end if
267 
268  end do
269  end do
270 
271  if (any(nb_used == af_no_box)) stop "a3_set_base: unresolved neighbors"
272 
273  ! Check if we have enough space, if not, increase space
274  if (n_boxes > size(tree%boxes(:))) then
275  call a3_resize_box_storage(tree, n_boxes)
276  end if
277 
278  ! Create coarser levels which are copies of lvl 1
279  do lvl = lbound(tree%lvls, 1), 1
280  deallocate(tree%lvls(lvl)%ids)
281  allocate(tree%lvls(lvl)%ids(n_boxes))
282 
283  call get_free_ids(tree, tree%lvls(lvl)%ids)
284  offset = tree%lvls(lvl)%ids(1) - 1
285 
286  do n = 1, n_boxes
287  id = tree%lvls(lvl)%ids(n)
288  ix = ix_list(:, n)
289  tree%boxes(id)%lvl = lvl
290  tree%boxes(id)%ix = ix
291  tree%boxes(id)%dr = tree%dr_base * 0.5_dp**(lvl-1)
292  tree%boxes(id)%r_min = tree%r_base + &
293  (ix - 1) * tree%dr_base * tree%n_cell
294  tree%boxes(id)%n_cell = tree%n_cell / (2**(1-lvl))
295  tree%boxes(id)%coord_t = tree%coord_t
296 
297  tree%boxes(id)%parent = af_no_box
298  tree%boxes(id)%children(:) = af_no_box ! Gets overwritten, see below
299 
300  ! Connectivity is the same for all lvls
301  where (nb_used(:, n) > af_no_box)
302  tree%boxes(id)%neighbors = nb_used(:, n) + offset
303  elsewhere
304  tree%boxes(id)%neighbors = nb_used(:, n)
305  end where
306 
307  tree%boxes(id)%neighbor_mat = af_no_box
308 
309  call a3_init_box(tree%boxes(id), tree%boxes(id)%n_cell, &
310  tree%n_var_cell, tree%n_var_face)
311  end do
312 
313  if (lvl == 1) then
314  deallocate(tree%lvls(lvl)%leaves)
315  allocate(tree%lvls(lvl)%leaves(n_boxes))
316  tree%lvls(lvl)%leaves = tree%lvls(lvl)%ids
317  else
318  deallocate(tree%lvls(lvl)%parents)
319  allocate(tree%lvls(lvl)%parents(n_boxes))
320  tree%lvls(lvl)%parents = tree%lvls(lvl)%ids
321  end if
322 
323  if (lvl > lbound(tree%lvls, 1)) then
324  tree%boxes(tree%lvls(lvl-1)%ids)%children(1) = &
325  tree%lvls(lvl)%ids
326  tree%boxes(tree%lvls(lvl)%ids)%parent = &
327  tree%lvls(lvl-1)%ids
328  end if
329  end do
330 
331  ! Set the diagonal neighbors
332  do lvl = lbound(tree%lvls, 1), 1
333  do n = 1, n_boxes
334  id = tree%lvls(lvl)%ids(n)
335  call set_neighbor_mat(tree%boxes, id)
336  end do
337  end do
338 
339  tree%highest_lvl = 1
340  tree%ready = .true.
341 
342  end subroutine a3_set_base
343 
344  subroutine set_neighbor_mat(boxes, id)
345  type(box3_t), intent(inout) :: boxes(:)
346  integer, intent(in) :: id
347  integer :: i, j, k, directions(3), ndir, nb_id
348 
349  do k = -1, 1; do j = -1, 1; do i = -1, 1
350  if (boxes(id)%neighbor_mat(i, j, k) /= af_no_box) cycle
351 
352  ndir = 0
353 
354  if (i /= 0) then
355  ndir = ndir + 1
356  directions(ndir) = a3_neighb_lowx + (i + 1)/2
357  end if
358 
359  if (j /= 0) then
360  ndir = ndir + 1
361  directions(ndir) = a3_neighb_lowy + (j + 1)/2
362  end if
363 
364  if (k /= 0) then
365  ndir = ndir + 1
366  directions(ndir) = a3_neighb_lowz + (k + 1)/2
367  end if
368 
369  boxes(id)%neighbor_mat(i, j, k) = &
370  find_neighb_id(boxes, id, directions(1:ndir))
371 
372  nb_id = boxes(id)%neighbor_mat(i, j, k)
373 
374  if (nb_id > af_no_box) then
375  boxes(nb_id)%neighbor_mat(-i, -j, -k) = id
376  end if
377  end do; end do; end do
378 
379  end subroutine set_neighbor_mat
380 
383  pure function find_neighb_id(boxes, id, nbs) result(nb_id)
384  type(box3_t), intent(in) :: boxes(:)
385  integer, intent(in) :: id
386  integer, intent(in) :: nbs(:) ! List of neighbor directions
387  integer :: i, j, k, nb, nb_id
388  integer :: nbs_perm(size(nbs))
389 
390  nb_id = id
391 
392  if (size(nbs) == 0) then
393  return
394  else
395  do i = 1, size(nbs)
396  nb_id = id
397 
398  ! Check if path exists starting from nbs(i)
399  do j = 1, size(nbs)
400  ! k starts at i and runs over the neighbors
401  k = 1 + mod(i + j - 2, size(nbs))
402  nb = nbs(k)
403 
404  nb_id = boxes(nb_id)%neighbors(nb)
405  if (nb_id <= af_no_box) exit
406  end do
407 
408  if (nb_id > af_no_box) exit ! Found it
409  end do
410  end if
411 
412  ! For a corner neighbor in 3D, try again using the permuted neighbor list to
413  ! covers all paths
414  if (size(nbs) == 3 .and. nb_id <= af_no_box) then
415  nbs_perm = nbs([2,1,3])
416 
417  do i = 1, size(nbs)
418  nb_id = id
419 
420  do j = 1, size(nbs)
421  k = 1 + mod(i + j - 2, size(nbs))
422  nb = nbs(k)
423 
424  nb_id = boxes(nb_id)%neighbors(nb)
425  if (nb_id <= af_no_box) exit
426  end do
427 
428  if (nb_id > af_no_box) exit ! Found it
429  end do
430  end if
431 
432  end function find_neighb_id
433 
435  subroutine a3_tidy_up(tree, max_hole_frac, n_clean_min)
437  type(a3_t), intent(inout) :: tree
439  real(dp), intent(in) :: max_hole_frac
441  integer, intent(in) :: n_clean_min
442  real(dp) :: hole_frac
443  integer :: n, lvl, id, n_clean, i, j, k
444  integer :: highest_id, n_used, n_stored, n_used_lvl
445  integer, allocatable :: ixs_sort(:), ixs_map(:)
446  type(box3_t), allocatable :: boxes_cpy(:)
447  integer(morton_k), allocatable :: mortons(:)
448 
449  if (.not. tree%ready) stop "Tree not ready"
450  if (max_hole_frac < 0) stop "a3_tidy_up: need max_hole_frac >= 0"
451  if (n_clean_min < 0) stop "a3_tidy_up: need n_clean_min >= 0"
452 
453  highest_id = tree%highest_id
454  n_used = a3_num_boxes_used(tree)
455  hole_frac = 1 - n_used / real(highest_id, dp)
456  n_clean = nint(hole_frac * highest_id)
457 
458  if (hole_frac > max_hole_frac .and. n_clean >= n_clean_min) then
459  allocate(boxes_cpy(n_used))
460  allocate(ixs_map(0:highest_id))
461  ixs_map(0) = 0
462  n_stored = 0
463 
464  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
465  n_used_lvl = size(tree%lvls(lvl)%ids)
466  allocate(mortons(n_used_lvl))
467  allocate(ixs_sort(n_used_lvl))
468 
469  do n = 1, n_used_lvl
470  id = tree%lvls(lvl)%ids(n)
471  ! Note the -1, since our indices start at 1
472  mortons(n) = morton_from_ix3(tree%boxes(id)%ix-1)
473  end do
474 
475  call morton_rank(mortons, ixs_sort)
476  tree%lvls(lvl)%ids = tree%lvls(lvl)%ids(ixs_sort)
477 
478  do n = 1, n_used_lvl
479  id = tree%lvls(lvl)%ids(n)
480  boxes_cpy(n_stored + n) = tree%boxes(id)
481  ixs_map(tree%lvls(lvl)%ids(n)) = n_stored + n
482  end do
483 
484  tree%lvls(lvl)%ids = [(n_stored+n, n=1,n_used_lvl)]
485  call set_leaves_parents(boxes_cpy, tree%lvls(lvl))
486  n_stored = n_stored + n_used_lvl
487  deallocate(mortons)
488  deallocate(ixs_sort)
489  end do
490 
491  ! Update id's to new indices
492  do n = 1, n_used
493  boxes_cpy(n)%parent = ixs_map(boxes_cpy(n)%parent)
494  boxes_cpy(n)%children = ixs_map(boxes_cpy(n)%children)
495  where (boxes_cpy(n)%neighbors > af_no_box)
496  boxes_cpy(n)%neighbors = ixs_map(boxes_cpy(n)%neighbors)
497  end where
498 
499  do k = -1, 1; do j = -1, 1; do i = -1, 1
500  if (boxes_cpy(n)%neighbor_mat(i, j, k) > af_no_box) then
501  boxes_cpy(n)%neighbor_mat(i, j, k) = &
502  ixs_map(boxes_cpy(n)%neighbor_mat(i, j, k))
503  end if
504  end do; end do; end do
505  end do
506 
507  tree%boxes(1:n_used) = boxes_cpy ! Copy ordered data
508  do n = n_used+1, highest_id
509  if (tree%boxes(n)%in_use) then
510  ! Remove moved data
511  call clear_box(tree%boxes(n))
512  end if
513  end do
514 
515  tree%highest_id = n_used
516  end if
517 
518  end subroutine a3_tidy_up
519 
521  subroutine set_leaves_parents(boxes, level)
522  type(box3_t), intent(in) :: boxes(:)
523  type(lvl_t), intent(inout) :: level
524  integer :: i, id, i_leaf, i_parent
525  integer :: n_parents, n_leaves
526 
527  n_parents = count(a3_has_children(boxes(level%ids)))
528  n_leaves = size(level%ids) - n_parents
529 
530  if (n_parents /= size(level%parents)) then
531  deallocate(level%parents)
532  allocate(level%parents(n_parents))
533  end if
534 
535  if (n_leaves /= size(level%leaves)) then
536  deallocate(level%leaves)
537  allocate(level%leaves(n_leaves))
538  end if
539 
540  i_leaf = 0
541  i_parent = 0
542  do i = 1, size(level%ids)
543  id = level%ids(i)
544  if (a3_has_children(boxes(id))) then
545  i_parent = i_parent + 1
546  level%parents(i_parent) = id
547  else
548  i_leaf = i_leaf + 1
549  level%leaves(i_leaf) = id
550  end if
551  end do
552  end subroutine set_leaves_parents
553 
556  subroutine a3_init_box(box, n_cell, n_cc, n_fc)
557  type(box3_t), intent(inout) :: box
558  integer, intent(in) :: n_cell
559  integer, intent(in) :: n_cc
560  integer, intent(in) :: n_fc
561 
562  box%in_use = .true.
563 
564  allocate(box%cc(0:n_cell+1, 0:n_cell+1, 0:n_cell+1, n_cc))
565  allocate(box%fc(n_cell+1, n_cell+1, n_cell+1, 3, n_fc))
566 
567  ! Initialize to zero
568  box%cc = 0
569  box%fc = 0
570  end subroutine a3_init_box
571 
573  subroutine clear_box(box)
574  type(box3_t), intent(inout) :: box
575 
576  box%in_use = .false.
577 
578  deallocate(box%cc)
579  deallocate(box%fc)
580  end subroutine clear_box
581 
582  ! Set the neighbors of id (using their parent)
583  subroutine set_neighbs_3d(boxes, id)
584  type(box3_t), intent(inout) :: boxes(:)
585  integer, intent(in) :: id
586  integer :: nb, nb_id, i, j, k
587 
588  do k = -1, 1; do j = -1, 1; do i = -1, 1
589  if (boxes(id)%neighbor_mat(i, j, k) == af_no_box) then
590  nb_id = find_neighb_3d(boxes, id, [i, j, k])
591  if (nb_id > af_no_box) then
592  boxes(id)%neighbor_mat(i, j, k) = nb_id
593  boxes(nb_id)%neighbor_mat(-i, -j, -k) = id
594  end if
595  end if
596  end do; end do; end do
597 
598  do nb = 1, a3_num_neighbors
599  if (boxes(id)%neighbors(nb) == af_no_box) then
600  nb_id = boxes(id)%neighbor_mat(a3_neighb_dix(1, nb), &
601  a3_neighb_dix(2, nb), a3_neighb_dix(3, nb))
602  if (nb_id > af_no_box) then
603  boxes(id)%neighbors(nb) = nb_id
604  boxes(nb_id)%neighbors(a3_neighb_rev(nb)) = id
605  end if
606  end if
607  end do
608  end subroutine set_neighbs_3d
609 
611  function find_neighb_3d(boxes, id, dix) result(nb_id)
612  type(box3_t), intent(in) :: boxes(:)
613  integer, intent(in) :: id
614  integer, intent(in) :: dix(3)
615  integer :: nb_id, p_id, c_ix, dix_c(3)
616 
617  p_id = boxes(id)%parent
618  c_ix = a3_ix_to_ichild(boxes(id)%ix)
619 
620  ! Check if neighbor is in same direction as dix is (low/high). If so, use
621  ! neighbor of parent
622  where ((dix == -1) .eqv. a3_child_low(:, c_ix))
623  dix_c = dix
624  elsewhere
625  dix_c = 0
626  end where
627 
628  p_id = boxes(p_id)%neighbor_mat(dix_c(1), dix_c(2), dix_c(3))
629 
630  if (p_id <= af_no_box) then
631  nb_id = p_id
632  else
633  c_ix = a3_ix_to_ichild(boxes(id)%ix + dix)
634  nb_id = boxes(p_id)%children(c_ix)
635  end if
636  end function find_neighb_3d
637 
639  subroutine a3_resize_box_storage(tree, new_size)
640  type(a3_t), intent(inout) :: tree
641  integer, intent(in) :: new_size
642  type(box3_t), allocatable :: boxes_cpy(:)
643 
644  if (.not. tree%ready) stop "a3_resize_box_storage: Tree not ready"
645  if (new_size < tree%highest_id) &
646  stop "a3_resize_box_storage: Cannot shrink tree"
647 
648  ! Store boxes in larger array boxes_cpy
649  allocate(boxes_cpy(new_size))
650  boxes_cpy(1:tree%highest_id) = tree%boxes(1:tree%highest_id)
651 
652  ! Deallocate current storage
653  deallocate(tree%boxes)
654 
655  ! Use new array
656  call move_alloc(boxes_cpy, tree%boxes)
657  end subroutine a3_resize_box_storage
658 
665  subroutine a3_adjust_refinement(tree, ref_subr, ref_info, ref_buffer, keep_buffer)
666  type(a3_t), intent(inout) :: tree
667  procedure(a3_subr_ref) :: ref_subr
668  type(ref_info_t), intent(inout) :: ref_info
669  integer, intent(in), optional :: ref_buffer
670  logical, intent(in), optional :: keep_buffer
671  integer :: lvl, id, i, c_ids(a3_num_children)
672  integer :: highest_id_prev, highest_id_req
673  integer, allocatable :: ref_flags(:)
674  integer :: num_new_boxes, total_num_boxes
675  integer :: new_size, ref_buffer_val
676  logical :: keep_buffer_val
677 
678  if (.not. tree%ready) stop "Tree not ready"
679 
680  ref_buffer_val = 0 ! Default buffer width (in cells) around refinement
681  if (present(ref_buffer)) ref_buffer_val = ref_buffer
682  keep_buffer_val = .false. ! Do not use a buffer for 'keep refinement' flags
683  if (present(keep_buffer)) keep_buffer_val = keep_buffer
684 
685  if (ref_buffer_val < 0) &
686  error stop "a3_adjust_refinement: ref_buffer < 0"
687  if (ref_buffer_val > tree%n_cell) &
688  error stop "a3_adjust_refinement: ref_buffer > tree%n_cell"
689 
690  ! Check whether the boxes array contains many holes, and if this is the
691  ! case, reorganize the array
692  call a3_tidy_up(tree, 0.5_dp, 1000)
693 
694  highest_id_prev = tree%highest_id
695  allocate(ref_flags(highest_id_prev))
696 
697  ! Set refinement values for all boxes
698  call consistent_ref_flags(tree, ref_flags, ref_subr, &
699  ref_buffer_val, keep_buffer_val)
700 
701  ! Compute number of new boxes
702  num_new_boxes = a3_num_children * count(ref_flags == af_refine)
703 
704  ! Determine number of boxes that could be in use after resizing
705  total_num_boxes = a3_num_boxes_used(tree) + num_new_boxes
706 
707  if (total_num_boxes > tree%box_limit) then
708  print *, "a3_adjust_refinement: exceeding memory limit"
709  write(*, '(A,E12.2)') " memory_limit (GByte): ", &
710  tree%box_limit * 0.5_dp**30 * &
711  a3_box_bytes(tree%n_cell, tree%n_var_cell, tree%n_var_face)
712  print *, "memory_limit (boxes): ", tree%box_limit
713  print *, "required number of boxes: ", total_num_boxes
714  print *, "You can increase the memory limit in your call to a3_init"
715  print *, "by setting mem_limit_gb to a higher value (in GBytes)"
716  error stop
717  end if
718 
719  ! Resize box storage if required
720  highest_id_req = highest_id_prev + num_new_boxes
721  if (highest_id_req > size(tree%boxes)) then
722  new_size = 2 * highest_id_req ! Allocate some extra empty boxes
723  call a3_resize_box_storage(tree, new_size)
724  end if
725 
726  do lvl = 1, tree%lvl_limit-1
727  do i = 1, size(tree%lvls(lvl)%ids)
728  id = tree%lvls(lvl)%ids(i)
729 
730  if (id > highest_id_prev) then
731  cycle ! This is a newly added box
732  else if (ref_flags(id) == af_refine) then
733  ! Add children. First need to get num_children free id's
734  call get_free_ids(tree, c_ids)
735  call add_children(tree%boxes, id, c_ids, &
736  tree%n_var_cell, tree%n_var_face)
737  else if (ref_flags(id) == af_derefine) then
738  ! Remove children
739  call auto_restrict(tree, id)
740  call remove_children(tree%boxes, id)
741  end if
742  end do
743 
744  ! Update leaves / parents
745  call set_leaves_parents(tree%boxes, tree%lvls(lvl))
746 
747  ! Set next level ids to children of this level
748  call set_child_ids(tree%lvls(lvl)%parents, &
749  tree%lvls(lvl+1)%ids, tree%boxes)
750 
751  ! Update connectivity of new children (id > highest_id_prev)
752  do i = 1, size(tree%lvls(lvl+1)%ids)
753  id = tree%lvls(lvl+1)%ids(i)
754  if (id > highest_id_prev) call set_neighbs_3d(tree%boxes, id)
755  end do
756 
757  if (size(tree%lvls(lvl+1)%ids) == 0) exit
758  end do
759 
760  tree%highest_lvl = lvl
761 
762  ! We still have to update leaves and parents for the last level, which is
763  ! either lvl+1 or lvl_limit. Note that lvl+1 is empty now, but maybe it was
764  ! not not empty before, and that lvl_limit is skipped in the above loop.
765  lvl = min(lvl+1, tree%lvl_limit)
766  call set_leaves_parents(tree%boxes, tree%lvls(lvl))
767 
768  ! Set information about the refinement
769  call set_ref_info(tree, ref_flags, ref_info)
770 
771  call auto_prolong(tree, ref_info)
772 
773  end subroutine a3_adjust_refinement
774 
776  subroutine auto_restrict(tree, id)
777  type(a3_t), intent(inout) :: tree
778  integer, intent(in) :: id
779  integer :: iv, i_ch, ch_id
780 
781  if (.not. any(tree%has_cc_method(:))) return
782 
783  do iv = 1, tree%n_var_cell
784  if (tree%has_cc_method(iv)) then
785  do i_ch = 1, a3_num_children
786  ch_id = tree%boxes(id)%children(i_ch)
787  call tree%cc_methods(iv)%restrict(tree%boxes(ch_id), &
788  tree%boxes(id), iv)
789  end do
790  end if
791  end do
792  end subroutine auto_restrict
793 
795  subroutine auto_prolong(tree, ref_info)
796  use m_a3_ghostcell, only: a3_gc_box
797  type(a3_t), intent(inout) :: tree
798  type(ref_info_t), intent(in) :: ref_info
799  integer :: lvl, i, n, iv, id, p_id
800 
801  ! Skip this routine when it won't do anything
802  if (.not. any(tree%has_cc_method(:)) .or. ref_info%n_add == 0) then
803  return
804  end if
805 
806  !$omp parallel private(lvl, i, n, iv, id, p_id)
807  do lvl = 1, tree%highest_lvl
808  !$omp do
809  do i = 1, size(ref_info%lvls(lvl)%add)
810  id = ref_info%lvls(lvl)%add(i)
811  p_id = tree%boxes(id)%parent
812 
813  do n = 1, size(tree%cc_method_vars)
814  iv = tree%cc_method_vars(n)
815  call tree%cc_methods(iv)%prolong(tree%boxes(p_id), &
816  tree%boxes(id), iv)
817  end do
818  end do
819  !$omp end do
820 
821  !$omp do
822  do i = 1, size(ref_info%lvls(lvl)%add)
823  id = ref_info%lvls(lvl)%add(i)
824  call a3_gc_box(tree, id, [tree%cc_method_vars])
825  end do
826  !$omp end do
827  end do
828  !$omp end parallel
829  end subroutine auto_prolong
830 
832  subroutine set_ref_info(tree, ref_flags, ref_info)
833  type(a3_t), intent(in) :: tree
834  integer, intent(in) :: ref_flags(:)
835  type(ref_info_t), intent(inout) :: ref_info
836  integer :: id, lvl, n, n_ch
837  integer, allocatable :: ref_count(:), drf_count(:)
838 
839  n_ch = a3_num_children
840  ref_info%n_add = n_ch * count(ref_flags == af_refine)
841  ref_info%n_rm = n_ch * count(ref_flags == af_derefine)
842 
843  ! Use highest_lvl+1 here because this lvl might have been completely removed
844  if (allocated(ref_info%lvls)) deallocate(ref_info%lvls)
845  allocate(ref_info%lvls(tree%highest_lvl+1))
846  allocate(ref_count(tree%highest_lvl+1))
847  allocate(drf_count(tree%highest_lvl+1))
848 
849  ! Find the number of (de)refined boxes per level
850  ref_count = 0
851  drf_count = 0
852 
853  do id = 1, size(ref_flags)
854  lvl = tree%boxes(id)%lvl
855 
856  if (ref_flags(id) == af_refine) then
857  ref_count(lvl) = ref_count(lvl) + 1
858  else if (ref_flags(id) == af_derefine) then
859  drf_count(lvl) = drf_count(lvl) + 1
860  end if
861  end do
862 
863  ! Allocate storage per level
864  ! There can be no new children at level 1
865  allocate(ref_info%lvls(1)%add(0))
866  allocate(ref_info%lvls(1)%rm(0))
867 
868  do lvl = 2, tree%highest_lvl+1
869  n = ref_count(lvl-1) * n_ch
870  allocate(ref_info%lvls(lvl)%add(n))
871  n = drf_count(lvl-1) * n_ch
872  allocate(ref_info%lvls(lvl)%rm(n))
873  end do
874 
875  ! Set the added and removed id's per level, these are the children of the
876  ! (de)refined boxes
877  ref_count = 0
878  drf_count = 0
879 
880  do id = 1, size(ref_flags)
881  lvl = tree%boxes(id)%lvl
882 
883  if (ref_flags(id) == af_refine) then
884  ref_count(lvl) = ref_count(lvl) + 1
885  n = n_ch * (ref_count(lvl)-1) + 1
886  ref_info%lvls(lvl+1)%add(n:n+n_ch-1) = tree%boxes(id)%children
887  else if (ref_flags(id) == af_derefine) then
888  drf_count(lvl) = drf_count(lvl) + 1
889  n = n_ch * (drf_count(lvl)-1) + 1
890  ref_info%lvls(lvl+1)%rm(n:n+n_ch-1) = tree%boxes(id)%children
891  end if
892  end do
893  end subroutine set_ref_info
894 
897  subroutine get_free_ids(tree, ids)
898  type(a3_t), intent(inout) :: tree
899  integer, intent(out) :: ids(:)
900  integer :: i, highest_id_prev, n_ids
901 
902  n_ids = size(ids)
903  !$omp critical (crit_free_ids)
904  highest_id_prev = tree%highest_id
905  tree%highest_id = tree%highest_id + n_ids
906  !$omp end critical (crit_free_ids)
907 
908  ids = [(highest_id_prev + i, i=1,n_ids)]
909  end subroutine get_free_ids
910 
916  subroutine consistent_ref_flags(tree, ref_flags, ref_subr, &
917  ref_buffer, keep_buffer)
918  use omp_lib, only: omp_get_max_threads, omp_get_thread_num
919  type(a3_t), intent(inout) :: tree
920  integer, intent(inout) :: ref_flags(:)
921  procedure(a3_subr_ref) :: ref_subr
922  integer, intent(in) :: ref_buffer
923  logical, intent(in) :: keep_buffer
924 
925  integer :: lvl, i, i_ch, ch_id, id, c_ids(a3_num_children)
926  integer :: nb, p_id, nb_id, p_nb_id
927  integer :: lvl_limit, thread_id
928  integer, allocatable :: tmp_flags(:, :)
929  integer :: cell_flags(tree%n_cell, tree%n_cell, tree%n_cell)
930  integer, parameter :: unset_flag = -huge(1)
931 
932  ! Set refinement flags for each thread individually, because we sometimes
933  ! modify the refinement flags of neighbors
934  allocate(tmp_flags(size(ref_flags), omp_get_max_threads()))
935 
936  lvl_limit = tree%lvl_limit
937  tmp_flags(:, :) = unset_flag
938 
939  ! Set refinement flags on all leaves and their immediate parents (on other
940  ! boxes the flags would not matter)
941 
942  !$omp parallel private(lvl, i, id, p_id, cell_flags, thread_id, i_ch, ch_id)
943  thread_id = omp_get_thread_num() + 1
944 
945  do lvl = 1, tree%highest_lvl
946  !$omp do
947  do i = 1, size(tree%lvls(lvl)%leaves)
948  id = tree%lvls(lvl)%leaves(i)
949 
950  call ref_subr(tree%boxes(id), cell_flags)
951  call cell_to_ref_flags(cell_flags, tree%n_cell, &
952  tmp_flags(:, thread_id), tree, id, ref_buffer, keep_buffer)
953 
954  ! If the parent exists, and this is the first child which is itself
955  ! not refined, set refinement flags for the parent
956  if (tree%boxes(id)%lvl > 1) then
957  p_id = tree%boxes(id)%parent
958  do i_ch = 1, a3_ix_to_ichild(tree%boxes(id)%ix)-1
959  ch_id = tree%boxes(p_id)%children(i_ch)
960  if (.not. a3_has_children(tree%boxes(ch_id))) exit
961  end do
962 
963  if (i_ch == a3_ix_to_ichild(tree%boxes(id)%ix)) then
964  ! This is the first child which is itself not refined
965  call ref_subr(tree%boxes(p_id), cell_flags)
966  call cell_to_ref_flags(cell_flags, tree%n_cell, &
967  tmp_flags(:, thread_id), tree, p_id, ref_buffer, keep_buffer)
968  end if
969  end if
970  end do
971  !$omp end do
972  end do
973  !$omp end parallel
974 
975  ! Take the highest value over the threads
976  do i = 1, size(ref_flags)
977  ref_flags(i) = maxval(tmp_flags(i, :))
978  if (ref_flags(i) == unset_flag) ref_flags(i) = af_keep_ref
979  end do
980 
981  if (maxval(ref_flags) > af_do_ref .or. minval(ref_flags) < af_rm_ref) &
982  stop "a3_adjust_refinement: invalid refinement flag given"
983 
984  ! Cannot refine beyond max level
985  do i = 1, size(tree%lvls(lvl_limit)%ids)
986  id = tree%lvls(lvl_limit)%ids(i)
987  if (ref_flags(id) == af_do_ref) ref_flags(id) = af_keep_ref
988  end do
989 
990  ! Ensure 2-1 balance
991  do lvl = tree%highest_lvl, 1, -1
992  do i = 1, size(tree%lvls(lvl)%leaves) ! We only check leaf tree%boxes
993  id = tree%lvls(lvl)%leaves(i)
994 
995  if (ref_flags(id) == af_do_ref .or. ref_flags(id) == af_refine) then
996  ref_flags(id) = af_refine ! Mark for actual refinement
997 
998  ! Ensure we will have the necessary neighbors
999  do nb = 1, a3_num_neighbors
1000  nb_id = tree%boxes(id)%neighbors(nb)
1001  if (nb_id == af_no_box) then
1002  ! Mark the parent containing neighbor for refinement
1003  p_id = tree%boxes(id)%parent
1004  p_nb_id = tree%boxes(p_id)%neighbors(nb)
1005  ref_flags(p_nb_id) = af_refine
1006  end if
1007  end do
1008 
1009  else if (ref_flags(id) == af_rm_ref) then
1010  ! Ensure we do not remove a required neighbor
1011  do nb = 1, a3_num_neighbors
1012  nb_id = tree%boxes(id)%neighbors(nb)
1013  if (nb_id > af_no_box) then
1014  if (a3_has_children(tree%boxes(nb_id)) .or. &
1015  ref_flags(nb_id) > af_keep_ref) then
1016  ref_flags(id) = af_keep_ref
1017  exit
1018  end if
1019  end if
1020  end do
1021  end if
1022 
1023  end do
1024  end do
1025 
1026  ! Make the (de)refinement flags consistent for blocks with children. Also
1027  ! ensure that at most one level can be removed at a time.
1028  do lvl = tree%highest_lvl-1, 1, -1
1029  do i = 1, size(tree%lvls(lvl)%parents)
1030  id = tree%lvls(lvl)%parents(i)
1031 
1032  ! Can only remove children if they are all marked for
1033  ! derefinement, and the box itself not for refinement.
1034  c_ids = tree%boxes(id)%children
1035  if (all(ref_flags(c_ids) == af_rm_ref) .and. &
1036  ref_flags(id) <= af_keep_ref) then
1037  ref_flags(id) = af_derefine
1038  else
1039  ref_flags(id) = af_keep_ref
1040  end if
1041  end do
1042  end do
1043 
1044  end subroutine consistent_ref_flags
1045 
1049  subroutine cell_to_ref_flags(cell_flags, nc, ref_flags, tree, id, &
1050  ref_buffer, keep_buffer)
1051  use m_a3_utils, only: a3_get_loc
1052  integer, intent(in) :: nc
1053  integer, intent(in) :: cell_flags(nc, nc, nc)
1054  integer, intent(inout) :: ref_flags(:)
1055  type(a3_t), intent(in) :: tree
1056  integer, intent(in) :: id
1057  integer, intent(in) :: ref_buffer
1058  logical, intent(in) :: keep_buffer
1059  integer :: ix0(3), ix1(3), i, j, k, nb_id
1060 
1061  if (minval(cell_flags) < af_rm_ref .or. &
1062  maxval(cell_flags) > af_do_ref) then
1063  error stop "Error: invalid cell flags given"
1064  end if
1065 
1066  ! Check whether the box needs to be refined or keep its refinement
1067  if (any(cell_flags == af_do_ref)) then
1068  ref_flags(id) = af_do_ref
1069  else if (any(cell_flags == af_keep_ref)) then
1070  ref_flags(id) = max(ref_flags(id), af_keep_ref)
1071  else ! All flags are af_rm_ref
1072  ref_flags(id) = max(ref_flags(id), af_rm_ref)
1073  end if
1074 
1075  if (ref_buffer <= 0) return ! No need to check neighbors
1076 
1077  ! Check whether neighbors also require refinement, which happens when cells
1078  ! close to the neighbor are flagged.
1079  do k = -1, 1; do j = -1, 1; do i = -1, 1
1080  if (all([i, j, k] == 0)) cycle
1081 
1082  nb_id = tree%boxes(id)%neighbor_mat(i, j, k)
1083 
1084  ! Skip neighbors that are not there
1085  if (nb_id <= af_no_box) cycle
1086 
1087  ! Compute index range relevant for neighbor
1088  ix0 = 1
1089  ix1 = nc
1090  where ([i, j, k] == 1)
1091  ix0 = nc - ref_buffer + 1
1092  ix1 = nc
1093  elsewhere ([i, j, k] == -1)
1094  ix0 = 1
1095  ix1 = ref_buffer
1096  end where
1097 
1098  if (any(cell_flags(ix0(1):ix1(1), ix0(2):ix1(2), &
1099  ix0(3):ix1(3)) == af_do_ref)) then
1100  ref_flags(nb_id) = af_do_ref
1101  else if (keep_buffer .and. any(cell_flags(ix0(1):ix1(1), &
1102  ix0(2):ix1(2), ix0(3):ix1(3)) == af_keep_ref)) then
1103  ref_flags(nb_id) = max(ref_flags(nb_id), af_keep_ref)
1104  end if
1105  end do; end do; end do
1106 
1107  end subroutine cell_to_ref_flags
1108 
1110  subroutine remove_children(boxes, id)
1111  type(box3_t), intent(inout) :: boxes(:)
1112  integer, intent(in) :: id
1113  integer :: ic, c_id, nb_id, nb_rev, nb, i, j, k
1114 
1115  do ic = 1, a3_num_children
1116  c_id = boxes(id)%children(ic)
1117 
1118  ! Remove from neighbors
1119  do nb = 1, a3_num_neighbors
1120  nb_id = boxes(c_id)%neighbors(nb)
1121  if (nb_id > af_no_box) then
1122  nb_rev = a3_neighb_rev(nb)
1123  boxes(nb_id)%neighbors(nb_rev) = af_no_box
1124  end if
1125  end do
1126 
1127  do k = -1, 1; do j = -1, 1; do i = -1, 1
1128  nb_id = boxes(c_id)%neighbor_mat(i, j, k)
1129  if (nb_id > af_no_box) then
1130  boxes(nb_id)%neighbor_mat(-i, -j, -k) = af_no_box
1131  end if
1132  end do; end do; end do
1133 
1134  call clear_box(boxes(c_id))
1135  end do
1136 
1137  boxes(id)%children = af_no_box
1138  end subroutine remove_children
1139 
1141  subroutine add_children(boxes, id, c_ids, n_cc, n_fc)
1142  type(box3_t), intent(inout) :: boxes(:)
1143  integer, intent(in) :: id
1144  integer, intent(in) :: c_ids(a3_num_children)
1145  integer, intent(in) :: n_cc
1146  integer, intent(in) :: n_fc
1147  integer :: i, nb, child_nb(2**(3-1))
1148  integer :: c_id, c_ix_base(3), dix(3)
1149 
1150  boxes(id)%children = c_ids
1151  c_ix_base = 2 * boxes(id)%ix - 1
1152 
1153  do i = 1, a3_num_children
1154  c_id = c_ids(i)
1155  boxes(c_id)%ix = c_ix_base + a3_child_dix(:,i)
1156  boxes(c_id)%lvl = boxes(id)%lvl+1
1157  boxes(c_id)%parent = id
1158  boxes(c_id)%tag = af_init_tag
1159  boxes(c_id)%children = af_no_box
1160  boxes(c_id)%neighbors = af_no_box
1161  boxes(c_id)%neighbor_mat = af_no_box
1162  boxes(c_id)%neighbor_mat(0, 0, 0) = c_id
1163  boxes(c_id)%n_cell = boxes(id)%n_cell
1164  boxes(c_id)%coord_t = boxes(id)%coord_t
1165  boxes(c_id)%dr = 0.5_dp * boxes(id)%dr
1166  boxes(c_id)%r_min = boxes(id)%r_min + 0.5_dp * boxes(id)%dr * &
1167  a3_child_dix(:,i) * boxes(id)%n_cell
1168 
1169  call a3_init_box(boxes(c_id), boxes(id)%n_cell, n_cc, n_fc)
1170  end do
1171 
1172  ! Set boundary conditions at children
1173  do nb = 1, a3_num_neighbors
1174  if (boxes(id)%neighbors(nb) < af_no_box) then
1175  child_nb = c_ids(a3_child_adj_nb(:, nb)) ! Neighboring children
1176  boxes(child_nb)%neighbors(nb) = boxes(id)%neighbors(nb)
1177  dix = a3_neighb_dix(:, nb)
1178  boxes(child_nb)%neighbor_mat(dix(1), dix(2), dix(3)) = &
1179  boxes(id)%neighbors(nb)
1180  end if
1181  end do
1182  end subroutine add_children
1183 
1186  subroutine set_child_ids(p_ids, c_ids, boxes)
1187  integer, intent(in) :: p_ids(:)
1188  integer, allocatable, intent(inout) :: c_ids(:)
1189  type(box3_t), intent(in) :: boxes(:)
1190  integer :: i, i0, i1, n_children
1191 
1192  n_children = a3_num_children * size(p_ids)
1193  if (n_children /= size(c_ids)) then
1194  deallocate(c_ids)
1195  allocate(c_ids(n_children))
1196  end if
1197 
1198  do i = 1, size(p_ids)
1199  i1 = i * a3_num_children
1200  i0 = i1 - a3_num_children + 1
1201  c_ids(i0:i1) = boxes(p_ids(i))%children
1202  end do
1203  end subroutine set_child_ids
1204 
1206  subroutine a3_consistent_fluxes(tree, f_ixs)
1207  type(a3_t), intent(inout) :: tree
1208  integer, intent(in) :: f_ixs(:)
1209  integer :: lvl, i, id, nb, nb_id
1210 
1211  if (.not. tree%ready) stop "Tree not ready"
1212  !$omp parallel private(lvl, i, id, nb, nb_id)
1213  do lvl = lbound(tree%lvls, 1), tree%highest_lvl-1
1214  !$omp do
1215  do i = 1, size(tree%lvls(lvl)%parents)
1216  id = tree%lvls(lvl)%parents(i)
1217  do nb = 1, a3_num_neighbors
1218  nb_id = tree%boxes(id)%neighbors(nb)
1219 
1220  ! If the neighbor exists and has no children, set flux
1221  if (nb_id > af_no_box) then
1222  if (.not. a3_has_children(tree%boxes(nb_id))) then
1223  call flux_from_children(tree%boxes, id, nb, f_ixs)
1224  end if
1225  end if
1226  end do
1227  end do
1228  !$omp end do
1229  end do
1230  !$omp end parallel
1231  end subroutine a3_consistent_fluxes
1232 
1235  subroutine flux_from_children(boxes, id, nb, f_ixs)
1236  type(box3_t), intent(inout) :: boxes(:)
1237  integer, intent(in) :: id
1238  integer, intent(in) :: nb
1239  integer, intent(in) :: f_ixs(:)
1240  integer :: nc, nch, c_id, i_ch, i, ic, d
1241  integer :: n_chnb, nb_id, i_nb, ioff(3)
1242 
1243  nc = boxes(id)%n_cell
1244  nch = ishft(nc, -1) ! nc/2
1245  d = a3_neighb_dim(nb)
1246  n_chnb = 2**(3-1)
1247  nb_id = boxes(id)%neighbors(nb)
1248 
1249  if (a3_neighb_low(nb)) then
1250  i = 1
1251  i_nb = nc+1
1252  else
1253  i = nc+1
1254  i_nb = 1
1255  end if
1256 
1257  select case (d)
1258  case (1)
1259  do ic = 1, n_chnb
1260  i_ch = a3_child_adj_nb(ic, nb)
1261  c_id = boxes(id)%children(i_ch)
1262  ioff = nch*a3_child_dix(:, i_ch)
1263  boxes(nb_id)%fc(i_nb, ioff(2)+1:ioff(2)+nch, &
1264  ioff(3)+1:ioff(3)+nch, 1, f_ixs) = 0.25_dp * ( &
1265  boxes(c_id)%fc(i, 1:nc:2, 1:nc:2, 1, f_ixs) + &
1266  boxes(c_id)%fc(i, 2:nc:2, 1:nc:2, 1, f_ixs) + &
1267  boxes(c_id)%fc(i, 1:nc:2, 2:nc:2, 1, f_ixs) + &
1268  boxes(c_id)%fc(i, 2:nc:2, 2:nc:2, 1, f_ixs))
1269  end do
1270  case (2)
1271  do ic = 1, n_chnb
1272  i_ch = a3_child_adj_nb(ic, nb)
1273  c_id = boxes(id)%children(i_ch)
1274  ioff = nch*a3_child_dix(:, i_ch)
1275  boxes(nb_id)%fc(ioff(1)+1:ioff(1)+nch, i_nb, &
1276  ioff(3)+1:ioff(3)+nch, 2, f_ixs) = 0.25_dp * ( &
1277  boxes(c_id)%fc(1:nc:2, i, 1:nc:2, 2, f_ixs) + &
1278  boxes(c_id)%fc(2:nc:2, i, 1:nc:2, 2, f_ixs) + &
1279  boxes(c_id)%fc(1:nc:2, i, 2:nc:2, 2, f_ixs) + &
1280  boxes(c_id)%fc(2:nc:2, i, 2:nc:2, 2, f_ixs))
1281  end do
1282  case (3)
1283  do ic = 1, n_chnb
1284  i_ch = a3_child_adj_nb(ic, nb)
1285  c_id = boxes(id)%children(i_ch)
1286  ioff = nch*a3_child_dix(:, i_ch)
1287  boxes(nb_id)%fc(ioff(1)+1:ioff(1)+nch, &
1288  ioff(2)+1:ioff(2)+nch, i_nb, 3, f_ixs) = 0.25_dp * ( &
1289  boxes(c_id)%fc(1:nc:2, 1:nc:2, i, 3, f_ixs) + &
1290  boxes(c_id)%fc(2:nc:2, 1:nc:2, i, 3, f_ixs) + &
1291  boxes(c_id)%fc(1:nc:2, 2:nc:2, i, 3, f_ixs) + &
1292  boxes(c_id)%fc(2:nc:2, 2:nc:2, i, 3, f_ixs))
1293  end do
1294  end select
1295  end subroutine flux_from_children
1296 
1297 end module m_a3_core
This module contains all kinds of different &#39;helper&#39; routines for Afivo. If the number of routines fo...
Definition: m_a3_utils.f90:5
Subroutine for restriction.
Definition: m_a3_types.f90:246
This module contains routines for restriction: going from fine to coarse variables.
To fill ghost cells near physical boundaries.
Definition: m_a3_types.f90:216
Subroutine for prolongation.
Definition: m_a3_types.f90:236
This module contains methods to convert indices to morton numbers.
Definition: m_morton.f90:7
The basic building block of afivo: a box with cell-centered and face centered data, and information about its position, neighbors, children etc.
Definition: m_a3_types.f90:109
Type which stores all the boxes and levels, as well as some information about the number of boxes...
Definition: m_a3_types.f90:131
To fill ghost cells near refinement boundaries.
Definition: m_a3_types.f90:207
This module contains the core routines of Afivo, namely those that deal with initializing and changin...
Definition: m_a3_core.f90:4
This module contains the routines related to prolongation: going from coarse to fine variables...
Definition: m_a3_prolong.f90:4
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
This module contains the basic types and constants that are used in the 3-dimensional version of Afiv...
Definition: m_a3_types.f90:5