Afivo  0.3
m_a3_utils.f90
1 
5 module m_a3_utils
6  use m_a3_types
7 
8  implicit none
9  private
10 
11  ! Public subroutines
12  public :: a3_loop_box
13  public :: a3_loop_box_arg
14  public :: a3_loop_boxes
15  public :: a3_loop_boxes_arg
16  public :: a3_tree_clear_cc
17  public :: a3_box_clear_cc
18  public :: a3_tree_clear_ghostcells
19  public :: a3_box_clear_ghostcells
20  public :: a3_box_add_cc
21  public :: a3_box_sub_cc
22  public :: a3_tree_times_cc
23  public :: a3_tree_apply
24  public :: a3_box_times_cc
25  public :: a3_box_lincomb_cc
26  public :: a3_box_copy_cc_to
27  public :: a3_box_copy_cc
28  public :: a3_boxes_copy_cc
29  public :: a3_tree_copy_cc
30  public :: a3_reduction
31  public :: a3_reduction_vec
32  public :: a3_reduction_loc
33  public :: a3_tree_max_cc
34  public :: a3_tree_min_cc
35  public :: a3_tree_max_fc
36  public :: a3_tree_sum_cc
37  public :: a3_box_copy_fc
38  public :: a3_boxes_copy_fc
39  public :: a3_tree_copy_fc
40 
41  ! Public functions
42  public :: a3_get_id_at
43  public :: a3_get_loc
44  public :: a3_r_loc
45  public :: a3_r_inside
46  public :: a3_n_cell
47 
48 contains
49 
51  subroutine a3_loop_box(tree, my_procedure, leaves_only)
52  type(a3_t), intent(inout) :: tree
53  procedure(a3_subr) :: my_procedure
54  logical, intent(in), optional :: leaves_only
55  logical :: leaves
56  integer :: lvl, i, id
57 
58  leaves = .false.; if (present(leaves_only)) leaves = leaves_only
59  if (.not. tree%ready) stop "a3_loop_box: set_base has not been called"
60 
61  !$omp parallel private(lvl, i, id)
62  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
63  if (leaves) then
64  !$omp do
65  do i = 1, size(tree%lvls(lvl)%leaves)
66  id = tree%lvls(lvl)%leaves(i)
67  call my_procedure(tree%boxes(id))
68  end do
69  !$omp end do
70  else
71  !$omp do
72  do i = 1, size(tree%lvls(lvl)%ids)
73  id = tree%lvls(lvl)%ids(i)
74  call my_procedure(tree%boxes(id))
75  end do
76  !$omp end do
77  end if
78  end do
79  !$omp end parallel
80  end subroutine a3_loop_box
81 
83  subroutine a3_loop_box_arg(tree, my_procedure, rarg, leaves_only)
84  type(a3_t), intent(inout) :: tree
85  procedure(a3_subr_arg) :: my_procedure
86  real(dp), intent(in) :: rarg(:)
87  logical, intent(in), optional :: leaves_only
88  logical :: leaves
89  integer :: lvl, i, id
90 
91  if (.not. tree%ready) stop "Tree not ready"
92  leaves = .false.; if (present(leaves_only)) leaves = leaves_only
93 
94  !$omp parallel private(lvl, i, id)
95  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
96  if (leaves) then
97  !$omp do
98  do i = 1, size(tree%lvls(lvl)%leaves)
99  id = tree%lvls(lvl)%leaves(i)
100  call my_procedure(tree%boxes(id), rarg)
101  end do
102  !$omp end do
103  else
104  !$omp do
105  do i = 1, size(tree%lvls(lvl)%ids)
106  id = tree%lvls(lvl)%ids(i)
107  call my_procedure(tree%boxes(id), rarg)
108  end do
109  !$omp end do
110  end if
111  end do
112  !$omp end parallel
113  end subroutine a3_loop_box_arg
114 
116  subroutine a3_loop_boxes(tree, my_procedure, leaves_only)
117  type(a3_t), intent(inout) :: tree
118  procedure(a3_subr_boxes) :: my_procedure
119  logical, intent(in), optional :: leaves_only
120  logical :: leaves
121  integer :: lvl, i, id
122 
123  if (.not. tree%ready) stop "Tree not ready"
124  leaves = .false.; if (present(leaves_only)) leaves = leaves_only
125 
126  !$omp parallel private(lvl, i, id)
127  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
128  if (leaves) then
129  !$omp do
130  do i = 1, size(tree%lvls(lvl)%leaves)
131  id = tree%lvls(lvl)%leaves(i)
132  call my_procedure(tree%boxes, id)
133  end do
134  !$omp end do
135  else
136  !$omp do
137  do i = 1, size(tree%lvls(lvl)%ids)
138  id = tree%lvls(lvl)%ids(i)
139  call my_procedure(tree%boxes, id)
140  end do
141  !$omp end do
142  end if
143  end do
144  !$omp end parallel
145  end subroutine a3_loop_boxes
146 
148  subroutine a3_loop_boxes_arg(tree, my_procedure, rarg, leaves_only)
149  type(a3_t), intent(inout) :: tree
150  procedure(a3_subr_boxes_arg) :: my_procedure
151  real(dp), intent(in) :: rarg(:)
152  logical, intent(in), optional :: leaves_only
153  logical :: leaves
154  integer :: lvl, i, id
155 
156  if (.not. tree%ready) stop "Tree not ready"
157  leaves = .false.; if (present(leaves_only)) leaves = leaves_only
158 
159  !$omp parallel private(lvl, i, id)
160  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
161  if (leaves) then
162  !$omp do
163  do i = 1, size(tree%lvls(lvl)%leaves)
164  id = tree%lvls(lvl)%leaves(i)
165  call my_procedure(tree%boxes, id, rarg)
166  end do
167  !$omp end do
168  else
169  !$omp do
170  do i = 1, size(tree%lvls(lvl)%ids)
171  id = tree%lvls(lvl)%ids(i)
172  call my_procedure(tree%boxes, id, rarg)
173  end do
174  !$omp end do
175  end if
176  end do
177  !$omp end parallel
178  end subroutine a3_loop_boxes_arg
179 
181  pure function a3_r_inside(box, r, d) result(inside)
182  type(box3_t), intent(in) :: box
183  real(dp), intent(in) :: r(3)
184  real(dp), intent(in), optional :: d
185  real(dp) :: r_max(3)
186  logical :: inside
187 
188  r_max = box%r_min + box%dr * box%n_cell
189  if (present(d)) then
190  inside = all(r+d >= box%r_min) .and. all(r-d <= r_max)
191  else
192  inside = all(r >= box%r_min) .and. all(r <= r_max)
193  end if
194  end function a3_r_inside
195 
199  pure function a3_get_id_at(tree, rr, highest_lvl, guess) result(id)
200  type(a3_t), intent(in) :: tree
201  real(dp), intent(in) :: rr(3)
202  integer, intent(in), optional :: highest_lvl
203  integer, intent(in), optional :: guess
204  integer :: id
205 
206  integer :: i, i_ch, lvl_max
207 
208  if (present(guess)) then
209  id = guess
210  if (id > 0 .and. id < tree%highest_id) then
211  if (tree%boxes(id)%in_use .and. &
212  a3_r_inside(tree%boxes(id), rr)) then
213 
214  ! Jump into potential children for as long as possible
215  do while (a3_has_children(tree%boxes(id)))
216  i_ch = child_that_contains(tree%boxes(id), rr)
217  id = tree%boxes(id)%children(i_ch)
218  end do
219 
220  ! return ! Exit routine
221  end if
222  end if
223  end if
224 
225  lvl_max = tree%lvl_limit
226  if (present(highest_lvl)) lvl_max = highest_lvl
227 
228  ! Find lvl 1 box that includes rr
229  do i = 1, size(tree%lvls(1)%ids)
230  id = tree%lvls(1)%ids(i)
231  if (a3_r_inside(tree%boxes(id), rr)) exit
232  end do
233 
234  if (i > size(tree%lvls(1)%ids)) then
235  ! Not inside any box
236  id = -1
237  else
238  ! Jump into children for as long as possible
239  do
240  if (tree%boxes(id)%lvl >= lvl_max .or. &
241  .not. a3_has_children(tree%boxes(id))) exit
242  i_ch = child_that_contains(tree%boxes(id), rr)
243  id = tree%boxes(id)%children(i_ch)
244  end do
245  end if
246 
247  end function a3_get_id_at
248 
252  pure function a3_get_loc(tree, rr, highest_lvl, guess) result(loc)
253  type(a3_t), intent(in) :: tree
254  real(dp), intent(in) :: rr(3)
255  integer, intent(in), optional :: highest_lvl
257  integer, intent(in), optional :: guess
258  type(a3_loc_t) :: loc
259 
260  loc%id = a3_get_id_at(tree, rr, highest_lvl, guess)
261 
262  if (loc%id == -1) then
263  loc%ix = -1
264  else
265  loc%ix = a3_cc_ix(tree%boxes(loc%id), rr)
266 
267  ! Fix indices for points exactly on the boundaries of a box (which could
268  ! get a ghost cell index)
269  where (loc%ix < 1) loc%ix = 1
270  where (loc%ix > tree%n_cell) loc%ix = tree%n_cell
271  end if
272  end function a3_get_loc
273 
275  pure function child_that_contains(box, rr) result(i_ch)
276  type(box3_t), intent(in) :: box
277  real(dp), intent(in) :: rr(3)
278  integer :: i_ch
279  real(dp) :: center(3)
280 
281  i_ch = 1
282  center = box%r_min + box%dr * ishft(box%n_cell, -1)
283 
284  if (rr(1) > center(1)) i_ch = i_ch + 1
285  if (rr(2) > center(2)) i_ch = i_ch + 2
286  if (rr(3) > center(3)) i_ch = i_ch + 4
287  end function child_that_contains
288 
290  pure function a3_r_loc(tree, loc) result(r)
291  type(a3_t), intent(in) :: tree
292  type(a3_loc_t), intent(in) :: loc
293  real(dp) :: r(3)
294  r = tree%boxes(loc%id)%r_min + &
295  (loc%ix-0.5_dp) * tree%boxes(loc%id)%dr
296  end function a3_r_loc
297 
298  subroutine a3_tree_clear_cc(tree, iv)
299  type(a3_t), intent(inout) :: tree
300  integer, intent(in) :: iv
301  integer :: lvl, i, id
302 
303  !$omp parallel private(lvl, i, id)
304  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
305  !$omp do
306  do i = 1, size(tree%lvls(lvl)%ids)
307  id = tree%lvls(lvl)%ids(i)
308  call a3_box_clear_cc(tree%boxes(id), iv)
309  end do
310  !$omp end do
311  end do
312  !$omp end parallel
313  end subroutine a3_tree_clear_cc
314 
316  subroutine a3_box_clear_cc(box, iv)
317  type(box3_t), intent(inout) :: box
318  integer, intent(in) :: iv
319  box%cc(:,:,:, iv) = 0
320  end subroutine a3_box_clear_cc
321 
322  subroutine a3_tree_clear_ghostcells(tree, iv)
323  type(a3_t), intent(inout) :: tree
324  integer, intent(in) :: iv
325  integer :: lvl, i, id
326 
327  !$omp parallel private(lvl, i, id)
328  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
329  !$omp do
330  do i = 1, size(tree%lvls(lvl)%ids)
331  id = tree%lvls(lvl)%ids(i)
332  call a3_box_clear_ghostcells(tree%boxes(id), iv)
333  end do
334  !$omp end do
335  end do
336  !$omp end parallel
337  end subroutine a3_tree_clear_ghostcells
338 
340  subroutine a3_box_clear_ghostcells(box, iv)
341  type(box3_t), intent(inout) :: box
342  integer, intent(in) :: iv
343  integer :: nc
344 
345  nc = box%n_cell
346 
347  box%cc(0, :, :, iv) = 0
348  box%cc(nc+1, :, :, iv) = 0
349  box%cc(:, 0, :, iv) = 0
350  box%cc(:, nc+1, :, iv) = 0
351  box%cc(:, :, 0, iv) = 0
352  box%cc(:, :, nc+1, iv) = 0
353  end subroutine a3_box_clear_ghostcells
354 
356  subroutine a3_box_add_cc(box, iv_from, iv_to)
357  type(box3_t), intent(inout) :: box
358  integer, intent(in) :: iv_from, iv_to
359  box%cc(:,:,:, iv_to) = box%cc(:,:,:, iv_to) + box%cc(:,:,:, iv_from)
360  end subroutine a3_box_add_cc
361 
363  subroutine a3_box_sub_cc(box, iv_from, iv_to)
364  type(box3_t), intent(inout) :: box
365  integer, intent(in) :: iv_from, iv_to
366  box%cc(:,:,:, iv_to) = box%cc(:,:,:, iv_to) - box%cc(:,:,:, iv_from)
367  end subroutine a3_box_sub_cc
368 
371  subroutine a3_tree_apply(tree, iv_a, iv_b, op, eps)
372  type(a3_t), intent(inout) :: tree
373  integer, intent(in) :: iv_a, iv_b
374  character(len=*), intent(in) :: op
375  real(dp), intent(in), optional :: eps
376  integer :: lvl, i, id
377  real(dp) :: use_eps
378 
379  use_eps = sqrt(tiny(1.0_dp))
380  if (present(eps)) use_eps = eps
381 
382  !$omp parallel private(lvl, i, id)
383  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
384  !$omp do
385  do i = 1, size(tree%lvls(lvl)%ids)
386  id = tree%lvls(lvl)%ids(i)
387  select case (op)
388  case ('+')
389  tree%boxes(id)%cc(:, :, :, iv_a) = &
390  tree%boxes(id)%cc(:, :, :, iv_a) + &
391  tree%boxes(id)%cc(:, :, :, iv_b)
392  case ('*')
393  tree%boxes(id)%cc(:, :, :, iv_a) = &
394  tree%boxes(id)%cc(:, :, :, iv_a) * &
395  tree%boxes(id)%cc(:, :, :, iv_b)
396  case ('/')
397  tree%boxes(id)%cc(:, :, :, iv_a) = &
398  tree%boxes(id)%cc(:, :, :, iv_a) / &
399  max(tree%boxes(id)%cc(:, :, :, iv_b), eps)
400  case default
401  error stop "a3_tree_apply: unknown operand"
402  end select
403  end do
404  !$omp end do
405  end do
406  !$omp end parallel
407  end subroutine a3_tree_apply
408 
409  subroutine a3_tree_times_cc(tree, ivs, facs)
410  type(a3_t), intent(inout) :: tree
411  integer, intent(in) :: ivs(:)
412  real(dp), intent(in) :: facs(:)
413  integer :: lvl, i, id, n
414 
415  if (size(ivs) /= size(facs)) &
416  error stop "a3_times_cc: invalid array size"
417 
418  !$omp parallel private(lvl, i, id, n)
419  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
420  !$omp do
421  do i = 1, size(tree%lvls(lvl)%ids)
422  id = tree%lvls(lvl)%ids(i)
423  do n = 1, size(ivs)
424  call a3_box_times_cc(tree%boxes(id), facs(n), ivs(n))
425  end do
426  end do
427  !$omp end do
428  end do
429  !$omp end parallel
430  end subroutine a3_tree_times_cc
431 
433  subroutine a3_box_times_cc(box, a, iv)
434  type(box3_t), intent(inout) :: box
435  real(dp), intent(in) :: a
436  integer, intent(in) :: iv
437  box%cc(:,:,:, iv) = a * box%cc(:,:,:, iv)
438  end subroutine a3_box_times_cc
439 
441  subroutine a3_box_lincomb_cc(box, a, iv_a, b, iv_b)
442  type(box3_t), intent(inout) :: box
443  real(dp), intent(in) :: a, b
444  integer, intent(in) :: iv_a, iv_b
445  box%cc(:,:,:, iv_b) = a * box%cc(:,:,:, iv_a) + b * box%cc(:,:,:, iv_b)
446  end subroutine a3_box_lincomb_cc
447 
449  subroutine a3_box_copy_cc_to(box_from, iv_from, box_to, iv_to)
450  type(box3_t), intent(in) :: box_from
451  type(box3_t), intent(inout) :: box_to
452  integer, intent(in) :: iv_from, iv_to
453  box_to%cc(:,:,:, iv_to) = box_from%cc(:,:,:, iv_from)
454  end subroutine a3_box_copy_cc_to
455 
457  subroutine a3_box_copy_cc(box, iv_from, iv_to)
458  type(box3_t), intent(inout) :: box
459  integer, intent(in) :: iv_from, iv_to
460  box%cc(:,:,:, iv_to) = box%cc(:,:,:, iv_from)
461  end subroutine a3_box_copy_cc
462 
464  subroutine a3_boxes_copy_cc(boxes, ids, iv_from, iv_to)
465  type(box3_t), intent(inout) :: boxes(:)
466  integer, intent(in) :: ids(:), iv_from, iv_to
467  integer :: i
468 
469  !$omp parallel do
470  do i = 1, size(ids)
471  call a3_box_copy_cc(boxes(ids(i)), iv_from, iv_to)
472  end do
473  !$omp end parallel do
474  end subroutine a3_boxes_copy_cc
475 
477  subroutine a3_tree_copy_cc(tree, iv_from, iv_to)
478  type(a3_t), intent(inout) :: tree
479  integer, intent(in) :: iv_from, iv_to
480  integer :: lvl
481 
482  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
483  call a3_boxes_copy_cc(tree%boxes, tree%lvls(lvl)%ids, iv_from, iv_to)
484  end do
485  end subroutine a3_tree_copy_cc
486 
488  subroutine a3_reduction(tree, box_func, reduction, init_val, out_val)
489  type(a3_t), intent(in) :: tree
490  real(dp), intent(in) :: init_val
491  real(dp), intent(out) :: out_val
492  real(dp) :: tmp, my_val
493  integer :: i, id, lvl
494 
495  interface
496 
497  real(dp) function box_func(box)
498  import
499  type(box3_t), intent(in) :: box
500  end function box_func
501 
503  real(dp) function reduction(a, b)
504  import
505  real(dp), intent(in) :: a, b
506  end function reduction
507  end interface
508 
509  if (.not. tree%ready) stop "Tree not ready"
510  out_val = init_val
511  my_val = init_val
512 
513  !$omp parallel private(lvl, i, id, tmp) firstprivate(my_val)
514  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
515  !$omp do
516  do i = 1, size(tree%lvls(lvl)%leaves)
517  id = tree%lvls(lvl)%leaves(i)
518  tmp = box_func(tree%boxes(id))
519  my_val = reduction(tmp, my_val)
520  end do
521  !$omp end do
522  end do
523 
524  !$omp critical
525  out_val = reduction(my_val, out_val)
526  !$omp end critical
527  !$omp end parallel
528  end subroutine a3_reduction
529 
531  subroutine a3_reduction_vec(tree, box_func, reduction, init_val, &
532  out_val, n_vals)
533  type(a3_t), intent(in) :: tree
534  integer, intent(in) :: n_vals
535  real(dp), intent(in) :: init_val(n_vals)
536  real(dp), intent(out) :: out_val(n_vals)
537  real(dp) :: tmp(n_vals), my_val(n_vals)
538  integer :: i, id, lvl
539 
540  interface
541 
542  function box_func(box, n_vals) result(vec)
543  import
544  type(box3_t), intent(in) :: box
545  integer, intent(in) :: n_vals
546  real(dp) :: vec(n_vals)
547  end function box_func
548 
550  function reduction(vec_1, vec_2, n_vals) result(vec)
551  import
552  integer, intent(in) :: n_vals
553  real(dp), intent(in) :: vec_1(n_vals), vec_2(n_vals)
554  real(dp) :: vec(n_vals)
555  end function reduction
556  end interface
557 
558  if (.not. tree%ready) stop "Tree not ready"
559  out_val = init_val
560  my_val = init_val
561 
562  !$omp parallel private(lvl, i, id, tmp) firstprivate(my_val)
563  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
564  !$omp do
565  do i = 1, size(tree%lvls(lvl)%leaves)
566  id = tree%lvls(lvl)%leaves(i)
567  tmp = box_func(tree%boxes(id), n_vals)
568  my_val = reduction(tmp, my_val, n_vals)
569  end do
570  !$omp end do
571  end do
572 
573  !$omp critical
574  out_val = reduction(my_val, out_val, n_vals)
575  !$omp end critical
576  !$omp end parallel
577  end subroutine a3_reduction_vec
578 
581  subroutine a3_reduction_loc(tree, iv, box_subr, reduction, &
582  init_val, out_val, out_loc)
583  type(a3_t), intent(in) :: tree
584  integer, intent(in) :: iv
585  real(dp), intent(in) :: init_val
586  real(dp), intent(out) :: out_val
587  type(a3_loc_t), intent(out) :: out_loc
588  real(dp) :: tmp, new_val, my_val
589  integer :: i, id, lvl, tmp_ix(3)
590  type(a3_loc_t) :: my_loc
591 
592  interface
593 
594  subroutine box_subr(box, iv, val, ix)
595  import
596  type(box3_t), intent(in) :: box
597  integer, intent(in) :: iv
598  real(dp), intent(out) :: val
599  integer, intent(out) :: ix(3)
600  end subroutine box_subr
601 
603  real(dp) function reduction(a, b)
604  import
605  real(dp), intent(in) :: a, b
606  end function reduction
607  end interface
608 
609  if (.not. tree%ready) stop "Tree not ready"
610  out_val = init_val
611  my_val = init_val
612  my_loc%id = -1
613  my_loc%ix = -1
614 
615  !$omp parallel private(lvl, i, id, tmp, tmp_ix, new_val) &
616  !$omp firstprivate(my_val, my_loc)
617  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
618  !$omp do
619  do i = 1, size(tree%lvls(lvl)%leaves)
620  id = tree%lvls(lvl)%leaves(i)
621  call box_subr(tree%boxes(id), iv, tmp, tmp_ix)
622  new_val = reduction(tmp, my_val)
623  if (abs(new_val - my_val) > 0) then
624  my_loc%id = id
625  my_loc%ix = tmp_ix
626  my_val = tmp
627  end if
628  end do
629  !$omp end do
630  end do
631 
632  !$omp critical
633  new_val = reduction(my_val, out_val)
634  if (abs(new_val - out_val) > 0) then
635  out_loc%id = my_loc%id
636  out_loc%ix = my_loc%ix
637  out_val = my_val
638  end if
639  !$omp end critical
640  !$omp end parallel
641  end subroutine a3_reduction_loc
642 
645  subroutine a3_tree_max_cc(tree, iv, cc_max, loc)
646  type(a3_t), intent(in) :: tree
647  integer, intent(in) :: iv
648  real(dp), intent(out) :: cc_max
650  type(a3_loc_t), intent(out), optional :: loc
651  type(a3_loc_t) :: tmp_loc
652 
653  call a3_reduction_loc(tree, iv, box_max_cc, reduce_max, &
654  -huge(1.0_dp)/10, cc_max, tmp_loc)
655  if (present(loc)) loc = tmp_loc
656  end subroutine a3_tree_max_cc
657 
660  subroutine a3_tree_min_cc(tree, iv, cc_min, loc)
661  type(a3_t), intent(in) :: tree
662  integer, intent(in) :: iv
663  real(dp), intent(out) :: cc_min
665  type(a3_loc_t), intent(out), optional :: loc
666  type(a3_loc_t) :: tmp_loc
667 
668  call a3_reduction_loc(tree, iv, box_min_cc, reduce_min, &
669  huge(1.0_dp)/10, cc_min, tmp_loc)
670 
671  if (present(loc)) loc = tmp_loc
672  end subroutine a3_tree_min_cc
673 
675  subroutine a3_tree_max_fc(tree, dim, iv, fc_max, loc)
676  type(a3_t), intent(in) :: tree
677  integer, intent(in) :: dim
678  integer, intent(in) :: iv
679  real(dp), intent(out) :: fc_max
681  type(a3_loc_t), intent(out), optional :: loc
682  type(a3_loc_t) :: tmp_loc
683  integer :: dim_iv
684 
685  ! Encode dim and iv in a single variable
686  dim_iv = (dim-1) * tree%n_var_face + iv - 1
687 
688  call a3_reduction_loc(tree, dim_iv, box_max_fc, reduce_max, &
689  -huge(1.0_dp)/10, fc_max, tmp_loc)
690  if (present(loc)) loc = tmp_loc
691  end subroutine a3_tree_max_fc
692 
693  subroutine box_max_cc(box, iv, val, ix)
694  type(box3_t), intent(in) :: box
695  integer, intent(in) :: iv
696  real(dp), intent(out) :: val
697  integer, intent(out) :: ix(3)
698  integer :: nc
699 
700  nc = box%n_cell
701  ix = maxloc(box%cc(1:nc, 1:nc, 1:nc, iv))
702  val = box%cc(ix(1), ix(2), ix(3), iv)
703  end subroutine box_max_cc
704 
705  subroutine box_min_cc(box, iv, val, ix)
706  type(box3_t), intent(in) :: box
707  integer, intent(in) :: iv
708  real(dp), intent(out) :: val
709  integer, intent(out) :: ix(3)
710  integer :: nc
711 
712  nc = box%n_cell
713  ix = minloc(box%cc(1:nc, 1:nc, 1:nc, iv))
714  val = box%cc(ix(1), ix(2), ix(3), iv)
715  end subroutine box_min_cc
716 
717  subroutine box_max_fc(box, dim_iv, val, ix)
718  type(box3_t), intent(in) :: box
719  integer, intent(in) :: dim_iv
720  real(dp), intent(out) :: val
721  integer, intent(out) :: ix(3)
722  integer :: dim, iv, n_fc, nc, dix(3)
723 
724  nc = box%n_cell
725  dix(:) = 0
726 
727  n_fc = size(box%fc, 5)
728  ! Decode dim and iv
729  dim = dim_iv / n_fc + 1
730  dix(dim) = 1
731  iv = dim_iv - (dim-1) * n_fc + 1
732  ix = maxloc(box%fc(1:nc+dix(1), 1:nc+dix(2), 1:nc+dix(3), dim, iv))
733  val = box%fc(ix(1), ix(2), ix(3), dim, iv)
734  end subroutine box_max_fc
735 
736  real(dp) function reduce_max(a, b)
737  real(dp), intent(in) :: a, b
738  reduce_max = max(a, b)
739  end function reduce_max
740 
741  real(dp) function reduce_min(a, b)
742  real(dp), intent(in) :: a, b
743  reduce_min = min(a, b)
744  end function reduce_min
745 
748  subroutine a3_tree_sum_cc(tree, iv, cc_sum)
749  type(a3_t), intent(in) :: tree
750  integer, intent(in) :: iv
751  real(dp), intent(out) :: cc_sum
752  real(dp) :: tmp, my_sum, fac
753  integer :: i, id, lvl, nc
754 
755  if (.not. tree%ready) stop "Tree not ready"
756  my_sum = 0
757 
758  !$omp parallel reduction(+: my_sum) private(lvl, i, id, nc, tmp, fac)
759  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
760  fac = a3_lvl_dr(tree, lvl)**3
761 
762  !$omp do
763  do i = 1, size(tree%lvls(lvl)%leaves)
764  id = tree%lvls(lvl)%leaves(i)
765  nc = tree%boxes(id)%n_cell
766  tmp = sum(tree%boxes(id)%cc(1:nc, 1:nc, 1:nc, iv))
767  my_sum = my_sum + fac * tmp
768  end do
769  !$omp end do
770  end do
771  !$omp end parallel
772 
773  cc_sum = my_sum
774 
775  end subroutine a3_tree_sum_cc
776 
778  subroutine a3_box_copy_fc(box, iv_from, iv_to)
779  type(box3_t), intent(inout) :: box
780  integer, intent(in) :: iv_from
781  integer, intent(in) :: iv_to
782  box%fc(:,:,:,:, iv_to) = box%fc(:,:,:,:, iv_from)
783  end subroutine a3_box_copy_fc
784 
786  subroutine a3_boxes_copy_fc(boxes, ids, iv_from, iv_to)
787  type(box3_t), intent(inout) :: boxes(:)
788  integer, intent(in) :: ids(:)
789  integer, intent(in) :: iv_from
790  integer, intent(in) :: iv_to
791  integer :: i
792 
793  !$omp parallel do
794  do i = 1, size(ids)
795  call a3_box_copy_fc(boxes(ids(i)), iv_from, iv_to)
796  end do
797  !$omp end parallel do
798  end subroutine a3_boxes_copy_fc
799 
801  subroutine a3_tree_copy_fc(tree, iv_from, iv_to)
802  type(a3_t), intent(inout) :: tree
803  integer, intent(in) :: iv_from
804  integer, intent(in) :: iv_to
805  integer :: lvl
806 
807  if (.not. tree%ready) stop "Tree not ready"
808  do lvl = lbound(tree%lvls, 1), tree%highest_lvl
809  call a3_boxes_copy_fc(tree%boxes, tree%lvls(lvl)%ids, iv_from, iv_to)
810  end do
811  end subroutine a3_tree_copy_fc
812 
816  pure function a3_n_cell(tree, lvl) result(n_cell)
817  type(a3_t), intent(in) :: tree
818  integer, intent(in) :: lvl
819  integer :: n_cell
820 
821  if (lvl >= 1) then
822  n_cell = tree%n_cell
823  else
824  n_cell = tree%n_cell / (2**(1-lvl))
825  end if
826  end function a3_n_cell
827 
828 end module m_a3_utils
Type specifying the location of a cell.
Definition: m_a3_types.f90:162
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 that gets a list of boxes, an id and an array of reals.
Definition: m_a3_types.f90:199
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
Subroutine that gets a box.
Definition: m_a3_types.f90:179
Subroutine that gets a box and an array of reals.
Definition: m_a3_types.f90:185
Subroutine that gets a list of boxes and a box id.
Definition: m_a3_types.f90:192
This module contains the basic types and constants that are used in the 3-dimensional version of Afiv...
Definition: m_a3_types.f90:5