1 |
!! Module neighborLists |
2 |
!! Impliments verlet neighbor lists for force modules. |
3 |
!! Automagically expands neighbor list if size too small |
4 |
!! up to maxAllocations times. If after maxAllocations we try to |
5 |
!! expand the neighbor list, we get an error message and quit. |
6 |
!! @author Charles F. Vardeman II |
7 |
!! @author Matthew Meineke |
8 |
!! @author J. Daniel Gezelter |
9 |
!! @version $Id: neighborLists.F90,v 1.7 2003-10-30 14:11:28 gezelter Exp $, |
10 |
|
11 |
module neighborLists |
12 |
|
13 |
use definitions |
14 |
#ifdef IS_MPI |
15 |
use mpiSimulation |
16 |
#endif |
17 |
|
18 |
implicit none |
19 |
PRIVATE |
20 |
|
21 |
!--------------MODULE VARIABLES----------------------> |
22 |
!! Parameter for size > # of long range particles neighbor list |
23 |
!! should be. |
24 |
integer, parameter :: listMultiplier = 80 |
25 |
!! Maximum number of times we should reallocate neighbor list. |
26 |
integer, parameter :: maxAllocations = 5 |
27 |
!! Number of times we have allocated the neighbor list. |
28 |
integer, save :: nAllocations = 0 |
29 |
!! Pointer array to location in list for atom i. |
30 |
integer, dimension(:),public, pointer :: point => null() |
31 |
!! Neighbor list for atom i. |
32 |
integer, dimension(:),public, pointer :: list => null() |
33 |
!! Position array of previous positions for check. Allocated first time |
34 |
!! into saveNeighborList. |
35 |
real( kind = dp ), dimension(:,:), allocatable, save :: q0 |
36 |
!! Current list size |
37 |
integer, save :: listSize |
38 |
!--------------MODULE ACCESS--------------------------> |
39 |
public :: expandNeighborList |
40 |
public :: checkNeighborList |
41 |
public :: saveNeighborList |
42 |
public :: getNeighborListSize |
43 |
|
44 |
contains |
45 |
|
46 |
|
47 |
subroutine expandNeighborList(natoms, error) |
48 |
integer, intent(out) :: error |
49 |
integer :: natoms |
50 |
integer :: alloc_error |
51 |
integer :: oldSize = 0 |
52 |
integer :: newSize = 0 |
53 |
integer :: i |
54 |
integer, dimension(:), pointer :: newList => null() |
55 |
error = 0 |
56 |
|
57 |
|
58 |
!! First time through we should allocate point and list. |
59 |
!! If one is associated and one is not, something is wrong |
60 |
!! and return a error. |
61 |
|
62 |
#ifdef IS_MPI !! // MPI |
63 |
if (.not. associated(point) .and. & |
64 |
.not. associated(list) ) then |
65 |
allocate(point(getNRow(plan_row)),stat=alloc_error) |
66 |
if (alloc_error /= 0) then |
67 |
write(default_error,*) & |
68 |
"ExpandNeighborLists: Error in allocating MPI point" |
69 |
error = -1 |
70 |
return |
71 |
end if |
72 |
allocate(list(listMultiplier * getNCol(plan_col)),stat=alloc_error) |
73 |
if (alloc_error /= 0) then |
74 |
write(default_error,*) & |
75 |
"ExpandNeighborLists: Error in allocating MPI list" |
76 |
error = -1 |
77 |
return |
78 |
end if |
79 |
listSize = size(list) |
80 |
nAllocations = nAllocations + 1 |
81 |
return |
82 |
end if |
83 |
#else !! // NONMPI |
84 |
if (.not. associated(point) .and. & |
85 |
.not. associated(list) ) then |
86 |
allocate(point(natoms),stat=alloc_error) |
87 |
if (alloc_error /= 0) then |
88 |
write(default_error,*) & |
89 |
"ExpandNeighborLists: Error in allocating point" |
90 |
error = -1 |
91 |
return |
92 |
end if |
93 |
allocate(list(listMultiplier * natoms),stat=alloc_error) |
94 |
if (alloc_error /= 0) then |
95 |
write(default_error,*) & |
96 |
"ExpandNeighborLists: Error in allocating list" |
97 |
error = -1 |
98 |
return |
99 |
end if |
100 |
listSize = size(list) |
101 |
nAllocations = nAllocations + 1 |
102 |
return |
103 |
end if |
104 |
#endif !! //MPI |
105 |
|
106 |
! Expand the neighbor list |
107 |
|
108 |
! Check to see if we have exceeded the maximum number of allocations. |
109 |
if (nAllocations > maxAllocations) then |
110 |
write(default_error,*) & |
111 |
"ExpandNeighborList: exceeded maximum number of re-allocations" |
112 |
error = -1 |
113 |
return |
114 |
else !! Expand the list. |
115 |
oldSize = size(list) |
116 |
|
117 |
|
118 |
#ifdef IS_MPI !! MPI |
119 |
newSize = listMultiplier * getNCol(plan_col) + oldSize |
120 |
#else |
121 |
newSize = listMultiplier * natoms + oldSize |
122 |
#endif !! MPI |
123 |
|
124 |
|
125 |
|
126 |
allocate(newList(newSize), stat=alloc_error) |
127 |
if (alloc_error /= 0) then |
128 |
write(*,*) "Error allocating new neighborlist" |
129 |
error = -1 |
130 |
return |
131 |
end if |
132 |
|
133 |
!! Copy old list to new list |
134 |
do i = 1, oldSize |
135 |
newList(i) = list(i) |
136 |
end do |
137 |
!! Free old list |
138 |
if(associated(list)) deallocate(list,stat=alloc_error) |
139 |
if (alloc_error /= 0) then |
140 |
error = -1 |
141 |
return |
142 |
end if |
143 |
|
144 |
!! Point list at new list |
145 |
|
146 |
list => newList |
147 |
end if |
148 |
|
149 |
nAllocations = nAllocations + 1 |
150 |
listSize = size(list) |
151 |
end subroutine expandNeighborList |
152 |
|
153 |
!! checks to see if any long range particle has moved |
154 |
!! through the neighbor list skin thickness. |
155 |
subroutine checkNeighborList(natoms, q, listSkin, update_nlist) |
156 |
integer :: natoms |
157 |
real(kind = dp), intent(in) :: listSkin |
158 |
real( kind = dp ), dimension(:,:) :: q |
159 |
integer :: i |
160 |
real( kind = DP ) :: dispmx |
161 |
logical, intent(out) :: update_nlist |
162 |
real( kind = DP ) :: dispmx_tmp |
163 |
integer :: nlocal |
164 |
|
165 |
|
166 |
nlocal = natoms |
167 |
dispmx = 0.0E0_DP |
168 |
!! calculate the largest displacement of any atom in any direction |
169 |
|
170 |
!! If we have changed the particle idents, then we need to update |
171 |
if (.not. allocated(q0)) then |
172 |
update_nlist = .true. |
173 |
return |
174 |
end if |
175 |
|
176 |
if (size(q0,2) /= nlocal) then |
177 |
update_nlist = .true. |
178 |
return |
179 |
end if |
180 |
|
181 |
|
182 |
#ifdef MPI |
183 |
|
184 |
dispmx_tmp = 0.0E0_DP |
185 |
do i = 1, nlocal |
186 |
dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx_tmp ) |
187 |
dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx_tmp ) |
188 |
dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx_tmp ) |
189 |
end do |
190 |
call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, & |
191 |
mpi_max,mpi_comm_world,mpi_err) |
192 |
|
193 |
#else |
194 |
|
195 |
dispmx = 0.0_DP |
196 |
do i = 1, nlocal |
197 |
dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx ) |
198 |
dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx ) |
199 |
dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx ) |
200 |
end do |
201 |
|
202 |
#endif |
203 |
|
204 |
|
205 |
!! a conservative test of list skin crossings |
206 |
dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx) |
207 |
|
208 |
update_nlist = (dispmx.gt.listSkin) |
209 |
|
210 |
end subroutine checkNeighborList |
211 |
|
212 |
|
213 |
!! Saves neighbor list for comparison in check. |
214 |
!! Save_neighborList will work even if the number of |
215 |
!! local atoms has changed. |
216 |
subroutine saveNeighborList(natoms, q) |
217 |
|
218 |
integer :: natoms |
219 |
real(kind = dp ), dimension(3,natoms), intent(in) :: q |
220 |
integer :: list_size |
221 |
|
222 |
|
223 |
|
224 |
!! get size of list |
225 |
list_size = natoms |
226 |
|
227 |
if (.not. allocated(q0)) then |
228 |
allocate(q0(3,list_size)) |
229 |
else if( list_size > size(q0,2)) then |
230 |
deallocate(q0) |
231 |
allocate(q0(3,list_size)) |
232 |
endif |
233 |
q0 = q |
234 |
end subroutine saveNeighborList |
235 |
|
236 |
|
237 |
function getNeighborListSize() result(returnListSize) |
238 |
integer :: returnListSize |
239 |
returnListSize = listSize |
240 |
end function getNeighborListSize |
241 |
|
242 |
end module neighborLists |