10 public :: a3_prolong_copy_from
11 public :: a3_prolong_copy
12 public :: a3_prolong_linear_from
13 public :: a3_prolong_sparse
14 public :: a3_prolong_linear
18 public :: a3_prolong_limit
19 public :: a3_prolong_linear_cons
24 subroutine a3_prolong_copy_from(boxes, id, iv, iv_to, add)
25 type(
box3_t),
intent(inout) :: boxes(:)
26 integer,
intent(in) :: id
27 integer,
intent(in) :: iv
28 integer,
intent(in),
optional :: iv_to
29 logical,
intent(in),
optional :: add
32 do i_c = 1, a3_num_children
33 c_id = boxes(id)%children(i_c)
34 if (c_id == af_no_box) cycle
35 call a3_prolong_copy(boxes(id), boxes(c_id), iv, iv_to=iv_to, add=add)
37 end subroutine a3_prolong_copy_from
40 subroutine a3_prolong_copy(box_p, box_c, iv, iv_to, low, high, add)
42 type(
box3_t),
intent(inout) :: box_c
43 integer,
intent(in) :: iv
44 integer,
intent(in),
optional :: iv_to
45 integer,
intent(in),
optional :: low(3)
46 integer,
intent(in),
optional :: high(3)
47 logical,
intent(in),
optional :: add
49 integer :: nc, ix_offset(3), ivc
50 integer :: i, j, i_c1, j_c1, lo(3), hi(3)
54 add_to = .false.;
if (
present(add)) add_to = add
55 ivc = iv;
if (
present(iv_to)) ivc = iv_to
56 lo = 1;
if (
present(low)) lo = low
57 hi = nc;
if (
present(high)) hi = high
60 ix_offset = a3_get_child_offset(box_c)
64 k_c1 = ix_offset(3) + ishft(k+1, -1)
66 j_c1 = ix_offset(2) + ishft(j+1, -1)
68 i_c1 = ix_offset(1) + ishft(i+1, -1)
69 box_c%cc(i, j, k, ivc) = box_c%cc(i, j, k, ivc) + &
70 box_p%cc(i_c1, j_c1, k_c1, iv)
76 k_c1 = ix_offset(3) + ishft(k+1, -1)
78 j_c1 = ix_offset(2) + ishft(j+1, -1)
80 i_c1 = ix_offset(1) + ishft(i+1, -1)
81 box_c%cc(i, j, k, ivc) = box_p%cc(i_c1, j_c1, k_c1, iv)
86 end subroutine a3_prolong_copy
90 subroutine a3_prolong_linear_from(boxes, id, iv, iv_to, add)
91 type(
box3_t),
intent(inout) :: boxes(:)
92 integer,
intent(in) :: id
93 integer,
intent(in) :: iv
94 integer,
intent(in),
optional :: iv_to
95 logical,
intent(in),
optional :: add
98 do i_c = 1, a3_num_children
99 c_id = boxes(id)%children(i_c)
100 if (c_id == af_no_box) cycle
101 call a3_prolong_linear(boxes(id), boxes(c_id), iv, iv_to, add)
103 end subroutine a3_prolong_linear_from
108 subroutine a3_prolong_sparse(box_p, box_c, iv, iv_to, add)
110 type(
box3_t),
intent(inout) :: box_c
111 integer,
intent(in) :: iv
112 integer,
intent(in),
optional :: iv_to
113 logical,
intent(in),
optional :: add
114 integer :: hnc, nc, ix_offset(3), ivc
115 integer :: i, j, i_c, i_f, j_c, j_f
116 real(dp) :: f0, flx, fhx, fly, fhy
119 integer :: k, k_c, k_f
122 hnc = ishft(box_c%n_cell, -1)
123 ix_offset = a3_get_child_offset(box_c)
124 add_to = .false.;
if (
present(add)) add_to = add
125 ivc = iv;
if (
present(iv_to)) ivc = iv_to
127 if (.not. add_to)
then 128 box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
132 k_c = k + ix_offset(3)
135 j_c = j + ix_offset(2)
138 i_c = i + ix_offset(1)
141 f0 = 0.25_dp * box_p%cc(i_c, j_c, k_c, iv)
142 flx = 0.25_dp * box_p%cc(i_c-1, j_c, k_c, iv)
143 fhx = 0.25_dp * box_p%cc(i_c+1, j_c, k_c, iv)
144 fly = 0.25_dp * box_p%cc(i_c, j_c-1, k_c, iv)
145 fhy = 0.25_dp * box_p%cc(i_c, j_c+1, k_c, iv)
146 flz = 0.25_dp * box_p%cc(i_c, j_c, k_c-1, iv)
147 fhz = 0.25_dp * box_p%cc(i_c, j_c, k_c+1, iv)
149 box_c%cc(i_f, j_f, k_f, ivc) = f0 + flx + &
150 fly + flz + box_c%cc(i_f, j_f, k_f, ivc)
151 box_c%cc(i_f+1, j_f, k_f, ivc) = f0 + fhx + &
152 fly + flz + box_c%cc(i_f+1, j_f, k_f, ivc)
153 box_c%cc(i_f, j_f+1, k_f, ivc) = f0 + flx + &
154 fhy + flz + box_c%cc(i_f, j_f+1, k_f, ivc)
155 box_c%cc(i_f+1, j_f+1, k_f, ivc) = f0 + fhx + &
156 fhy + flz + box_c%cc(i_f+1, j_f+1, k_f, ivc)
157 box_c%cc(i_f, j_f, k_f+1, ivc) = f0 + flx + &
158 fly + fhz + box_c%cc(i_f, j_f, k_f+1, ivc)
159 box_c%cc(i_f+1, j_f, k_f+1, ivc) = f0 + fhx + &
160 fly + fhz + box_c%cc(i_f+1, j_f, k_f+1, ivc)
161 box_c%cc(i_f, j_f+1, k_f+1, ivc) = f0 + flx + &
162 fhy + fhz + box_c%cc(i_f, j_f+1, k_f+1, ivc)
163 box_c%cc(i_f+1, j_f+1, k_f+1, ivc) = f0 + fhx + &
164 fhy + fhz + box_c%cc(i_f+1, j_f+1, k_f+1, ivc)
168 end subroutine a3_prolong_sparse
173 subroutine a3_prolong_limit(box_p, box_c, iv, iv_to, add)
175 type(
box3_t),
intent(inout) :: box_c
176 integer,
intent(in) :: iv
177 integer,
intent(in),
optional :: iv_to
178 logical,
intent(in),
optional :: add
179 integer :: hnc, nc, ix_offset(3), ivc
180 integer :: i, j, i_c, i_f, j_c, j_f
181 real(dp) :: f0, fxL, fxR, fyL, fyR, fx, fy
183 real(dp) :: fzL, fzR, fz
184 integer :: k, k_c, k_f
187 hnc = ishft(box_c%n_cell, -1)
188 ix_offset = a3_get_child_offset(box_c)
189 add_to = .false.;
if (
present(add)) add_to = add
190 ivc = iv;
if (
present(iv_to)) ivc = iv_to
192 if (.not. add_to)
then 193 box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
197 k_c = k + ix_offset(3)
200 j_c = j + ix_offset(2)
203 i_c = i + ix_offset(1)
206 f0 = box_p%cc(i_c, j_c, k_c, iv)
207 fxl = f0 - box_p%cc(i_c-1, j_c, k_c, iv)
208 fxr = box_p%cc(i_c+1, j_c, k_c, iv) - f0
209 fyl = f0 - box_p%cc(i_c, j_c-1, k_c, iv)
210 fyr = box_p%cc(i_c, j_c+1, k_c, iv) - f0
211 fzl = f0 - box_p%cc(i_c, j_c, k_c-1, iv)
212 fzr = box_p%cc(i_c, j_c, k_c+1, iv) - f0
214 fx = 0.25_dp * limit_slope(fxl, fxr)
215 fy = 0.25_dp * limit_slope(fyl, fyr)
216 fz = 0.25_dp * limit_slope(fzl, fzr)
218 box_c%cc(i_f, j_f, k_f, ivc) = f0 - fx - &
219 fy - fz + box_c%cc(i_f, j_f, k_f, ivc)
220 box_c%cc(i_f+1, j_f, k_f, ivc) = f0 + fx - &
221 fy - fz + box_c%cc(i_f+1, j_f, k_f, ivc)
222 box_c%cc(i_f, j_f+1, k_f, ivc) = f0 - fx + &
223 fy - fz + box_c%cc(i_f, j_f+1, k_f, ivc)
224 box_c%cc(i_f+1, j_f+1, k_f, ivc) = f0 + fx + &
225 fy - fz + box_c%cc(i_f+1, j_f+1, k_f, ivc)
226 box_c%cc(i_f, j_f, k_f+1, ivc) = f0 - fx - &
227 fy + fz + box_c%cc(i_f, j_f, k_f+1, ivc)
228 box_c%cc(i_f+1, j_f, k_f+1, ivc) = f0 + fx - &
229 fy + fz + box_c%cc(i_f+1, j_f, k_f+1, ivc)
230 box_c%cc(i_f, j_f+1, k_f+1, ivc) = f0 - fx + &
231 fy + fz + box_c%cc(i_f, j_f+1, k_f+1, ivc)
232 box_c%cc(i_f+1, j_f+1, k_f+1, ivc) = f0 + fx + &
233 fy + fz + box_c%cc(i_f+1, j_f+1, k_f+1, ivc)
240 elemental function limit_slope(ll, rr)
result(slope)
241 real(dp),
intent(in) :: ll, rr
244 if (ll * rr < 0)
then 246 else if (ll * ll < rr * rr)
then 251 end function limit_slope
253 end subroutine a3_prolong_limit
257 subroutine a3_prolong_linear_cons(box_p, box_c, iv, iv_to, add)
259 type(
box3_t),
intent(inout) :: box_c
260 integer,
intent(in) :: iv
261 integer,
intent(in),
optional :: iv_to
262 logical,
intent(in),
optional :: add
263 integer :: hnc, nc, ix_offset(3), ivc
264 integer :: i, j, i_c, i_f, j_c, j_f
265 real(dp) :: f0, fx, fy
268 integer :: k, k_c, k_f
271 hnc = ishft(box_c%n_cell, -1)
272 ix_offset = a3_get_child_offset(box_c)
273 add_to = .false.;
if (
present(add)) add_to = add
274 ivc = iv;
if (
present(iv_to)) ivc = iv_to
276 if (.not. add_to)
then 277 box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
281 k_c = k + ix_offset(3)
284 j_c = j + ix_offset(2)
287 i_c = i + ix_offset(1)
290 f0 = box_p%cc(i_c, j_c, k_c, iv)
291 fx = 0.125_dp * (box_p%cc(i_c+1, j_c, k_c, iv) - &
292 box_p%cc(i_c-1, j_c, k_c, iv))
293 fy = 0.125_dp * (box_p%cc(i_c, j_c+1, k_c, iv) - &
294 box_p%cc(i_c, j_c-1, k_c, iv))
295 fz = 0.125_dp * (box_p%cc(i_c, j_c, k_c+1, iv) - &
296 box_p%cc(i_c, j_c, k_c-1, iv))
298 box_c%cc(i_f, j_f, k_f, ivc) = f0 - fx - &
299 fy - fz + box_c%cc(i_f, j_f, k_f, ivc)
300 box_c%cc(i_f+1, j_f, k_f, ivc) = f0 + fx - &
301 fy - fz + box_c%cc(i_f+1, j_f, k_f, ivc)
302 box_c%cc(i_f, j_f+1, k_f, ivc) = f0 - fx + &
303 fy - fz + box_c%cc(i_f, j_f+1, k_f, ivc)
304 box_c%cc(i_f+1, j_f+1, k_f, ivc) = f0 + fx + &
305 fy - fz + box_c%cc(i_f+1, j_f+1, k_f, ivc)
306 box_c%cc(i_f, j_f, k_f+1, ivc) = f0 - fx - &
307 fy + fz + box_c%cc(i_f, j_f, k_f+1, ivc)
308 box_c%cc(i_f+1, j_f, k_f+1, ivc) = f0 + fx - &
309 fy + fz + box_c%cc(i_f+1, j_f, k_f+1, ivc)
310 box_c%cc(i_f, j_f+1, k_f+1, ivc) = f0 - fx + &
311 fy + fz + box_c%cc(i_f, j_f+1, k_f+1, ivc)
312 box_c%cc(i_f+1, j_f+1, k_f+1, ivc) = f0 + fx + &
313 fy + fz + box_c%cc(i_f+1, j_f+1, k_f+1, ivc)
318 end subroutine a3_prolong_linear_cons
321 subroutine a3_prolong_linear(box_p, box_c, iv, iv_to, add)
323 type(
box3_t),
intent(inout) :: box_c
324 integer,
intent(in) :: iv
325 integer,
intent(in),
optional :: iv_to
326 logical,
intent(in),
optional :: add
327 integer :: hnc, nc, ix_offset(3), ivc
328 integer :: i, j, i_c, i_f, j_c, j_f
330 real(dp) :: f000, f00l, f0l0, f0ll, fl00, fl0l, fll0
331 real(dp) :: flll, f00h, f0h0, f0hh, fh00, fh0h, fhh0
332 real(dp) :: fhhh, f0lh, f0hl, fl0h, fh0l, flh0, fhl0
333 real(dp) :: fllh, flhl, fhll, fhhl, fhlh, flhh
334 real(dp),
parameter :: f1 = 1/64.0_dp, f3=3/64.0_dp, f9=9/64.0_dp
335 real(dp),
parameter :: f27 = 27/64.0_dp
336 integer :: k, k_c, k_f
339 hnc = ishft(box_c%n_cell, -1)
340 ix_offset = a3_get_child_offset(box_c)
341 add_to = .false.;
if (
present(add)) add_to = add
342 ivc = iv;
if (
present(iv_to)) ivc = iv_to
344 if (.not. add_to)
then 345 box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
349 k_c = k + ix_offset(3)
352 j_c = j + ix_offset(2)
355 i_c = i + ix_offset(1)
358 f000 = f27 * box_p%cc(i_c, j_c, k_c, iv)
360 f00l = f9 * box_p%cc(i_c, j_c, k_c-1, iv)
361 f0l0 = f9 * box_p%cc(i_c, j_c-1, k_c, iv)
362 f0ll = f3 * box_p%cc(i_c, j_c-1, k_c-1, iv)
363 fl00 = f9 * box_p%cc(i_c-1, j_c, k_c, iv)
364 fl0l = f3 * box_p%cc(i_c-1, j_c, k_c-1, iv)
365 fll0 = f3 * box_p%cc(i_c-1, j_c-1, k_c, iv)
366 flll = f1 * box_p%cc(i_c-1, j_c-1, k_c-1, iv)
368 f00h = f9 * box_p%cc(i_c, j_c, k_c+1, iv)
369 f0h0 = f9 * box_p%cc(i_c, j_c+1, k_c, iv)
370 f0hh = f3 * box_p%cc(i_c, j_c+1, k_c+1, iv)
371 fh00 = f9 * box_p%cc(i_c+1, j_c, k_c, iv)
372 fh0h = f3 * box_p%cc(i_c+1, j_c, k_c+1, iv)
373 fhh0 = f3 * box_p%cc(i_c+1, j_c+1, k_c, iv)
374 fhhh = f1 * box_p%cc(i_c+1, j_c+1, k_c+1, iv)
376 fl0h = f3 * box_p%cc(i_c-1, j_c, k_c+1, iv)
377 fh0l = f3 * box_p%cc(i_c+1, j_c, k_c-1, iv)
378 flh0 = f3 * box_p%cc(i_c-1, j_c+1, k_c, iv)
379 fhl0 = f3 * box_p%cc(i_c+1, j_c-1, k_c, iv)
380 f0lh = f3 * box_p%cc(i_c, j_c-1, k_c+1, iv)
381 f0hl = f3 * box_p%cc(i_c, j_c+1, k_c-1, iv)
383 fllh = f1 * box_p%cc(i_c-1, j_c-1, k_c+1, iv)
384 flhl = f1 * box_p%cc(i_c-1, j_c+1, k_c-1, iv)
385 fhll = f1 * box_p%cc(i_c+1, j_c-1, k_c-1, iv)
387 fhhl = f1 * box_p%cc(i_c+1, j_c+1, k_c-1, iv)
388 fhlh = f1 * box_p%cc(i_c+1, j_c-1, k_c+1, iv)
389 flhh = f1 * box_p%cc(i_c-1, j_c+1, k_c+1, iv)
391 box_c%cc(i_f, j_f, k_f, ivc) = f000 + fl00 + &
392 f0l0 + f00l + fll0 + fl0l + f0ll + flll + &
393 box_c%cc(i_f, j_f, k_f, ivc)
394 box_c%cc(i_f+1, j_f, k_f, ivc) = f000 + fh00 + &
395 f0l0 + f00l + fhl0 + fh0l + f0ll + fhll + &
396 box_c%cc(i_f+1, j_f, k_f, ivc)
397 box_c%cc(i_f, j_f+1, k_f, ivc) = f000 + fl00 + &
398 f0h0 + f00l + flh0 + fl0l + f0hl + flhl + &
399 box_c%cc(i_f, j_f+1, k_f, ivc)
400 box_c%cc(i_f+1, j_f+1, k_f, ivc) = f000 + fh00 + &
401 f0h0 + f00l + fhh0 + fh0l + f0hl + fhhl + &
402 box_c%cc(i_f+1, j_f+1, k_f, ivc)
403 box_c%cc(i_f, j_f, k_f+1, ivc) = f000 + fl00 + &
404 f0l0 + f00h + fll0 + fl0h + f0lh + fllh + &
405 box_c%cc(i_f, j_f, k_f+1, ivc)
406 box_c%cc(i_f+1, j_f, k_f+1, ivc) = f000 + fh00 + &
407 f0l0 + f00h + fhl0 + fh0h + f0lh + fhlh + &
408 box_c%cc(i_f+1, j_f, k_f+1, ivc)
409 box_c%cc(i_f, j_f+1, k_f+1, ivc) = f000 + fl00 + &
410 f0h0 + f00h + flh0 + fl0h + f0hh + flhh + &
411 box_c%cc(i_f, j_f+1, k_f+1, ivc)
412 box_c%cc(i_f+1, j_f+1, k_f+1, ivc) = f000 + fh00 + &
413 f0h0 + f00h + fhh0 + fh0h + f0hh + fhhh + &
414 box_c%cc(i_f+1, j_f+1, k_f+1, ivc)
418 end subroutine a3_prolong_linear
The basic building block of afivo: a box with cell-centered and face centered data, and information about its position, neighbors, children etc.
This module contains the routines related to prolongation: going from coarse to fine variables...
This module contains the basic types and constants that are used in the 3-dimensional version of Afiv...