Afivo  0.3
m_a3_prolong.f90
1 
5  use m_a3_types
6 
7  implicit none
8  private
9 
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
15  ! public :: a3_prolong_quadratic_from
16  ! public :: a3_prolong_quadratic
17 
18  public :: a3_prolong_limit
19  public :: a3_prolong_linear_cons
20 
21 contains
22 
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
30  integer :: i_c, c_id
31 
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)
36  end do
37  end subroutine a3_prolong_copy_from
38 
40  subroutine a3_prolong_copy(box_p, box_c, iv, iv_to, low, high, add)
41  type(box3_t), intent(in) :: box_p
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
48  logical :: add_to
49  integer :: nc, ix_offset(3), ivc
50  integer :: i, j, i_c1, j_c1, lo(3), hi(3)
51  integer :: k, k_c1
52 
53  nc = box_c%n_cell
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
58 
59  ! Offset of child w.r.t. parent
60  ix_offset = a3_get_child_offset(box_c)
61 
62  if (add_to) then
63  do k = lo(3), hi(3)
64  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
65  do j = lo(2), hi(2)
66  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
67  do i = lo(1), hi(1)
68  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
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)
71  end do
72  end do
73  end do
74  else
75  do k = lo(3), hi(3)
76  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
77  do j = lo(2), hi(2)
78  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
79  do i = lo(1), hi(1)
80  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
81  box_c%cc(i, j, k, ivc) = box_p%cc(i_c1, j_c1, k_c1, iv)
82  end do
83  end do
84  end do
85  end if
86  end subroutine a3_prolong_copy
87 
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
96  integer :: i_c, c_id
97 
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)
102  end do
103  end subroutine a3_prolong_linear_from
104 
108  subroutine a3_prolong_sparse(box_p, box_c, iv, iv_to, add)
109  type(box3_t), intent(in) :: box_p
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
117  logical :: add_to
118  real(dp) :: flz, fhz
119  integer :: k, k_c, k_f
120 
121  nc = box_c%n_cell
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
126 
127  if (.not. add_to) then
128  box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
129  end if
130 
131  do k = 1, hnc
132  k_c = k + ix_offset(3)
133  k_f = 2 * k - 1
134  do j = 1, hnc
135  j_c = j + ix_offset(2)
136  j_f = 2 * j - 1
137  do i = 1, hnc
138  i_c = i + ix_offset(1)
139  i_f = 2 * i - 1
140 
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)
148 
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)
165  end do
166  end do
167  end do
168  end subroutine a3_prolong_sparse
169 
170  ! Prolong with a limited slope in the coarse cells, which takes the minimum of
171  ! the left/right slopes if they have the same sign (else it is zero). This
172  ! procedure is conservative.
173  subroutine a3_prolong_limit(box_p, box_c, iv, iv_to, add)
174  type(box3_t), intent(in) :: box_p
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
182  logical :: add_to
183  real(dp) :: fzL, fzR, fz
184  integer :: k, k_c, k_f
185 
186  nc = box_c%n_cell
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
191 
192  if (.not. add_to) then
193  box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
194  end if
195 
196  do k = 1, hnc
197  k_c = k + ix_offset(3)
198  k_f = 2 * k - 1
199  do j = 1, hnc
200  j_c = j + ix_offset(2)
201  j_f = 2 * j - 1
202  do i = 1, hnc
203  i_c = i + ix_offset(1)
204  i_f = 2 * i - 1
205 
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
213 
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)
217 
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)
234  end do
235  end do
236  end do
237  contains
238 
239  ! Take minimum of two slopes if they have the same sign, else take zero
240  elemental function limit_slope(ll, rr) result(slope)
241  real(dp), intent(in) :: ll, rr
242  real(dp) :: slope
243 
244  if (ll * rr < 0) then
245  slope = 0.0_dp
246  else if (ll * ll < rr * rr) then
247  slope = ll
248  else
249  slope = rr
250  end if
251  end function limit_slope
252 
253  end subroutine a3_prolong_limit
254 
255  ! Prolong with a linear (unlimited) slope in the coarse cells, which can
256  ! result in negative densities. This procedure is conservative.
257  subroutine a3_prolong_linear_cons(box_p, box_c, iv, iv_to, add)
258  type(box3_t), intent(in) :: box_p
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
266  logical :: add_to
267  real(dp) :: fz
268  integer :: k, k_c, k_f
269 
270  nc = box_c%n_cell
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
275 
276  if (.not. add_to) then
277  box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
278  end if
279 
280  do k = 1, hnc
281  k_c = k + ix_offset(3)
282  k_f = 2 * k - 1
283  do j = 1, hnc
284  j_c = j + ix_offset(2)
285  j_f = 2 * j - 1
286  do i = 1, hnc
287  i_c = i + ix_offset(1)
288  i_f = 2 * i - 1
289 
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))
297 
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)
314  end do
315  end do
316  end do
317 
318  end subroutine a3_prolong_linear_cons
319 
321  subroutine a3_prolong_linear(box_p, box_c, iv, iv_to, add)
322  type(box3_t), intent(in) :: box_p
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
329  logical :: add_to
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
337 
338  nc = box_c%n_cell
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
343 
344  if (.not. add_to) then
345  box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
346  end if
347 
348  do k = 1, hnc
349  k_c = k + ix_offset(3)
350  k_f = 2 * k - 1
351  do j = 1, hnc
352  j_c = j + ix_offset(2)
353  j_f = 2 * j - 1
354  do i = 1, hnc
355  i_c = i + ix_offset(1)
356  i_f = 2 * i - 1
357 
358  f000 = f27 * box_p%cc(i_c, j_c, k_c, iv)
359 
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)
367 
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)
375 
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)
382 
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)
386 
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)
390 
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)
415  end do
416  end do
417  end do
418  end subroutine a3_prolong_linear
419 
420  ! !> Quadratic prolongation to children. We use stencils that do not require
421  ! !> corner ghost cells.
422  ! subroutine a3_prolong_quadratic_from(boxes, id, iv, iv_to, add)
423  ! type(box3_t), intent(inout) :: boxes(:) !< List of all boxes
424  ! integer, intent(in) :: id !< Box whose children we will fill
425  ! integer, intent(in) :: iv !< Variable that is filled
426  ! integer, intent(in), optional :: iv_to !< Destination variable
427  ! logical, intent(in), optional :: add !< Add to old values
428  ! integer :: i_c, c_id
429 
430  ! do i_c = 1, a3_num_children
431  ! c_id = boxes(id)%children(i_c)
432  ! if (c_id == af_no_box) cycle
433  ! call a3_prolong_quadratic(boxes(id), boxes(c_id), iv, iv_to, add)
434  ! end do
435  ! end subroutine a3_prolong_quadratic_from
436 
437 ! !> Prolongation to a child (from parent) using quadratic interpolation (using
438 ! !> a local Taylor approximation)
439 ! !> \todo Mixed derivatives in 3D version
440 ! subroutine a3_prolong_quadratic(box_p, box_c, iv, iv_to, add)
441 ! type(box3_t), intent(in) :: box_p !< Parent box
442 ! type(box3_t), intent(inout) :: box_c !< Child box
443 ! integer, intent(in) :: iv !< Variable to fill
444 ! integer, intent(in), optional :: iv_to !< Destination variable
445 ! logical, intent(in), optional :: add !< Add to old values
446 ! logical :: add_to
447 ! integer :: hnc, ix_offset(3)
448 ! integer :: i, j, nc, ivc
449 ! integer :: i_c, i_f, j_c, j_f
450 ! real(dp) :: f0, fx, fy, fxx, fyy, f2
451 ! #if 3 == 2
452 ! real(dp) :: fxy(2**3)
453 ! #elif 3 == 3
454 ! real(dp) :: fz, fzz
455 ! integer :: k, k_c, k_f
456 ! #endif
457 
458 ! nc = box_c%n_cell
459 ! hnc = ishft(box_c%n_cell, -1)
460 ! ix_offset = a3_get_child_offset(box_c)
461 ! add_to = .false.; if (present(add)) add_to = add
462 ! ivc = iv; if (present(iv_to)) ivc = iv_to
463 
464 ! if (.not. add_to) then
465 ! #if 3 == 2
466 ! box_c%cc(1:nc, 1:nc, ivc) = 0
467 ! #elif 3 == 3
468 ! box_c%cc(1:nc, 1:nc, 1:nc, ivc) = 0
469 ! #endif
470 ! end if
471 
472 ! #if 3 == 2
473 ! do j = 1, hnc
474 ! j_c = j + ix_offset(2)
475 ! j_f = 2 * j - 1
476 ! do i = 1, hnc
477 ! i_c = i + ix_offset(1)
478 ! i_f = 2 * i - 1
479 
480 ! f0 = box_p%cc(i_c, j_c, iv)
481 ! fx = 0.125_dp * (box_p%cc(i_c+1, j_c, iv) - &
482 ! box_p%cc(i_c-1, j_c, iv))
483 ! fy = 0.125_dp * (box_p%cc(i_c, j_c+1, iv) - &
484 ! box_p%cc(i_c, j_c-1, iv))
485 ! fxx = 0.03125_dp * (box_p%cc(i_c-1, j_c, iv) - &
486 ! 2 * f0 + box_p%cc(i_c+1, j_c, iv))
487 ! fyy = 0.03125_dp * (box_p%cc(i_c, j_c-1, iv) - &
488 ! 2 * f0 + box_p%cc(i_c, j_c+1, iv))
489 ! f2 = fxx + fyy
490 
491 ! fxy(1) = 0.0625_dp * (box_p%cc(i_c-1, j_c-1, iv) + f0 - &
492 ! box_p%cc(i_c-1, j_c, iv) - box_p%cc(i_c, j_c-1, iv))
493 ! fxy(2) = 0.0625_dp * (box_p%cc(i_c+1, j_c-1, iv) + f0 - &
494 ! box_p%cc(i_c+1, j_c, iv) - box_p%cc(i_c, j_c-1, iv))
495 ! fxy(3) = 0.0625_dp * (box_p%cc(i_c-1, j_c+1, iv) + f0 - &
496 ! box_p%cc(i_c-1, j_c, iv) - box_p%cc(i_c, j_c+1, iv))
497 ! fxy(4) = 0.0625_dp * (box_p%cc(i_c+1, j_c+1, iv) + f0 - &
498 ! box_p%cc(i_c+1, j_c, iv) - box_p%cc(i_c, j_c+1, iv))
499 
500 ! box_c%cc(i_f, j_f, ivc) = f0 - fx - fy + f2 + fxy(1) + &
501 ! box_c%cc(i_f, j_f, ivc)
502 ! box_c%cc(i_f+1, j_f, ivc) = f0 + fx - fy + f2 + fxy(2) + &
503 ! box_c%cc(i_f+1, j_f, ivc)
504 ! box_c%cc(i_f, j_f+1, ivc) = f0 - fx + fy + f2 + fxy(3) + &
505 ! box_c%cc(i_f, j_f+1, ivc)
506 ! box_c%cc(i_f+1, j_f+1, ivc) = f0 + fx + fy + f2 + fxy(4) + &
507 ! box_c%cc(i_f+1, j_f+1, ivc)
508 ! end do
509 ! end do
510 ! #elif 3 == 3
511 ! do k = 1, hnc
512 ! k_c = k + ix_offset(3)
513 ! k_f = 2 * k - 1
514 ! do j = 1, hnc
515 ! j_c = j + ix_offset(2)
516 ! j_f = 2 * j - 1
517 ! do i = 1, hnc
518 ! i_c = i + ix_offset(1)
519 ! i_f = 2 * i - 1
520 
521 ! f0 = box_p%cc(i_c, j_c, k_c, iv)
522 ! fx = 0.125_dp * (box_p%cc(i_c+1, j_c, k_c, iv) - &
523 ! box_p%cc(i_c-1, j_c, k_c, iv))
524 ! fy = 0.125_dp * (box_p%cc(i_c, j_c+1, k_c, iv) - &
525 ! box_p%cc(i_c, j_c-1, k_c, iv))
526 ! fz = 0.125_dp * (box_p%cc(i_c, j_c, k_c+1, iv) - &
527 ! box_p%cc(i_c, j_c, k_c-1, iv))
528 ! fxx = 0.03125_dp * (box_p%cc(i_c-1, j_c, k_c, iv) - &
529 ! 2 * f0 + box_p%cc(i_c+1, j_c, k_c, iv))
530 ! fyy = 0.03125_dp * (box_p%cc(i_c, j_c-1, k_c, iv) - &
531 ! 2 * f0 + box_p%cc(i_c, j_c+1, k_c, iv))
532 ! fzz = 0.03125_dp * (box_p%cc(i_c, j_c, k_c-1, iv) - &
533 ! 2 * f0 + box_p%cc(i_c, j_c, k_c+1, iv))
534 ! f2 = fxx + fyy + fzz
535 
536 ! box_c%cc(i_f, j_f, k_f, ivc) = f0 - fx - fy - fz + f2 + &
537 ! box_c%cc(i_f, j_f, k_f, ivc)
538 ! box_c%cc(i_f+1, j_f, k_f, ivc) = f0 + fx - fy - fz + f2 + &
539 ! box_c%cc(i_f+1, j_f, k_f, ivc)
540 ! box_c%cc(i_f, j_f+1, k_f, ivc) = f0 - fx + fy - fz + f2 + &
541 ! box_c%cc(i_f, j_f+1, k_f, ivc)
542 ! box_c%cc(i_f+1, j_f+1, k_f, ivc) = f0 + fx + fy - fz + f2 + &
543 ! box_c%cc(i_f+1, j_f+1, k_f, ivc)
544 ! box_c%cc(i_f, j_f, k_f+1, ivc) = f0 - fx - fy + fz + f2 + &
545 ! box_c%cc(i_f, j_f, k_f+1, ivc)
546 ! box_c%cc(i_f+1, j_f, k_f+1, ivc) = f0 + fx - fy + fz + f2 + &
547 ! box_c%cc(i_f+1, j_f, k_f+1, ivc)
548 ! box_c%cc(i_f, j_f+1, k_f+1, ivc) = f0 - fx + fy + fz + f2 + &
549 ! box_c%cc(i_f, j_f+1, k_f+1, ivc)
550 ! box_c%cc(i_f+1, j_f+1, k_f+1, ivc) = f0 + fx + fy + fz + f2 + &
551 ! box_c%cc(i_f+1, j_f+1, k_f+1, ivc)
552 ! end do
553 ! end do
554 ! end do
555 ! #endif
556 ! end subroutine a3_prolong_quadratic
557 
558 end module m_a3_prolong
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
This module contains the routines related to prolongation: going from coarse to fine variables...
Definition: m_a3_prolong.f90:4
This module contains the basic types and constants that are used in the 3-dimensional version of Afiv...
Definition: m_a3_types.f90:5