Afivo  0.3
m_af_particles.f90
1 #include "cpp_macros.h"
2 
10  use m_af_types
11 
12  implicit none
13  private
14 
15  abstract interface
16 
17  subroutine subr_particle_id(ix, id)
18  import
19  integer, intent(in) :: ix
20  integer, intent(out) :: id
21  end subroutine subr_particle_id
22 
24  subroutine subr_particle_rw(ix, r, w)
25  import
26  integer, intent(in) :: ix
27  real(dp), intent(out) :: r(NDIM)
28  real(dp), intent(out) :: w
29  end subroutine subr_particle_rw
30  end interface
31 
32  public :: af_particles_to_grid
33 
34 contains
35 
39  subroutine af_particles_to_grid(tree, iv, n_particles, get_id, get_rw, &
40  order, density, fill_gc, iv_tmp, offset_particles)
41  use m_af_restrict, only: af_restrict_tree
42  use m_af_ghostcell, only: af_gc_tree
43  use m_af_utils, only: af_get_id_at, af_tree_clear_cc, af_tree_clear_ghostcells
44  type(af_t), intent(inout) :: tree
45  integer, intent(in) :: iv
46  integer, intent(in) :: n_particles
47  procedure(subr_particle_id) :: get_id
48  procedure(subr_particle_rw) :: get_rw
49  integer, intent(in) :: order
51  logical, intent(in), optional :: density
53  logical, intent(in), optional :: fill_gc
57  integer, intent(in), optional :: iv_tmp
59  integer, intent(in), optional :: offset_particles
60 
61  integer :: n, m, p_offset
62  integer :: current_thread, current_work
63  integer :: threads_left, work_left
64  integer, allocatable :: ids(:)
65  integer, allocatable :: npart_per_box(:)
66  integer, allocatable :: box_threads(:)
67  integer, allocatable :: threads(:)
68  real(dp) :: r(ndim), weight
69  logical :: as_density
70  logical :: fill_gc_at_end
71  logical :: use_tmp_var
72 
73  allocate(ids(n_particles))
74  allocate(npart_per_box(-1:tree%highest_id))
75  allocate(box_threads(tree%highest_id))
76  allocate(threads(n_particles))
77 
78  npart_per_box(:) = 0
79  as_density = .true.
80  if (present(density)) as_density = density
81  fill_gc_at_end = .true.
82  if (present(fill_gc)) fill_gc_at_end = fill_gc
83  use_tmp_var = .false.
84  if (present(iv_tmp)) then
85  if (iv_tmp > 0) use_tmp_var = .true.
86  end if
87  p_offset = 0
88  if (present(offset_particles)) p_offset = offset_particles
89 
90  if (use_tmp_var .and. .not. as_density) &
91  error stop "Use iv_tmp only for density = .true."
92 
93  !$omp parallel do reduction(+:npart_per_box)
94  do n = 1, n_particles
95  call get_id(p_offset + n, ids(n))
96  npart_per_box(ids(n)) = npart_per_box(ids(n)) + 1
97  end do
98  !$omp end parallel do
99 
100  if (sum(npart_per_box(-1:0)) > 0) then
101  print *, "af_particles_to_grid: some are outside domain"
102  m = 0
103  do n = 1, n_particles
104  if (ids(n) <= af_no_box) then
105  call get_rw(p_offset + n, r, weight)
106  print *, n, r
107  m = m + 1
108  end if
109  if (m > 10) then
110  print *, "..."
111  exit
112  end if
113  end do
114  error stop "af_particles_to_grid: some are outside domain"
115  end if
116 
117  threads_left = af_get_max_threads()
118  current_thread = 0
119  current_work = 0
120  work_left = n_particles
121 
122  do m = 1, tree%highest_id
123  box_threads(m) = current_thread
124  current_work = current_work + npart_per_box(m)
125 
126  if (current_work > work_left/threads_left) then
127  current_thread = current_thread + 1
128  threads_left = threads_left - 1
129  work_left = work_left - current_work
130  current_work = 0
131  end if
132  end do
133 
134  !$omp parallel do
135  do n = 1, n_particles
136  threads(n) = box_threads(ids(n))
137  end do
138  !$omp end parallel do
139 
140  ! Set density to zero in ghost cells
141  call af_tree_clear_ghostcells(tree, iv)
142 
143  ! Set density to zero in all cells
144  if (use_tmp_var) call af_tree_clear_cc(tree, iv_tmp)
145 
146  select case (order)
147  case (0)
148  if (use_tmp_var) then
149  call particles_to_grid_0(tree, iv_tmp, get_rw, ids, &
150  threads, n_particles, .false., p_offset)
151  call add_as_density(tree, iv_tmp, iv)
152  else
153  call particles_to_grid_0(tree, iv, get_rw, ids, &
154  threads, n_particles, as_density, p_offset)
155  end if
156  case (1)
157  if (use_tmp_var) then
158  call particles_to_grid_1(tree, iv_tmp, get_rw, ids, &
159  threads, n_particles, .false., p_offset)
160  call tree_add_from_ghostcells(tree, iv_tmp)
161  call add_as_density(tree, iv_tmp, iv)
162  else
163  call particles_to_grid_1(tree, iv, get_rw, ids, &
164  threads, n_particles, as_density, p_offset)
165  call tree_add_from_ghostcells(tree, iv)
166  end if
167  case default
168  error stop "af_particles_to_grid: Invalid interpolation order"
169  end select
170 
171  call af_restrict_tree(tree, [iv])
172 
173  if (fill_gc_at_end) then
174  if (.not. tree%has_cc_method(iv)) then
175  print *, "Variable with index ", iv, "has no cc_method"
176  print *, "do this with call af_set_cc_methods(tree, iv, ...)"
177  error stop "af_particles_to_grid: no ghost cell method defined"
178  else
179  call af_gc_tree(tree, [iv])
180  end if
181  end if
182  end subroutine af_particles_to_grid
183 
184  subroutine particles_to_grid_0(tree, iv, get_rw, ids, &
185  threads, n_particles, density, p_offset)
186  use omp_lib
187  type(af_t), intent(inout) :: tree
188  integer, intent(in) :: iv
189  integer, intent(in) :: n_particles
190  procedure(subr_particle_rw) :: get_rw
191  integer, intent(in) :: ids(n_particles)
192  integer, intent(in) :: threads(n_particles)
193  logical, intent(in) :: density
194  integer, intent(in) :: p_offset
195  integer :: n, thread_id, ix(NDIM)
196  real(dp) :: r(NDIM), weight, inv_volume
197 
198  !$omp parallel private(n, thread_id, ix, inv_volume, r, weight)
199  thread_id = omp_get_thread_num()
200 
201  do n = 1, n_particles
202  if (threads(n) /= thread_id) cycle
203  ! Handle this particle
204  call get_rw(p_offset + n, r, weight)
205 
206  ix = af_cc_ix(tree%boxes(ids(n)), r)
207 
208  ! Fix indices for points exactly on the boundaries of a box (which could
209  ! get a ghost cell index)
210  where (ix < 1) ix = 1
211  where (ix > tree%n_cell) ix = tree%n_cell
212 
213  if (density) then
214 #if NDIM == 2
215  if (tree%coord_t == af_cyl) then
216  inv_volume = 1 / af_cyl_volume_cc(tree%boxes(ids(n)), ix(1))
217  else
218  ! Cartesian
219  inv_volume = 1 / product(tree%boxes(ids(n))%dr)
220  end if
221 #else
222  inv_volume = 1 / product(tree%boxes(ids(n))%dr)
223 #endif
224 
225  tree%boxes(ids(n))%cc(dindex(ix), iv) = &
226  tree%boxes(ids(n))%cc(dindex(ix), iv) + &
227  weight * inv_volume
228  else
229  tree%boxes(ids(n))%cc(dindex(ix), iv) = &
230  tree%boxes(ids(n))%cc(dindex(ix), iv) + &
231  weight
232  end if
233  end do
234  !$omp end parallel
235  end subroutine particles_to_grid_0
236 
239  subroutine particles_to_grid_1(tree, iv, get_rw, ids, &
240  threads, n_particles, density, p_offset)
241  use omp_lib
242  type(af_t), intent(inout) :: tree
243  integer, intent(in) :: iv
244  integer, intent(in) :: n_particles
245  procedure(subr_particle_rw) :: get_rw
246  integer, intent(in) :: ids(n_particles)
247  integer, intent(in) :: threads(n_particles)
248  logical, intent(in) :: density
249  integer, intent(in) :: p_offset
250  real(dp) :: tmp(NDIM), inv_dr(NDIM)
251  real(dp) :: wu(NDIM), wl(NDIM), w(DTIMES(2))
252  real(dp) :: inv_volume, r(NDIM), weight
253  integer :: id, ix(NDIM), n, thread_id
254 
255  if (tree%coord_t == af_cyl .and. density) &
256  error stop "For cylindrical coordinates, use iv_tmp"
257 
258  !$omp parallel private(n, inv_dr, tmp, thread_id, ix, id, wu, wl, &
259  !$omp& w, inv_volume, r, weight)
260  thread_id = omp_get_thread_num()
261 
262  do n = 1, n_particles
263  if (threads(n) /= thread_id) cycle
264 
265  call get_rw(p_offset + n, r, weight)
266 
267  id = ids(n)
268  inv_dr = 1.0_dp/tree%boxes(id)%dr
269  tmp = (r - tree%boxes(id)%r_min) * inv_dr + 0.5_dp
270  ix = floor(tmp)
271  wu = tmp - ix
272  wl = 1 - wu
273 
274 #if NDIM == 1
275  w(1) = wl(1)
276  w(2) = wu(1)
277 #elif NDIM == 2
278  w(:, 1) = [wl(1), wu(1)] * wl(2)
279  w(:, 2) = [wl(1), wu(1)] * wu(2)
280 #elif NDIM == 3
281  w(:, 1, 1) = [wl(1), wu(1)] * wl(2) * wl(3)
282  w(:, 2, 1) = [wl(1), wu(1)] * wu(2) * wl(3)
283  w(:, 1, 2) = [wl(1), wu(1)] * wl(2) * wu(3)
284  w(:, 2, 2) = [wl(1), wu(1)] * wu(2) * wu(3)
285 #endif
286 
287  ! Linear interpolation
288  if (density) then
289  inv_volume = 1 / product(tree%boxes(ids(n))%dr)
290 #if NDIM == 1
291  tree%boxes(id)%cc(ix(1):ix(1)+1, iv) = &
292  tree%boxes(id)%cc(ix(1):ix(1)+1, iv) + &
293  w * inv_volume * weight
294 #elif NDIM == 2
295  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) = &
296  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) + &
297  w * inv_volume * weight
298 #elif NDIM == 3
299  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, ix(3):ix(3)+1, iv) = &
300  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, ix(3):ix(3)+1, iv) + &
301  w * inv_volume * weight
302 #endif
303  else
304 #if NDIM == 1
305  tree%boxes(id)%cc(ix(1):ix(1)+1, iv) = &
306  tree%boxes(id)%cc(ix(1):ix(1)+1, iv) + w * weight
307 #elif NDIM == 2
308  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) = &
309  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) + &
310  w * weight
311 #elif NDIM == 3
312  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, ix(3):ix(3)+1, iv) = &
313  tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, ix(3):ix(3)+1, iv) + &
314  w * weight
315 #endif
316  end if
317  end do
318  !$omp end parallel
319 
320  end subroutine particles_to_grid_1
321 
322  subroutine tree_add_from_ghostcells(tree, iv)
323  type(af_t), intent(inout) :: tree
324  integer, intent(in) :: iv
325  integer :: lvl, i, id
326 
327  !$omp parallel private(lvl, i, id)
328  do lvl = 1, tree%highest_lvl
329  !$omp do
330  do i = 1, size(tree%lvls(lvl)%leaves)
331  id = tree%lvls(lvl)%leaves(i)
332  call add_from_ghostcells(tree%boxes, id, iv)
333  end do
334  !$omp end do nowait
335  end do
336  !$omp end parallel
337  end subroutine tree_add_from_ghostcells
338 
339  subroutine add_from_ghostcells(boxes, id, iv)
340  type(box_t), intent(inout) :: boxes(:)
341  integer, intent(in) :: id
342  integer, intent(in) :: iv
343  integer :: IJK, i0(NDIM), i1(NDIM)
344  integer :: n0(NDIM), n1(NDIM), nb_id, nc
345  logical :: copy_own
346 
347  nc = boxes(id)%n_cell
348 
349  do kji_do(-1,1)
350  if (all([ijk] == 0)) cycle
351 
352  nb_id = boxes(id)%neighbor_mat(ijk)
353 
354  copy_own = .false.
355 
356  if (nb_id <= af_no_box) then
357  copy_own = .true.
358  else if (nb_id > af_no_box) then
359  if (af_has_children(boxes(nb_id))) copy_own = .true.
360  end if
361 
362  if (copy_own) then
363  ! Physical boundary
364  i0 = 1
365  i1 = nc
366  n0 = 1
367  n1 = nc
368 
369  where ([ijk] == 1)
370  i0 = nc
371  n0 = nc+1
372  n1 = nc+1
373  elsewhere ([ijk] == -1)
374  i1 = 1
375  n0 = 0
376  n1 = 0
377  end where
378 
379  boxes(id)%cc(dslice(i0, i1), iv) = &
380  boxes(id)%cc(dslice(i0, i1), iv) + &
381  boxes(id)%cc(dslice(n0, n1), iv)
382  else
383  i0 = 1
384  i1 = nc
385  n0 = 1
386  n1 = nc
387  where ([ijk] == 1)
388  i0 = nc
389  n0 = 0
390  n1 = 0
391  elsewhere ([ijk] == -1)
392  i1 = 1
393  n0 = nc+1
394  n1 = nc+1
395  end where
396 
397  boxes(id)%cc(dslice(i0, i1), iv) = &
398  boxes(id)%cc(dslice(i0, i1), iv) + &
399  boxes(nb_id)%cc(dslice(n0, n1), iv)
400  end if
401  end do; close_do
402  end subroutine add_from_ghostcells
403 
405  subroutine add_as_density(tree, iv_from, iv_to)
406  type(af_t), intent(inout) :: tree
407  integer, intent(in) :: iv_from
408  integer, intent(in) :: iv_to
409  integer :: lvl, i, id
410 
411  !$omp parallel private(lvl, i, id)
412  do lvl = 1, tree%highest_lvl
413  !$omp do
414  do i = 1, size(tree%lvls(lvl)%leaves)
415  id = tree%lvls(lvl)%leaves(i)
416  call add_as_density_box(tree%boxes(id), iv_from, iv_to)
417  end do
418  !$omp end do nowait
419  end do
420  !$omp end parallel
421  end subroutine add_as_density
422 
424  subroutine add_as_density_box(box, iv_from, iv_to)
425  type(box_t), intent(inout) :: box
426  integer, intent(in) :: iv_from
427  integer, intent(in) :: iv_to
428  real(dp) :: inv_volume
429 #if NDIM == 2
430  integer :: i
431  real(dp), parameter :: twopi = 2 * acos(-1.0_dp)
432  real(dp) :: radius, inv_cyl
433 #endif
434 
435  inv_volume = 1.0_dp / product(box%dr)
436 
437 #if NDIM == 2
438  if (box%coord_t == af_cyl) then
439  do i = 0, box%n_cell+1
440  ! abs() accounts for ghost cell on other side of r = 0
441  radius = abs(af_cyl_radius_cc(box, i))
442  inv_cyl = inv_volume / (twopi * radius)
443  box%cc(i, :, iv_to) = box%cc(i, :, iv_to) + &
444  box%cc(i, :, iv_from) * inv_cyl
445  end do
446  else
447  box%cc(dtimes(:), iv_to) = box%cc(dtimes(:), iv_to) + &
448  box%cc(dtimes(:), iv_from) * inv_volume
449  end if
450 #else
451  box%cc(dtimes(:), iv_to) = box%cc(dtimes(:), iv_to) + &
452  box%cc(dtimes(:), iv_from) * inv_volume
453 #endif
454  end subroutine add_as_density_box
455 
456 end module m_af_particles
To get a particle's coordinates and weight.
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
This module contains routines related to , which can interpolate 'to' the grid and 'from' the grid (u...
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:4
This module contains all kinds of different 'helper' routines for Afivo. If the number of routines fo...
Definition: m_af_utils.f90:5
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