Afivo  0.3
m_af_stencil.f90
1 #include "cpp_macros.h"
2 
4  use m_af_types
5 
6  implicit none
7  private
8 
10  integer, parameter :: stencil_list_initial_size = 10
11 
12  integer, parameter, public :: stencil_constant = 1
13  integer, parameter, public :: stencil_variable = 2
14  integer, parameter, public :: stencil_sparse = 3
15 
17  integer, parameter, private :: num_shapes = 5
18 
20  integer, parameter, public :: af_stencil_357 = 1
21 
23  integer, parameter, public :: af_stencil_p234 = 2
24 
26  integer, parameter, public :: af_stencil_p248 = 3
27 
29  integer, parameter, public :: af_stencil_246 = 4
30 
32  integer, parameter, public :: af_stencil_mask = 5
33 
35  integer, parameter, public :: af_stencil_sizes(num_shapes) = &
36  [2*ndim+1, ndim+1, 2**ndim, 2*ndim, 1]
37 
38  abstract interface
39 
40  subroutine af_subr_stencil(box, stencil)
41  import
42  type(box_t), intent(in) :: box
43  type(stencil_t), intent(inout) :: stencil
44  end subroutine af_subr_stencil
45  end interface
46 
47  public :: af_stencil_print_info
48  public :: af_stencil_index
49  public :: af_stencil_prepare_store
50  public :: af_stencil_remove
51  public :: af_stencil_store_box
52  public :: af_stencil_check_box
53  public :: af_stencil_allocate_coeff
54  public :: af_stencil_store
55  public :: af_stencil_apply
56  public :: af_stencil_apply_box
57  public :: af_stencil_gsrb_box
58  public :: af_stencil_prolong_box
59  public :: af_stencil_get_box
60  public :: af_stencil_try_constant
61 
62 contains
63 
65  subroutine af_stencil_print_info(tree)
66  type(af_t), intent(in) :: tree
67  integer :: id, n_stencils_stored
68  integer :: n_boxes_with_stencils
69  integer :: n_constant_stored, n_variable_stored
70  integer :: n, n_stencils, n_sparse_stored
71 
72  n_stencils_stored = 0
73  n_constant_stored = 0
74  n_variable_stored = 0
75  n_sparse_stored = 0
76  n_boxes_with_stencils = 0
77 
78  do id = 1, tree%highest_id
79  if (.not. tree%boxes(id)%in_use) cycle
80 
81  n_stencils = tree%boxes(id)%n_stencils
82  if (n_stencils > 0) then
83  n_boxes_with_stencils = n_boxes_with_stencils + 1
84  n_stencils_stored = n_stencils_stored + n_stencils
85 
86  do n = 1, n_stencils
87  select case (tree%boxes(id)%stencils(n)%stype)
88  case (stencil_constant)
89  n_constant_stored = n_constant_stored + 1
90  case (stencil_variable)
91  n_variable_stored = n_variable_stored + 1
92  case (stencil_sparse)
93  n_sparse_stored = n_sparse_stored + 1
94  end select
95  end do
96  end if
97  end do
98 
99  write(*, '(A)') ''
100  write(*, '(A)') ' ** Information about stencils **'
101  write(*, '(A, I0)') ' #boxes with stencils: ', n_boxes_with_stencils
102  write(*, '(A, I0)') ' #stencils stored: ', n_stencils_stored
103  write(*, '(A, I0)') ' #constant stencils: ', n_constant_stored
104  write(*, '(A, I0)') ' #variable stencils: ', n_variable_stored
105  write(*, '(A, I0)') ' #sparse stencils: ', n_sparse_stored
106  end subroutine af_stencil_print_info
107 
109  pure integer function af_stencil_index(box, key)
110  type(box_t), intent(in) :: box
111  integer, intent(in) :: key
112  integer :: i
113 
114  af_stencil_index = af_stencil_none
115  do i = 1, box%n_stencils
116  if (box%stencils(i)%key == key) then
117  af_stencil_index = i
118  exit
119  end if
120  end do
121  end function af_stencil_index
122 
124  subroutine af_stencil_store(tree, key, set_stencil)
125  type(af_t), intent(inout) :: tree
126  integer, intent(in) :: key
128  procedure(af_subr_stencil) :: set_stencil
129  integer :: lvl, i, id
130 
131  !$omp parallel private(lvl, i, id)
132  do lvl = 1, tree%highest_lvl
133  !$omp do
134  do i = 1, size(tree%lvls(lvl)%ids)
135  id = tree%lvls(lvl)%ids(i)
136  call af_stencil_store_box(tree%boxes(id), key, set_stencil)
137  end do
138  !$omp end do
139  end do
140  !$omp end parallel
141  end subroutine af_stencil_store
142 
144  subroutine af_stencil_prepare_store(box, key, ix)
145  type(box_t), intent(inout) :: box
146  integer, intent(in) :: key
147  integer, intent(out) :: ix
148  type(stencil_t), allocatable :: tmp(:)
149 
150  ix = af_stencil_index(box, key)
151 
152  if (ix == af_stencil_none) then
153  ! Allocate storage if necessary
154  if (.not. allocated(box%stencils)) then
155  allocate(box%stencils(stencil_list_initial_size))
156  else if (box%n_stencils == size(box%stencils)) then
157  allocate(tmp(2 * box%n_stencils))
158  tmp(1:box%n_stencils) = box%stencils
159  call move_alloc(tmp, box%stencils)
160  end if
161 
162  box%n_stencils = box%n_stencils + 1
163  ix = box%n_stencils
164  box%stencils(ix)%key = key
165  else
166  error stop "Stencil with this key has already been stored"
167  end if
168  end subroutine af_stencil_prepare_store
169 
171  subroutine af_stencil_allocate_coeff(stencil, nc, use_f, n_sparse)
172  type(stencil_t), intent(inout) :: stencil
174  integer, intent(in) :: nc
177  logical, intent(in), optional :: use_f
179  integer, intent(in), optional :: n_sparse
180  logical :: allocate_f
181  integer :: n_coeff
182 
183  allocate_f = .false.; if (present(use_f)) allocate_f = use_f
184  if (stencil%shape < 1 .or. stencil%shape > num_shapes) &
185  error stop "Unknown stencil shape"
186 
187  n_coeff = af_stencil_sizes(stencil%shape)
188 
189  select case (stencil%stype)
190  case (stencil_constant)
191  if (.not. allocated(stencil%c)) then
192  allocate(stencil%c(n_coeff))
193  else if (size(stencil%c) /= n_coeff) then
194  deallocate(stencil%c)
195  allocate(stencil%c(n_coeff))
196  end if
197  case (stencil_variable)
198  if (.not. allocated(stencil%v)) then
199  allocate(stencil%v(n_coeff, dtimes(nc)))
200  else if (size(stencil%v, 1) /= n_coeff) then
201  deallocate(stencil%v)
202  allocate(stencil%v(n_coeff, dtimes(nc)))
203  end if
204 
205  if (allocate_f .and. .not. allocated(stencil%f)) then
206  allocate(stencil%f(dtimes(nc)))
207  allocate(stencil%bc_correction(dtimes(nc)))
208  end if
209  case (stencil_sparse)
210  if (.not. present(n_sparse)) error stop "n_sparse required"
211 
212  if (.not. allocated(stencil%sparse_v)) then
213  allocate(stencil%sparse_ix(ndim, n_sparse))
214  allocate(stencil%sparse_v(n_coeff, n_sparse))
215  else if (any(size(stencil%sparse_v) /= [n_coeff, n_sparse])) then
216  deallocate(stencil%sparse_v)
217  deallocate(stencil%sparse_ix)
218  allocate(stencil%sparse_ix(ndim, n_sparse))
219  allocate(stencil%sparse_v(n_coeff, n_sparse))
220  end if
221  case default
222  error stop "Unknow stencil%stype"
223  end select
224  end subroutine af_stencil_allocate_coeff
225 
227  subroutine af_stencil_remove(tree, key)
228  type(af_t), intent(inout) :: tree
229  integer, intent(in) :: key
230  integer :: lvl, i, id
231 
232  !$omp parallel private(lvl, i, id)
233  do lvl = 1, tree%highest_lvl
234  !$omp do
235  do i = 1, size(tree%lvls(lvl)%ids)
236  id = tree%lvls(lvl)%ids(i)
237  call af_stencil_remove_box(tree%boxes(id), key)
238  end do
239  !$omp end do
240  end do
241  !$omp end parallel
242  end subroutine af_stencil_remove
243 
245  subroutine af_stencil_store_box(box, key, set_stencil)
246  type(box_t), intent(inout) :: box
247  integer, intent(in) :: key
249  procedure(af_subr_stencil) :: set_stencil
250  integer :: ix
251 
252  call af_stencil_prepare_store(box, key, ix)
253  call set_stencil(box, box%stencils(ix))
254  end subroutine af_stencil_store_box
255 
257  subroutine af_stencil_remove_box(box, key)
258  type(box_t), intent(inout) :: box
259  integer, intent(in) :: key
260  integer :: ix
261  type(stencil_t) :: empty_stencil
262 
263  ix = af_stencil_index(box, key)
264  if (ix == af_stencil_none) error stop "Stencil not present"
265 
266  ! Ensure there are no gaps in the stencil list
267  if (ix == box%n_stencils) then
268  box%stencils(ix) = empty_stencil
269  else
270  box%stencils(ix) = box%stencils(box%n_stencils)
271  box%stencils(box%n_stencils) = empty_stencil
272  end if
273 
274  box%n_stencils = box%n_stencils - 1
275  end subroutine af_stencil_remove_box
276 
278  subroutine af_stencil_check_box(box, key, ix)
279  type(box_t), intent(in) :: box
280  integer, intent(in) :: key
281  integer, intent(in) :: ix
282 
283  if (key == af_stencil_none) &
284  error stop "key == af_stencil_none"
285 
286  if (af_stencil_index(box, key) /= ix) &
287  error stop "Stencil with key not found at index"
288 
289  associate(stencil => box%stencils(ix))
290  if (stencil%shape < 1 .or. stencil%shape > num_shapes) &
291  error stop "Unknown stencil shape"
292  select case (stencil%stype)
293  case (stencil_constant)
294  if (.not. allocated(stencil%c)) &
295  error stop "stencil%c not allocated"
296  case (stencil_variable)
297  if (.not. allocated(stencil%v)) &
298  error stop "stencil%v not allocated"
299  case (stencil_sparse)
300  if (.not. allocated(stencil%sparse_ix)) &
301  error stop "stencil%sparse_ix not allocated"
302  if (.not. allocated(stencil%sparse_v)) &
303  error stop "stencil%sparse_v not allocated"
304  case default
305  error stop "Unknow stencil%stype"
306  end select
307  end associate
308  end subroutine af_stencil_check_box
309 
310  subroutine af_stencil_apply(tree, key, iv, i_out)
311  type(af_t), intent(inout) :: tree
312  integer, intent(in) :: key
313  integer, intent(in) :: iv
314  integer, intent(in) :: i_out
315 
316  integer :: lvl, i, id
317 
318  !$omp parallel private(lvl, i, id)
319  do lvl = 1, tree%highest_lvl
320  !$omp do
321  do i = 1, size(tree%lvls(lvl)%ids)
322  id = tree%lvls(lvl)%ids(i)
323  call af_stencil_apply_box(tree%boxes(id), key, iv, i_out)
324  end do
325  !$omp end do
326  end do
327  !$omp end parallel
328  end subroutine af_stencil_apply
329 
330  subroutine af_stencil_apply_box(box, key, iv, i_out)
331  type(box_t), intent(inout) :: box
332  integer, intent(in) :: key
333  integer, intent(in) :: iv
334  integer, intent(in) :: i_out
335  integer :: i
336 
337  i = af_stencil_index(box, key)
338  if (i == af_stencil_none) error stop "Stencil not stored"
339 
340  select case (box%stencils(i)%shape)
341  case (af_stencil_357)
342  call stencil_apply_357(box, box%stencils(i), iv, i_out)
343  case default
344  error stop "Unknown stencil"
345  end select
346  end subroutine af_stencil_apply_box
347 
349  subroutine stencil_apply_357(box, stencil, iv, i_out)
350  type(box_t), intent(inout) :: box
351  type(stencil_t), intent(in) :: stencil
352  integer, intent(in) :: iv
353  integer, intent(in) :: i_out
354  real(dp) :: c(2*NDIM+1)
355  integer :: IJK
356 #if NDIM == 2
357  real(dp) :: rfac(2, box%n_cell), c_cyl(2*NDIM+1)
358  real(dp) :: cc_cyl(2*NDIM+1, box%n_cell)
359 #endif
360 
361  if (iv == i_out) error stop "Cannot have iv == i_out"
362  if (stencil%stype == stencil_sparse) error stop "sparse not implemented"
363 
364  associate(cc => box%cc, nc => box%n_cell)
365 
366 #if NDIM == 1
367  if (stencil%stype == stencil_constant) then
368  c = stencil%c
369  do kji_do(1, nc)
370  cc(i, i_out) = &
371  c(1) * cc(i, iv) + c(2) * cc(i-1, iv) + c(3) * cc(i+1, iv)
372  end do; close_do
373  else
374  do kji_do(1, nc)
375  c = stencil%v(:, ijk)
376  cc(i, i_out) = &
377  c(1) * cc(i, iv) + c(2) * cc(i-1, iv) + c(3) * cc(i+1, iv)
378  end do; close_do
379  end if
380 #elif NDIM == 2
381  if (stencil%cylindrical_gradient) then
382  ! Correct for cylindrical coordinates, assuming the elements correspond
383  ! to a gradient operation
384  call af_cyl_flux_factors(box, rfac)
385 
386  if (stencil%stype == stencil_constant) then
387  c = stencil%c
388 
389  ! Pre-compute coefficients for each i-index
390  do i = 1, nc
391  cc_cyl(2:3, i) = rfac(1:2, i) * c(2:3)
392  cc_cyl(1, i) = c(1) - (cc_cyl(2, i) - c(2)) &
393  - (cc_cyl(3, i) - c(3))
394  cc_cyl(4:, i) = c(4:)
395  end do
396 
397  do kji_do(1, nc)
398  cc(i, j, i_out) = &
399  cc_cyl(1, i) * cc(i, j, iv) + &
400  cc_cyl(2, i) * cc(i-1, j, iv) + &
401  cc_cyl(3, i) * cc(i+1, j, iv) + &
402  cc_cyl(4, i) * cc(i, j-1, iv) + &
403  cc_cyl(5, i) * cc(i, j+1, iv)
404  end do; close_do
405  else
406  ! Variable stencil
407  do kji_do(1, nc)
408  c = stencil%v(:, ijk)
409  c_cyl(2:3) = rfac(1:2, i) * c(2:3)
410  c_cyl(1) = c(1) - (c_cyl(2) - c(2)) - (c_cyl(3) - c(3))
411  c_cyl(4:) = c(4:)
412  cc(i, j, i_out) = &
413  c_cyl(1) * cc(i, j, iv) + &
414  c_cyl(2) * cc(i-1, j, iv) + &
415  c_cyl(3) * cc(i+1, j, iv) + &
416  c_cyl(4) * cc(i, j-1, iv) + &
417  c_cyl(5) * cc(i, j+1, iv)
418  end do; close_do
419  end if
420  else
421  ! No cylindrical gradient correction
422  if (stencil%stype == stencil_constant) then
423  c = stencil%c
424  do kji_do(1, nc)
425  cc(i, j, i_out) = &
426  c(1) * cc(i, j, iv) + &
427  c(2) * cc(i-1, j, iv) + &
428  c(3) * cc(i+1, j, iv) + &
429  c(4) * cc(i, j-1, iv) + &
430  c(5) * cc(i, j+1, iv)
431  end do; close_do
432  else
433  do kji_do(1, nc)
434  c = stencil%v(:, ijk)
435  cc(i, j, i_out) = &
436  c(1) * cc(i, j, iv) + &
437  c(2) * cc(i-1, j, iv) + &
438  c(3) * cc(i+1, j, iv) + &
439  c(4) * cc(i, j-1, iv) + &
440  c(5) * cc(i, j+1, iv)
441  end do; close_do
442  end if
443  end if
444 #elif NDIM == 3
445  if (stencil%stype == stencil_constant) then
446  c = stencil%c
447  do kji_do(1, nc)
448  cc(i, j, k, i_out) = &
449  c(1) * cc(i, j, k, iv) + &
450  c(2) * cc(i-1, j, k, iv) + &
451  c(3) * cc(i+1, j, k, iv) + &
452  c(4) * cc(i, j-1, k, iv) + &
453  c(5) * cc(i, j+1, k, iv) + &
454  c(6) * cc(i, j, k-1, iv) + &
455  c(7) * cc(i, j, k+1, iv)
456  end do; close_do
457  else
458  do kji_do(1, nc)
459  c = stencil%v(:, ijk)
460  cc(i, j, k, i_out) = &
461  c(1) * cc(i, j, k, iv) + &
462  c(2) * cc(i-1, j, k, iv) + &
463  c(3) * cc(i+1, j, k, iv) + &
464  c(4) * cc(i, j-1, k, iv) + &
465  c(5) * cc(i, j+1, k, iv) + &
466  c(6) * cc(i, j, k-1, iv) + &
467  c(7) * cc(i, j, k+1, iv)
468  end do; close_do
469  end if
470 #endif
471 
472  if (allocated(stencil%bc_correction)) then
473  cc(dtimes(1:nc), i_out) = cc(dtimes(1:nc), i_out) - &
474  stencil%bc_correction
475  end if
476  end associate
477 
478  end subroutine stencil_apply_357
479 
480  ! subroutine stencil_correct_bc_357(box, stencil, iv, tmp)
481  ! type(box_t), intent(inout) :: box !< Operate on this box
482  ! type(stencil_t), intent(in) :: stencil !< Stencil
483  ! integer, intent(in) :: iv !< Input variable
484  ! real(dp), intent(inout) :: tmp(DTIMES(box%n_cell))
485  ! integer :: bc_ix, bc_type, nb, lo(NDIM), hi(NDIM)
486  ! integer :: nb_dim, olo(NDIM), ohi(NDIM)
487  ! real(dp) :: bc_val(box%n_cell**(NDIM-1))
488  ! real(dp) :: stencil_coeffs(box%n_cell**(NDIM-1))
489  ! real(dp) :: values_inside(box%n_cell**(NDIM-1))
490  ! real(dp) :: values_outside(box%n_cell**(NDIM-1))
491 
492  ! do bc_ix = 1, box%n_bc
493  ! nb = box%bc_index_to_nb(bc_ix)
494 
495  ! bc_val = box%bc_val(:, iv, bc_ix)
496  ! bc_type = box%bc_type(iv, bc_ix)
497 
498  ! ! Determine index range next to boundary
499  ! call af_get_index_bc_inside(nb, box%n_cell, 1, lo, hi)
500  ! call af_get_index_bc_outside(nb, box%n_cell, 1, olo, ohi)
501 
502  ! ! Get stencil coefficients near boundary
503  ! if (stencil%constant) then
504  ! stencil_coeffs = stencil%c(nb+1)
505  ! else
506  ! stencil_coeffs = pack(stencil%v(nb+1, &
507  ! DSLICE(lo, hi)), .true.)
508  ! end if
509 
510  ! values_inside = pack(box%cc(DSLICE(lo, hi), iv), .true.)
511  ! values_outside = pack(box%cc(DSLICE(olo, ohi), iv), .true.)
512 
513  ! select case (bc_type)
514  ! case (af_bc_dirichlet)
515  ! ! Dirichlet value at cell face, so compute gradient over h/2
516  ! ! E.g. 1 -2 1 becomes 0 -3 1 for a 1D Laplacian
517  ! tmp(DSLICE(lo, hi)) = tmp(DSLICE(lo, hi)) &
518  ! + reshape(2 * stencil_coeffs * bc_val - &
519  ! stencil_coeffs * (values_inside + values_outside), &
520  ! hi - lo + 1)
521  ! case (af_bc_neumann)
522  ! nb_dim = af_neighb_dim(nb)
523  ! ! E.g. 1 -2 1 becomes 0 -1 1 for a 1D Laplacian
524  ! tmp(DSLICE(lo, hi)) = tmp(DSLICE(lo, hi)) &
525  ! - reshape(stencil_coeffs * (values_outside - values_inside) &
526  ! - stencil_coeffs * af_neighb_high_pm(nb) * box%dr(nb_dim) * &
527  ! bc_val, hi - lo + 1)
528  ! case default
529  ! error stop "unsupported boundary condition"
530  ! end select
531  ! end do
532  ! end subroutine stencil_correct_bc_357
533 
534  subroutine af_stencil_prolong_box(box_p, box_c, key, iv, iv_to, add)
535  type(box_t), intent(in) :: box_p
536  type(box_t), intent(inout) :: box_c
537  integer, intent(in) :: key
538  integer, intent(in) :: iv
539  integer, intent(in) :: iv_to
540  logical, intent(in), optional :: add
541  logical :: add_to
542  integer :: i
543 
544  add_to = .false.
545  if (present(add)) add_to = add
546 
547  i = af_stencil_index(box_c, key)
548  if (i == af_stencil_none) error stop "Stencil not stored"
549 
550  select case (box_c%stencils(i)%shape)
551  case (af_stencil_p234)
552  call stencil_prolong_234(box_p, box_c, box_c%stencils(i), &
553  iv, iv_to, add_to)
554  case (af_stencil_p248)
555  call stencil_prolong_248(box_p, box_c, box_c%stencils(i), &
556  iv, iv_to, add_to)
557  case default
558  print *, "box_c%stencils(i)%shape: ", box_c%stencils(i)%shape
559  error stop "Unknown stencil"
560  end select
561  end subroutine af_stencil_prolong_box
562 
564  subroutine stencil_prolong_234(box_p, box_c, stencil, iv, iv_to, add)
565  type(box_t), intent(in) :: box_p
566  type(box_t), intent(inout) :: box_c
567  type(stencil_t), intent(in) :: stencil
568  integer, intent(in) :: iv
569  integer, intent(in) :: iv_to
570  logical, intent(in) :: add
571 
572  integer :: nc, ix_offset(NDIM), IJK
573  integer :: IJK_(c1), IJK_(c2)
574  real(dp) :: c(NDIM+1)
575 
576  nc = box_c%n_cell
577  ix_offset = af_get_child_offset(box_c)
578  if (stencil%stype == stencil_sparse) error stop "sparse not implemented"
579  if (.not. add) box_c%cc(dtimes(1:nc), iv_to) = 0
580 
581  ! In these loops, we calculate the closest coarse index (_c1), and the
582  ! one-but-closest (_c2). The fine cell lies in between.
583 #if NDIM == 1
584  if (stencil%stype == stencil_constant) then
585  c = stencil%c
586  do i = 1, nc
587  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
588  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
589  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
590  c(1) * box_p%cc(i_c1, iv) + &
591  c(2) * box_p%cc(i_c2, iv)
592  end do
593  else
594  do i = 1, nc
595  c = stencil%v(:, ijk)
596  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
597  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
598  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
599  c(1) * box_p%cc(i_c1, iv) + &
600  c(2) * box_p%cc(i_c2, iv)
601  end do
602  end if
603 #elif NDIM == 2
604  if (stencil%stype == stencil_constant) then
605  c = stencil%c
606  do j = 1, nc
607  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
608  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
609  do i = 1, nc
610  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
611  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
612  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
613  c(1) * box_p%cc(i_c1, j_c1, iv) + &
614  c(2) * box_p%cc(i_c2, j_c1, iv) + &
615  c(3) * box_p%cc(i_c1, j_c2, iv)
616  end do
617  end do
618  else
619  do j = 1, nc
620  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
621  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
622  do i = 1, nc
623  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
624  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
625  c = stencil%v(:, ijk)
626  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
627  c(1) * box_p%cc(i_c1, j_c1, iv) + &
628  c(2) * box_p%cc(i_c2, j_c1, iv) + &
629  c(3) * box_p%cc(i_c1, j_c2, iv)
630  end do
631  end do
632  end if
633 #elif NDIM == 3
634  if (stencil%stype == stencil_constant) then
635  c = stencil%c
636  do k = 1, nc
637  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
638  k_c2 = k_c1 + 1 - 2 * iand(k, 1) ! even: +1, odd: -1
639  do j = 1, nc
640  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
641  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
642  do i = 1, nc
643  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
644  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
645  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
646  c(1) * box_p%cc(i_c1, j_c1, k_c1, iv) + &
647  c(2) * box_p%cc(i_c2, j_c1, k_c1, iv) + &
648  c(3) * box_p%cc(i_c1, j_c2, k_c1, iv) + &
649  c(4) * box_p%cc(i_c1, j_c1, k_c2, iv)
650  end do
651  end do
652  end do
653  else
654  do k = 1, nc
655  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
656  k_c2 = k_c1 + 1 - 2 * iand(k, 1) ! even: +1, odd: -1
657  do j = 1, nc
658  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
659  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
660  do i = 1, nc
661  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
662  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
663  c = stencil%v(:, ijk)
664  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
665  c(1) * box_p%cc(i_c1, j_c1, k_c1, iv) + &
666  c(2) * box_p%cc(i_c2, j_c1, k_c1, iv) + &
667  c(3) * box_p%cc(i_c1, j_c2, k_c1, iv) + &
668  c(4) * box_p%cc(i_c1, j_c1, k_c2, iv)
669  end do
670  end do
671  end do
672  end if
673 #endif
674  end subroutine stencil_prolong_234
675 
677  subroutine stencil_prolong_248(box_p, box_c, stencil, iv, iv_to, add)
678  type(box_t), intent(in) :: box_p
679  type(box_t), intent(inout) :: box_c
680  type(stencil_t), intent(in) :: stencil
681  integer, intent(in) :: iv
682  integer, intent(in) :: iv_to
683  logical, intent(in) :: add
684 
685  integer :: nc, ix_offset(NDIM), IJK
686  integer :: IJK_(c1), IJK_(c2)
687  real(dp) :: c(2**NDIM)
688 
689  nc = box_c%n_cell
690  ix_offset = af_get_child_offset(box_c)
691  if (stencil%stype == stencil_sparse) error stop "sparse not implemented"
692  if (.not. add) box_c%cc(dtimes(1:nc), iv_to) = 0
693 
694  ! In these loops, we calculate the closest coarse index (_c1), and the
695  ! one-but-closest (_c2). The fine cell lies in between.
696 #if NDIM == 1
697  if (stencil%stype == stencil_constant) then
698  c = stencil%c
699  do i = 1, nc
700  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
701  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
702  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
703  c(1) * box_p%cc(i_c1, iv) + &
704  c(2) * box_p%cc(i_c2, iv)
705  end do
706  else
707  do i = 1, nc
708  c = stencil%v(:, ijk)
709  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
710  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
711  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
712  c(1) * box_p%cc(i_c1, iv) + &
713  c(2) * box_p%cc(i_c2, iv)
714  end do
715  end if
716 #elif NDIM == 2
717  if (stencil%stype == stencil_constant) then
718  c = stencil%c
719  do j = 1, nc
720  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
721  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
722  do i = 1, nc
723  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
724  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
725  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
726  c(1) * box_p%cc(i_c1, j_c1, iv) + &
727  c(2) * box_p%cc(i_c2, j_c1, iv) + &
728  c(3) * box_p%cc(i_c1, j_c2, iv) + &
729  c(4) * box_p%cc(i_c2, j_c2, iv)
730  end do
731  end do
732  else
733  do j = 1, nc
734  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
735  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
736  do i = 1, nc
737  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
738  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
739  c = stencil%v(:, ijk)
740  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
741  c(1) * box_p%cc(i_c1, j_c1, iv) + &
742  c(2) * box_p%cc(i_c2, j_c1, iv) + &
743  c(3) * box_p%cc(i_c1, j_c2, iv) + &
744  c(4) * box_p%cc(i_c2, j_c2, iv)
745  end do
746  end do
747  end if
748 #elif NDIM == 3
749  if (stencil%stype == stencil_constant) then
750  c = stencil%c
751  do k = 1, nc
752  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
753  k_c2 = k_c1 + 1 - 2 * iand(k, 1) ! even: +1, odd: -1
754  do j = 1, nc
755  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
756  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
757  do i = 1, nc
758  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
759  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
760  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
761  c(1) * box_p%cc(i_c1, j_c1, k_c1, iv) + &
762  c(2) * box_p%cc(i_c2, j_c1, k_c1, iv) + &
763  c(3) * box_p%cc(i_c1, j_c2, k_c1, iv) + &
764  c(4) * box_p%cc(i_c2, j_c2, k_c1, iv) + &
765  c(5) * box_p%cc(i_c1, j_c1, k_c2, iv) + &
766  c(6) * box_p%cc(i_c2, j_c1, k_c2, iv) + &
767  c(7) * box_p%cc(i_c1, j_c2, k_c2, iv) + &
768  c(8) * box_p%cc(i_c2, j_c2, k_c2, iv)
769  end do
770  end do
771  end do
772  else
773  do k = 1, nc
774  k_c1 = ix_offset(3) + ishft(k+1, -1) ! (k+1)/2
775  k_c2 = k_c1 + 1 - 2 * iand(k, 1) ! even: +1, odd: -1
776  do j = 1, nc
777  j_c1 = ix_offset(2) + ishft(j+1, -1) ! (j+1)/2
778  j_c2 = j_c1 + 1 - 2 * iand(j, 1) ! even: +1, odd: -1
779  do i = 1, nc
780  i_c1 = ix_offset(1) + ishft(i+1, -1) ! (i+1)/2
781  i_c2 = i_c1 + 1 - 2 * iand(i, 1) ! even: +1, odd: -1
782  c = stencil%v(:, ijk)
783  box_c%cc(ijk, iv_to) = box_c%cc(ijk, iv_to) + &
784  c(1) * box_p%cc(i_c1, j_c1, k_c1, iv) + &
785  c(2) * box_p%cc(i_c2, j_c1, k_c1, iv) + &
786  c(3) * box_p%cc(i_c1, j_c2, k_c1, iv) + &
787  c(4) * box_p%cc(i_c2, j_c2, k_c1, iv) + &
788  c(5) * box_p%cc(i_c1, j_c1, k_c2, iv) + &
789  c(6) * box_p%cc(i_c2, j_c1, k_c2, iv) + &
790  c(7) * box_p%cc(i_c1, j_c2, k_c2, iv) + &
791  c(8) * box_p%cc(i_c2, j_c2, k_c2, iv)
792  end do
793  end do
794  end do
795  end if
796 #endif
797  end subroutine stencil_prolong_248
798 
800  subroutine af_stencil_gsrb_box(box, key, redblack, iv, i_rhs)
801  type(box_t), intent(inout) :: box
802  integer, intent(in) :: key
803  integer, intent(in) :: redblack
804  integer, intent(in) :: iv
805  integer, intent(in) :: i_rhs
806  integer :: i
807 
808  i = af_stencil_index(box, key)
809  if (i == af_stencil_none) error stop "Stencil not stored"
810 
811  select case (box%stencils(i)%shape)
812  case (af_stencil_357)
813  call stencil_gsrb_357(box, box%stencils(i), redblack, iv, i_rhs)
814  case default
815  error stop "Unknown stencil"
816  end select
817  end subroutine af_stencil_gsrb_box
818 
820  subroutine stencil_gsrb_357(box, stencil, redblack, iv, i_rhs)
821  type(box_t), intent(inout) :: box
822  type(stencil_t), intent(in) :: stencil
823  integer, intent(in) :: redblack
824  integer, intent(in) :: iv
825  integer, intent(in) :: i_rhs
826 
827  real(dp) :: c(2*NDIM+1), inv_c1
828  integer :: IJK, i0, nc
829 #if NDIM == 2
830  real(dp) :: rfac(2, box%n_cell), c_cyl(2*NDIM+1)
831  real(dp) :: cc_cyl(2*NDIM+1, box%n_cell), inv_cc1(box%n_cell)
832 #endif
833 
834  if (stencil%stype == stencil_sparse) error stop "sparse not implemented"
835  nc = box%n_cell
836 
837  associate(cc => box%cc, nc => box%n_cell)
838  if (allocated(stencil%bc_correction)) then
839  cc(dtimes(1:nc), i_rhs) = cc(dtimes(1:nc), i_rhs) + &
840  stencil%bc_correction
841  end if
842 
843 #if NDIM == 1
844  i0 = 2 - iand(redblack, 1)
845  if (stencil%stype == stencil_constant) then
846  c = stencil%c
847  inv_c1 = 1 / c(1)
848 
849  do i = i0, nc, 2
850  cc(ijk, iv) = (cc(ijk, i_rhs) &
851  - c(2) * cc(i-1, iv) &
852  - c(3) * cc(i+1, iv)) * inv_c1
853  end do
854  else
855  do i = i0, nc, 2
856  c = stencil%v(:, ijk)
857  cc(ijk, iv) = (cc(ijk, i_rhs) &
858  - c(2) * cc(i-1, iv) &
859  - c(3) * cc(i+1, iv)) / c(1)
860  end do
861  end if
862 #elif NDIM == 2
863  if (stencil%cylindrical_gradient) then
864  ! Correct for cylindrical coordinates, assuming the elements correspond
865  ! to a gradient operation
866  call af_cyl_flux_factors(box, rfac)
867 
868  if (stencil%stype == stencil_constant) then
869  c = stencil%c
870 
871  ! Pre-compute coefficients for each i-index
872  do i = 1, nc
873  cc_cyl(2:3, i) = rfac(1:2, i) * c(2:3)
874  cc_cyl(1, i) = c(1) - (cc_cyl(2, i) - c(2)) &
875  - (cc_cyl(3, i) - c(3))
876  cc_cyl(4:, i) = c(4:)
877  inv_cc1(i) = 1 / cc_cyl(1, i)
878  end do
879 
880  do j = 1, nc
881  i0 = 2 - iand(ieor(redblack, j), 1)
882  do i = i0, nc, 2
883  cc(ijk, iv) = (cc(ijk, i_rhs) &
884  - cc_cyl(2, i) * cc(i-1, j, iv) &
885  - cc_cyl(3, i) * cc(i+1, j, iv) &
886  - cc_cyl(4, i) * cc(i, j-1, iv) &
887  - cc_cyl(5, i) * cc(i, j+1, iv)) * inv_cc1(i)
888  end do
889  end do
890  else
891  ! Variable stencil
892  do j = 1, nc
893  i0 = 2 - iand(ieor(redblack, j), 1)
894  do i = i0, nc, 2
895  c = stencil%v(:, ijk)
896  c_cyl(2:3) = rfac(1:2, i) * c(2:3)
897  c_cyl(1) = c(1) - (c_cyl(2) - c(2)) - (c_cyl(3) - c(3))
898  c_cyl(4:) = c(4:)
899 
900  cc(ijk, iv) = (cc(ijk, i_rhs) &
901  - c_cyl(2) * cc(i-1, j, iv) &
902  - c_cyl(3) * cc(i+1, j, iv) &
903  - c_cyl(4) * cc(i, j-1, iv) &
904  - c_cyl(5) * cc(i, j+1, iv)) / c_cyl(1)
905  end do
906  end do
907  end if
908  else ! No cylindrical gradient correction
909  if (stencil%stype == stencil_constant) then
910  c = stencil%c
911  inv_c1 = 1 / c(1)
912 
913  do j = 1, nc
914  i0 = 2 - iand(ieor(redblack, j), 1)
915  do i = i0, nc, 2
916  cc(ijk, iv) = (cc(ijk, i_rhs) &
917  - c(2) * cc(i-1, j, iv) &
918  - c(3) * cc(i+1, j, iv) &
919  - c(4) * cc(i, j-1, iv) &
920  - c(5) * cc(i, j+1, iv)) * inv_c1
921  end do
922  end do
923  else
924  do j = 1, nc
925  i0 = 2 - iand(ieor(redblack, j), 1)
926  do i = i0, nc, 2
927  c = stencil%v(:, ijk)
928  cc(ijk, iv) = (cc(ijk, i_rhs) &
929  - c(2) * cc(i-1, j, iv) &
930  - c(3) * cc(i+1, j, iv) &
931  - c(4) * cc(i, j-1, iv) &
932  - c(5) * cc(i, j+1, iv)) / c(1)
933  end do
934  end do
935  end if
936  end if
937 #elif NDIM == 3
938  if (stencil%stype == stencil_constant) then
939  c = stencil%c
940  inv_c1 = 1 / c(1)
941 
942  do k = 1, nc
943  do j = 1, nc
944  i0 = 2 - iand(ieor(redblack, k+j), 1)
945  do i = i0, nc, 2
946  cc(ijk, iv) = (cc(ijk, i_rhs) &
947  - c(2) * cc(i-1, j, k, iv) &
948  - c(3) * cc(i+1, j, k, iv) &
949  - c(4) * cc(i, j-1, k, iv) &
950  - c(5) * cc(i, j+1, k, iv) &
951  - c(6) * cc(i, j, k-1, iv) &
952  - c(7) * cc(i, j, k+1, iv)) * inv_c1
953  end do
954  end do
955  end do
956  else
957  do k = 1, nc
958  do j = 1, nc
959  i0 = 2 - iand(ieor(redblack, k+j), 1)
960  do i = i0, nc, 2
961  c = stencil%v(:, ijk)
962  cc(ijk, iv) = (cc(ijk, i_rhs) &
963  - c(2) * cc(i-1, j, k, iv) &
964  - c(3) * cc(i+1, j, k, iv) &
965  - c(4) * cc(i, j-1, k, iv) &
966  - c(5) * cc(i, j+1, k, iv) &
967  - c(6) * cc(i, j, k-1, iv) &
968  - c(7) * cc(i, j, k+1, iv)) / c(1)
969  end do
970  end do
971  end do
972  end if
973 #endif
974 
975  if (allocated(stencil%bc_correction)) then
976  cc(dtimes(1:nc), i_rhs) = cc(dtimes(1:nc), i_rhs) - &
977  stencil%bc_correction
978  end if
979  end associate
980  end subroutine stencil_gsrb_357
981 
983  subroutine af_stencil_try_constant(box, ix, abs_tol, success)
984  type(box_t), intent(inout) :: box
985  integer, intent(in) :: ix
986  real(dp), intent(in) :: abs_tol
987  logical :: success
988  integer :: ijk, nc, n_coeff
989 
990  if (box%stencils(ix)%stype == stencil_sparse) &
991  error stop "sparse not implemented"
992  if (.not. allocated(box%stencils(ix)%v)) &
993  error stop "No variable stencil present"
994 
995  nc = box%n_cell
996  success = .false.
997 
998  ! Check if all coefficients are the same, otherwise return
999  do kji_do(1, nc)
1000  if (any(abs(box%stencils(ix)%v(:, ijk) - &
1001  box%stencils(ix)%v(:, dtimes(1))) > abs_tol)) return
1002  end do; close_do
1003 
1004  box%stencils(ix)%stype = stencil_constant
1005  n_coeff = size(box%stencils(ix)%v, 1)
1006  allocate(box%stencils(ix)%c(n_coeff))
1007  box%stencils(ix)%c = box%stencils(ix)%v(:, dtimes(1))
1008  deallocate(box%stencils(ix)%v)
1009  success = .true.
1010  end subroutine af_stencil_try_constant
1011 
1013  subroutine af_stencil_get_box(box, key, v)
1014  type(box_t), intent(inout) :: box
1015  integer, intent(in) :: key
1017  real(dp), intent(inout) :: v(:, dtimes(:))
1018  integer :: i
1019 
1020  i = af_stencil_index(box, key)
1021  if (i == af_stencil_none) error stop "Stencil not stored"
1022 
1023  select case (box%stencils(i)%shape)
1024  case (af_stencil_357)
1025  call stencil_get_357(box, box%stencils(i), v)
1026  case default
1027  error stop "Unknown stencil"
1028  end select
1029  end subroutine af_stencil_get_box
1030 
1032  subroutine stencil_get_357(box, stencil, v)
1033  type(box_t), intent(inout) :: box
1034  type(stencil_t), intent(in) :: stencil
1036  real(dp), intent(inout) :: v(:, DTIMES(:))
1037 #if NDIM == 2
1038  real(dp) :: c_cyl(size(v, 1)), rfac(2, box%n_cell)
1039 #endif
1040  integer :: IJK, nc
1041 
1042  nc = box%n_cell
1043 
1044  if (size(v, 1) /= af_stencil_sizes(stencil%shape)) &
1045  error stop "Argument v has wrong size"
1046 
1047  if (stencil%stype == stencil_constant) then
1048  do kji_do(1, nc)
1049  v(:, ijk) = stencil%c
1050  end do; close_do
1051  else
1052  v = stencil%v
1053  end if
1054 
1055 #if NDIM == 2
1056  if (stencil%cylindrical_gradient) then
1057  ! Correct for cylindrical coordinates, assuming the elements correspond
1058  ! to a gradient operation
1059  call af_cyl_flux_factors(box, rfac)
1060  do kji_do(1, nc)
1061  c_cyl(2:3) = rfac(1:2, i) * v(2:3, ijk)
1062  c_cyl(1) = v(1, ijk) - (c_cyl(2) - v(2, ijk)) &
1063  - (c_cyl(3) - v(3, ijk))
1064  c_cyl(4:) = v(4:, ijk)
1065  v(:, ijk) = c_cyl
1066  end do; close_do
1067  end if
1068 #endif
1069  end subroutine stencil_get_357
1070 
1071 end module m_af_stencil
Subroutine for setting a stencil on a box.
This module contains functionality for dealing with numerical stencils.
Definition: m_af_stencil.f90:3
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Definition: m_af_types.f90:4
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
Type for storing a numerical stencil for a box.
Definition: m_af_types.f90:260