-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathget_neighbour_lists.f90
318 lines (307 loc) · 12.6 KB
/
get_neighbour_lists.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
!------------------------------------------------------------------------------!
! NDSPMHD: A Smoothed Particle (Magneto)Hydrodynamics code for (astrophysical) !
! fluid dynamics simulations in 1, 2 and 3 spatial dimensions. !
! !
! (c) 2002-2015 Daniel Price !
! !
! http://users.monash.edu.au/~dprice/ndspmhd !
! daniel.price@monash.edu -or- dprice@cantab.net (forwards to current address) !
! !
! NDSPMHD comes with ABSOLUTELY NO WARRANTY. !
! This is free software; and you are welcome to redistribute !
! it under the terms of the GNU General Public License !
! (see LICENSE file for details) and the provision that !
! this notice remains intact. If you modify this file, please !
! note section 2a) of the GPLv2 states that: !
! !
! a) You must cause the modified files to carry prominent notices !
! stating that you changed the files and the date of any change. !
! !
! ChangeLog: !
!------------------------------------------------------------------------------!
module get_neighbour_lists
implicit none
contains
!!-----------------------------------------------------------------------
!! Using the link list setup, compiles the neighbour list for the
!! current cell (this list is common to all particles in the cell)
!!
!! the list is returned in 'listneigh' (length nneigh)
!! the neigbouring cells to the current cell are returned in 'neighcell'
!!
!!-----------------------------------------------------------------------
subroutine get_neighbour_list(icell,listneigh,nneigh)
use dimen_mhd
use debug
use loguns
use linklist
use options
use part, only:ntotal
!
!--define local variables
!
implicit none
integer, intent(IN) :: icell
integer, dimension(:), intent(OUT) :: listneigh
integer, intent(OUT) :: nneigh
integer, dimension(3**ndim) :: neighcell
integer :: nneighcell,ipart,ncellsxy
integer :: i,j,k
logical :: debugging
!
!--allow for tracing flow (removed for speed)
!
! IF (trace) WRITE(iprint,*) ' Entering subroutine get_neighbour_list'
debugging = .false.
if (idebug(1:4).EQ.'neig') debugging = .true.
!
!--if cell is empty, neighbouring cells are irrelevant, so return
!
if (ifirstincell(icell).le.0) then
if (debugging) write(iprint,*) 'icell:',icell,' cell empty: no neighbours'
!listneigh = 0 ! unnecessary
nneigh = 0
return
endif
!
!--work out which cells current cell should interact with
! same template is used for all cells as a block of empty cells surrounds those with particles
! (ie. so there is always a cell to the left, above or whatever)
!
nneighcell = 2
neighcell(1) = icell ! always interacts with itself
neighcell(2) = icell + 1
if (ndim.ge.2) then
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1) - 1 ! above, left
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1) ! above
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1) + 1 ! above, right
endif
if (ndim.ge.3) then ! next block
ncellsxy = ncellsx(1)*ncellsx(2)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - 1 ! next block, left
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy ! next block
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + 1 ! next block, right
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1) - 1 ! next block, above left
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1) ! next block, above
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1) + 1 ! next block, above right
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1) - 1 ! next block, below left
nneighcell = nneighcell + 1 ! next block, below
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1) ! next block, below
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1) + 1 ! next block, below right
endif
if (debugging) then
print*,' current cell = ',icell,' first particle = ',ifirstincell(icell)
print*,' number of neighbouring cells = ',nneighcell
print*, (neighcell(i),i=1,nneighcell)
read*
endif
!
!--construct list of neighbours for particles in the current cell
! using the link lists from the neighbouring cells
!
j=0
do k = 1,nneighcell ! construct list of neighbours
if (neighcell(k).GT.ncells) then
print*,' k, neighcell = ',k,neighcell(1:nneighcell),ncells
stop 'get_neighbour_list: error: cell > ncells'
endif
ipart = ifirstincell(neighcell(k))
! PRINT*,'neighbouring cell = ',neighcell(k)
if (ipart.NE.-1) then
do while (ll(ipart).NE.-1)
j = j + 1
if (j.GT.size(listneigh)) then ! should reallocate array here
write(iprint,*) 'getneigh: # neighbours > array size:',size(listneigh)
write(iprint,*) 'ipart = ',ipart
call quit
endif
listneigh(j) = ipart
ipart = ll(ipart)
!--make sure we have a real particle
if (ipart.gt.ntotal .or. ipart.lt.0) then
write(iprint,*) 'ERROR in neighbour lists: ipart > ntotal ',ipart
endif
! print*,'listneigh ',j,'= ',listneigh(j)
enddo
j = j + 1
!--make sure we have a real particle
if (ipart.gt.ntotal .or. ipart.lt.0) then
write(iprint,*) 'ERROR in neighbour lists: ipart = ',ipart
endif
if (j.gt.size(listneigh)) then
write(iprint,*) 'getneigh: # neighbours > array size:',size(listneigh)
write(iprint,*) 'j = ',j
call quit
endif
listneigh(j) = ipart
endif
! IF (j.GT.0) THEN
! print*,'listneigh ',j,'= ',listneigh(j)
! ELSE
! print*,' no neighbours'
! ENDIF
enddo
nneigh = j
! print*,'number of neighbours = ',nneigh
! read*
return
end subroutine get_neighbour_list
!!-----------------------------------------------------------------------
!! Using the link list setup, compiles the neighbour list for the
!! current cell (this list is common to all particles in the cell)
!!
!! the list is returned in 'listneigh' (length nneigh)
!! the neigbouring cells to the current cell are returned in 'neighcell'
!!
!!-----------------------------------------------------------------------
subroutine get_neighbour_list_partial(icell,listneigh,nneigh)
use dimen_mhd
use debug
use loguns
use linklist
use options
!
!--define local variables
!
implicit none
integer, intent(IN) :: icell
integer, dimension(:), intent(OUT) :: listneigh
integer, intent(OUT) :: nneigh
integer, dimension(3**ndim) :: neighcell
integer :: nneighcell,ipart,ncellsxy
integer :: i,j,k
logical :: debugging
!
!--allow for tracing flow (removed for speed)
!
! IF (trace) WRITE(iprint,*) ' Entering subroutine get_neighbour_list'
debugging = .false.
if (idebug(1:4).EQ.'neig') debugging = .true.
!
!--if cell is empty, neighbouring cells are irrelevant, so return
!
if (ifirstincell(icell).le.0) then
if (debugging) write(iprint,*) 'icell:',icell,' cell empty: no neighbours'
!listneigh = 0 unnecessary
nneigh = 0
return
endif
!
!--work out which cells current cell should interact with
!
nneighcell = 3
neighcell(1) = icell - 1 ! always interacts with itself
neighcell(2) = icell
neighcell(3) = icell + 1
if (ndim.ge.2) then ! add incrementally to stop compiler errors if ndim < 2
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsx(1) + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsx(1) + 1
if (ndim.ge.3) then
ncellsxy = ncellsx(1)*ncellsx(2)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + ncellsx(1) + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell + ncellsxy - ncellsx(1) + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy + ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy + ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy + ncellsx(1) + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy + 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy - ncellsx(1) - 1
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy - ncellsx(1)
nneighcell = nneighcell + 1
neighcell(nneighcell) = icell - ncellsxy - ncellsx(1) + 1
endif
endif
if (debugging) then
print*,' current cell = ',icell,' first particle = ',ifirstincell(icell)
print*,' number of neighbouring cells = ',nneighcell
print*, (neighcell(i),i=1,nneighcell)
read*
endif
!
!--construct list of neighbours for particles in the current cell
! using the link lists from the neighbouring cells
!
j=0
do k = 1,nneighcell ! construct list of neighbours
if (neighcell(k).GT.ncells) then
print*,' k=',k,' neighcell = ',neighcell(1:nneighcell),'ncells=',ncells
stop 'get_neighbour_list: error: cell > ncells'
endif
ipart = ifirstincell(neighcell(k))
! PRINT*,'neighbouring cell = ',neighcell(k)
if (ipart.NE.-1) then
do while (ll(ipart).NE.-1)
j = j + 1
if (j.GT.size(listneigh)) then ! should reallocate array here
write(iprint,*) 'getneighpartial: # neighbours > array size:',size(listneigh)
write(iprint,*) 'icell = ',icell,' ncells = ',nneighcell,'=',neighcell(1:nneighcell)
write(iprint,*) 'neighbouring cell = ',neighcell(k)
write(iprint,*) 'listneigh = ',listneigh(:)
write(iprint,*) 'ipart = ',ipart
call quit
endif
listneigh(j) = ipart
ipart = ll(ipart)
! print*,'listneigh ',j,'= ',listneigh(j)
enddo
j = j + 1
listneigh(j) = ipart
endif
! IF (j.GT.0) THEN
! print*,'listneigh ',j,'= ',listneigh(j)
! ELSE
! print*,' no neighbours'
! ENDIF
enddo
nneigh = j
! print*,'number of neighbours = ',nneigh
! read*
return
end subroutine get_neighbour_list_partial
end module get_neighbour_lists