Afivo 0.3
All Classes Namespaces Functions Variables Pages
m_af_restrict.f90
1!> This module contains routines for restriction: going from fine to coarse
2!> variables.
4#include "cpp_macros.h"
5 use m_af_types
6
7 implicit none
8 private
9
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
15 ! public :: af_restrict_box_face
16
17contains
18
19 !> Restrict the children of a box to the box (e.g., in 2D, average the values
20 !> at the four children to get the value for the parent)
21 subroutine af_restrict_to_box(boxes, id, ivs)
22 type(box_t), intent(inout) :: boxes(:) !< List of all the boxes
23 integer, intent(in) :: id !< Box whose children will be restricted to it
24 integer, intent(in) :: ivs(:) !< Variable to restrict
25 integer :: nc, i_c, c_id
26
27 nc = boxes(id)%n_cell
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)
32 end do
33 end subroutine af_restrict_to_box
34
35 !> Restrict the children of boxes ids(:) to them.
36 subroutine af_restrict_to_boxes(boxes, ids, ivs)
37 type(box_t), intent(inout) :: boxes(:) !< List of all the boxes
38 integer, intent(in) :: ids(:) !< Boxes whose children will be restricted to it
39 integer, intent(in) :: ivs(:) !< Variables to restrict
40 integer :: i
41
42 !$omp parallel do
43 do i = 1, size(ids)
44 call af_restrict_to_box(boxes, ids(i), ivs)
45 end do
46 !$omp end parallel do
47 end subroutine af_restrict_to_boxes
48
49 !> Restrict variables iv to all parent boxes, from the highest to the lowest level
50 subroutine af_restrict_tree(tree, ivs)
51 type(af_t), intent(inout) :: tree !< Tree to restrict on
52 integer, intent(in) :: ivs(:) !< Variables to restrict
53 integer :: lvl
54
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)
58 end do
59 end subroutine af_restrict_tree
60
61 !> Restriction of child box (box_c) to its parent (box_p)
62 subroutine af_restrict_box(box_c, box_p, ivs, use_geometry)
63 type(box_t), intent(in) :: box_c !< Child box to restrict
64 type(box_t), intent(inout) :: box_p !< Parent box to restrict to
65 integer, intent(in) :: ivs(:) !< Variable to restrict
66 logical, intent(in), optional :: use_geometry !< If set to false, don't use geometry
67 integer :: ijk, ijk_(c), ijk_(f), n, iv
68 integer :: hnc, ix_offset(ndim)
69 logical :: use_geom
70#if NDIM == 2
71 real(dp) :: w1, w2
72#endif
73
74 hnc = ishft(box_c%n_cell, -1) ! n_cell / 2
75 ix_offset = af_get_child_offset(box_c)
76
77 if (present(use_geometry)) then
78 use_geom = use_geometry
79 else
80 use_geom = .true.
81 end if
82
83 do n = 1, size(ivs)
84 iv = ivs(n)
85#if NDIM == 1
86 do i = 1, hnc
87 i_c = ix_offset(1) + i
88 i_f = 2 * i - 1
89 box_p%cc(i_c, iv) = 0.5_dp * &
90 sum(box_c%cc(i_f:i_f+1, iv))
91 end do
92#elif NDIM == 2
93 if (box_p%coord_t == af_cyl .and. use_geom) then
94 do j = 1, hnc
95 j_c = ix_offset(2) + j
96 j_f = 2 * j - 1
97 do i = 1, hnc
98 i_c = ix_offset(1) + i
99 i_f = 2 * i - 1
100
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)))
105 end do
106 end do
107 else
108 do j = 1, hnc
109 j_c = ix_offset(2) + j
110 j_f = 2 * j - 1
111 do i = 1, hnc
112 i_c = ix_offset(1) + i
113 i_f = 2 * i - 1
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))
116 end do
117 end do
118 endif
119#elif NDIM == 3
120 do k = 1, hnc
121 k_c = ix_offset(3) + k
122 k_f = 2 * k - 1
123 do j = 1, hnc
124 j_c = ix_offset(2) + j
125 j_f = 2 * j - 1
126 do i = 1, hnc
127 i_c = ix_offset(1) + i
128 i_f = 2 * i - 1
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))
131 end do
132 end do
133 end do
134#endif
135 end do
136 end subroutine af_restrict_box
137
138 !> Restrict only next to refinement boundaries, which which can be required
139 !> for filling coarse-grid ghost cells
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
144
145 !$omp parallel private(lvl, i, id, p_id)
146 do lvl = 1, tree%highest_lvl
147 !$omp do
148 do i = 1, size(tree%lvls(lvl)%leaves)
149 id = tree%lvls(lvl)%leaves(i)
150 p_id = tree%boxes(id)%parent
151
152 ! Only restrict near refinement boundaries
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)
156 end if
157 end do
158 !$omp end do
159 end do
160 !$omp end parallel
161 end subroutine af_restrict_ref_boundary
162
163! !> Restriction of face-centered variables from child to parent
164! subroutine af_restrict_box_face(box_c, box_p, ivf, ivf_to)
165! type(box_t), intent(in) :: box_c !< Child box to restrict
166! type(box_t), intent(inout) :: box_p !< Parent box to restrict to
167! integer, intent(in) :: ivf !< Face-variable to restrict
168! integer, intent(in), optional :: ivf_to !< Destination (if /= ivf)
169! integer :: i, j, i_f, j_f, i_c, j_c, i_dest
170! integer :: hnc, ix_offset(NDIM)
171! #if NDIM == 3
172! integer :: k, k_f, k_c
173! #endif
174
175! hnc = ishft(box_c%n_cell, -1) ! n_cell / 2
176! ix_offset = af_get_child_offset(box_c)
177
178! if (present(ivf_to)) then
179! i_dest = ivf_to
180! else
181! i_dest = ivf
182! end if
183
184! if (box_p%coord_t == af_cyl) &
185! stop "restrict_box_face not implemented for cylindrical case"
186
187! #if NDIM == 2
188! do j = 1, hnc
189! j_c = ix_offset(2) + j
190! j_f = 2 * j - 1
191! do i = 1, hnc
192! i_c = ix_offset(1) + i
193! i_f = 2 * i - 1
194
195! box_p%fc(i_c, j_c, 1, i_dest) = 0.5_dp * &
196! sum(box_c%fc(i_f, j_f:j_f+1, 1, ivf))
197! if (i == hnc) then
198! box_p%fc(i_c+1, j_c, 1, i_dest) = 0.5_dp * &
199! sum(box_c%fc(i_f+2, j_f:j_f+1, 1, ivf))
200! end if
201
202! box_p%fc(i_c, j_c, 2, i_dest) = 0.5_dp * &
203! sum(box_c%fc(i_f:i_f+1, j_f, 2, ivf))
204! if (j == hnc) then
205! box_p%fc(i_c, j_c+1, 2, i_dest) = 0.5_dp * &
206! sum(box_c%fc(i_f:i_f+1, j_f+2, 2, ivf))
207! end if
208
209! end do
210! end do
211! #elif NDIM == 3
212! if (box_p%coord_t == af_cyl) &
213! stop "restrict_box_face not implemented for 3D case"
214! #endif
215! end subroutine af_restrict_box_face
216
217end module m_af_restrict
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...
Definition: m_af_types.f90:3
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
Definition: m_af_types.f90:326
The basic building block of afivo: a box with cell-centered and face centered data,...
Definition: m_af_types.f90:286