11 public :: a2_particles_to_grid
18 subroutine a2_particles_to_grid(tree, iv, coords, weights, n_particles, &
19 order, id_guess, density, fill_gc)
22 use m_a2_utils, only: a2_get_id_at, a2_tree_clear_ghostcells
23 type(
a2_t),
intent(inout) :: tree
24 integer,
intent(in) :: iv
25 integer,
intent(in) :: n_particles
26 real(dp),
intent(in) :: coords(2, n_particles)
27 real(dp),
intent(in) :: weights(n_particles)
28 integer,
intent(in) :: order
30 integer,
intent(inout),
optional :: id_guess(n_particles)
32 logical,
intent(in),
optional :: density
34 logical,
intent(in),
optional :: fill_gc
37 integer :: current_thread, current_work
38 integer :: threads_left, work_left
39 integer,
allocatable :: ids(:)
40 integer,
allocatable :: npart_per_box(:)
41 integer,
allocatable :: box_threads(:)
42 integer,
allocatable :: threads(:)
44 logical :: fill_gc_at_end
46 allocate(ids(n_particles))
47 allocate(npart_per_box(-1:tree%highest_id))
48 allocate(box_threads(tree%highest_id))
49 allocate(threads(n_particles))
53 if (
present(density)) as_density = density
54 fill_gc_at_end = .true.
55 if (
present(fill_gc)) fill_gc_at_end = fill_gc
57 if (
present(id_guess))
then 60 ids(n) = a2_get_id_at(tree, coords(:, n), &
63 npart_per_box(ids(n)) = npart_per_box(ids(n)) + 1
69 ids(n) = a2_get_id_at(tree, coords(:, n))
70 npart_per_box(ids(n)) = npart_per_box(ids(n)) + 1
75 if (sum(npart_per_box(-1:0)) > 0)
then 76 print *,
"a2_particles_to_grid: some are outside domain" 79 if (ids(n) <= af_no_box)
then 80 print *, n, coords(:, n)
88 error stop
"a2_particles_to_grid: some are outside domain" 91 threads_left = af_get_max_threads()
94 work_left = n_particles
96 do m = 1, tree%highest_id
97 box_threads(m) = current_thread
98 current_work = current_work + npart_per_box(m)
100 if (current_work > work_left/threads_left)
then 101 current_thread = current_thread + 1
102 threads_left = threads_left - 1
103 work_left = work_left - current_work
109 do n = 1, n_particles
110 threads(n) = box_threads(ids(n))
115 call a2_tree_clear_ghostcells(tree, iv)
119 call particles_to_grid_0(tree, iv, coords, weights, ids, &
120 threads, n_particles, as_density)
122 call particles_to_grid_1(tree, iv, coords, weights, ids, &
123 threads, n_particles, as_density)
124 call tree_add_from_ghostcells(tree, iv)
126 error stop
"a2_particles_to_grid: Invalid interpolation order" 129 call a2_restrict_tree(tree, iv)
131 if (fill_gc_at_end)
then 132 if (.not. tree%has_cc_method(iv))
then 133 print *,
"Variable with index ", iv,
"has no cc_method" 134 print *,
"do this with call a2_set_cc_methods(tree, iv, ...)" 135 error stop
"a2_particles_to_grid: no ghost cell method defined" 137 call a2_gc_tree(tree, iv)
140 end subroutine a2_particles_to_grid
142 subroutine particles_to_grid_0(tree, iv, coords, weights, ids, &
143 threads, n_particles, density)
145 type(
a2_t),
intent(inout) :: tree
146 integer,
intent(in) :: iv
147 integer,
intent(in) :: n_particles
148 real(dp),
intent(in) :: coords(2, n_particles)
149 real(dp),
intent(in) :: weights(n_particles)
150 integer,
intent(in) :: ids(n_particles)
151 integer,
intent(in) :: threads(n_particles)
152 logical,
intent(in) :: density
153 integer :: n, thread_id, ix(2)
156 thread_id = omp_get_thread_num()
158 do n = 1, n_particles
159 if (threads(n) /= thread_id) cycle
161 ix = a2_cc_ix(tree%boxes(ids(n)), coords(:, n))
165 where (ix < 1) ix = 1
166 where (ix > tree%n_cell) ix = tree%n_cell
169 tree%boxes(ids(n))%cc(ix(1), ix(2), iv) = &
170 tree%boxes(ids(n))%cc(ix(1), ix(2), iv) + &
171 weights(n) / tree%boxes(ids(n))%dr**2
173 tree%boxes(ids(n))%cc(ix(1), ix(2), iv) = &
174 tree%boxes(ids(n))%cc(ix(1), ix(2), iv) + &
179 end subroutine particles_to_grid_0
183 subroutine particles_to_grid_1(tree, iv, coords, weights, ids, &
184 threads, n_particles, density)
186 type(
a2_t),
intent(inout) :: tree
187 integer,
intent(in) :: iv
188 integer,
intent(in) :: n_particles
189 real(dp),
intent(in) :: coords(2, n_particles)
190 real(dp),
intent(in) :: weights(n_particles)
191 integer,
intent(in) :: ids(n_particles)
192 integer,
intent(in) :: threads(n_particles)
193 logical,
intent(in) :: density
194 real(dp) :: tmp(2), inv_dr
195 real(dp) :: wu(2), wl(2), w(2, 2)
196 integer :: id, ix(2), n, thread_id
199 thread_id = omp_get_thread_num()
201 do n = 1, n_particles
202 if (threads(n) /= thread_id) cycle
205 inv_dr = 1.0_dp/tree%boxes(id)%dr
206 tmp = (coords(:, n) - tree%boxes(id)%r_min) * inv_dr + 0.5_dp
211 w(:, 1) = [wl(1), wu(1)] * wl(2)
212 w(:, 2) = [wl(1), wu(1)] * wu(2)
216 tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) = &
217 tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) + &
218 w * weights(n) / tree%boxes(id)%dr**2
220 tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) = &
221 tree%boxes(id)%cc(ix(1):ix(1)+1, ix(2):ix(2)+1, iv) + &
227 end subroutine particles_to_grid_1
229 subroutine tree_add_from_ghostcells(tree, iv)
230 type(
a2_t),
intent(inout) :: tree
231 integer,
intent(in) :: iv
232 integer :: lvl, i, id
235 do lvl = 1, tree%highest_lvl
237 do i = 1,
size(tree%lvls(lvl)%leaves)
238 id = tree%lvls(lvl)%leaves(i)
239 call add_from_ghostcells(tree%boxes, id, iv)
244 end subroutine tree_add_from_ghostcells
246 subroutine add_from_ghostcells(boxes, id, iv)
247 type(
box2_t),
intent(inout) :: boxes(:)
248 integer,
intent(in) :: id
249 integer,
intent(in) :: iv
250 integer :: i, j, i0(2), i1(2)
251 integer :: n0(2), n1(2), nb_id, nc
254 nc = boxes(id)%n_cell
256 do j = -1, 1;
do i = -1, 1
257 if (all([i, j] == 0)) cycle
259 nb_id = boxes(id)%neighbor_mat(i, j)
263 if (nb_id <= af_no_box)
then 265 else if (nb_id > af_no_box)
then 266 if (a2_has_children(boxes(nb_id))) copy_own = .true.
280 elsewhere ([i, j] == -1)
286 boxes(id)%cc(i0(1):i1(1), i0(2):i1(2), iv) = &
287 boxes(id)%cc(i0(1):i1(1), i0(2):i1(2), iv) + &
288 boxes(id)%cc(n0(1):n1(1), n0(2):n1(2), iv)
298 elsewhere ([i, j] == -1)
304 boxes(id)%cc(i0(1):i1(1), i0(2):i1(2), iv) = &
305 boxes(id)%cc(i0(1):i1(1), i0(2):i1(2), iv) + &
306 boxes(nb_id)%cc(n0(1):n1(1), n0(2):n1(2), iv)
309 end subroutine add_from_ghostcells
This module contains all kinds of different 'helper' routines for Afivo. If the number of routines fo...
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
Type which stores all the boxes and levels, as well as some information about the number of boxes...
This module contains the basic types and constants that are used in the 2-dimensional version of Afiv...
This module contains routines for restriction: going from fine to coarse variables.
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 routines related to , which can interpolate 'to' the grid and 'from' the grid (u...