4 #include "cpp_macros.h"
10 public :: af_restrict_to_box
11 public :: af_restrict_to_boxes
12 public :: af_restrict_tree
13 public :: af_restrict_box
14 public :: af_restrict_ref_boundary
21 subroutine af_restrict_to_box(boxes, id, ivs)
22 type(
box_t),
intent(inout) :: boxes(:)
23 integer,
intent(in) :: id
24 integer,
intent(in) :: ivs(:)
25 integer :: nc, i_c, c_id
28 do i_c = 1, af_num_children
29 c_id = boxes(id)%children(i_c)
30 if (c_id == af_no_box) cycle
31 call af_restrict_box(boxes(c_id), boxes(id), ivs)
33 end subroutine af_restrict_to_box
36 subroutine af_restrict_to_boxes(boxes, ids, ivs)
37 type(
box_t),
intent(inout) :: boxes(:)
38 integer,
intent(in) :: ids(:)
39 integer,
intent(in) :: ivs(:)
44 call af_restrict_to_box(boxes, ids(i), ivs)
47 end subroutine af_restrict_to_boxes
50 subroutine af_restrict_tree(tree, ivs)
51 type(
af_t),
intent(inout) :: tree
52 integer,
intent(in) :: ivs(:)
55 if (.not. tree%ready) stop
"Tree not ready"
56 do lvl = tree%highest_lvl-1, 1, -1
57 call af_restrict_to_boxes(tree%boxes, tree%lvls(lvl)%parents, ivs)
59 end subroutine af_restrict_tree
62 subroutine af_restrict_box(box_c, box_p, ivs, use_geometry)
63 type(
box_t),
intent(in) :: box_c
64 type(
box_t),
intent(inout) :: box_p
65 integer,
intent(in) :: ivs(:)
66 logical,
intent(in),
optional :: use_geometry
67 integer :: ijk, ijk_(c), ijk_(f), n, iv
68 integer :: hnc, ix_offset(ndim)
74 hnc = ishft(box_c%n_cell, -1)
75 ix_offset = af_get_child_offset(box_c)
77 if (
present(use_geometry))
then
78 use_geom = use_geometry
87 i_c = ix_offset(1) + i
89 box_p%cc(i_c, iv) = 0.5_dp * &
90 sum(box_c%cc(i_f:i_f+1, iv))
93 if (box_p%coord_t == af_cyl .and. use_geom)
then
95 j_c = ix_offset(2) + j
98 i_c = ix_offset(1) + i
101 call af_cyl_child_weights(box_p, i_c, w1, w2)
102 box_p%cc(i_c, j_c, iv) = 0.25_dp * (&
103 w1 * sum(box_c%cc(i_f, j_f:j_f+1, iv)) + &
104 w2 * sum(box_c%cc(i_f+1, j_f:j_f+1, iv)))
109 j_c = ix_offset(2) + j
112 i_c = ix_offset(1) + i
114 box_p%cc(i_c, j_c, iv) = 0.25_dp * &
115 sum(box_c%cc(i_f:i_f+1, j_f:j_f+1, iv))
121 k_c = ix_offset(3) + k
124 j_c = ix_offset(2) + j
127 i_c = ix_offset(1) + i
129 box_p%cc(i_c, j_c, k_c, iv) = 0.125_dp * &
130 sum(box_c%cc(i_f:i_f+1, j_f:j_f+1, k_f:k_f+1, iv))
136 end subroutine af_restrict_box
140 subroutine af_restrict_ref_boundary(tree, ivs)
141 type(
af_t),
intent(inout) :: tree
142 integer,
intent(in) :: ivs(:)
143 integer :: lvl, i, id, p_id
146 do lvl = 1, tree%highest_lvl
148 do i = 1,
size(tree%lvls(lvl)%leaves)
149 id = tree%lvls(lvl)%leaves(i)
150 p_id = tree%boxes(id)%parent
153 if (p_id > af_no_box .and. &
154 any(tree%boxes(id)%neighbors == af_no_box))
then
155 call af_restrict_box(tree%boxes(id), tree%boxes(p_id), ivs)
161 end subroutine af_restrict_ref_boundary
This module contains routines for restriction: going from fine to coarse variables.
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
The basic building block of afivo: a box with cell-centered and face centered data,...