1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
!! publication of scientific results based in part on use of the |
11 |
!! program. An acceptable form of acknowledgement is citation of |
12 |
!! the article in which the program was described (Matthew |
13 |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
!! |
18 |
!! 2. Redistributions of source code must retain the above copyright |
19 |
!! notice, this list of conditions and the following disclaimer. |
20 |
!! |
21 |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
!! notice, this list of conditions and the following disclaimer in the |
23 |
!! documentation and/or other materials provided with the |
24 |
!! distribution. |
25 |
!! |
26 |
!! This software is provided "AS IS," without a warranty of any |
27 |
!! kind. All express or implied conditions, representations and |
28 |
!! warranties, including any implied warranty of merchantability, |
29 |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
!! be liable for any damages suffered by licensee as a result of |
32 |
!! using, modifying or distributing the software or its |
33 |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
!! damages, however caused and regardless of the theory of liability, |
37 |
!! arising out of the use of or inability to use software, even if the |
38 |
!! University of Notre Dame has been advised of the possibility of |
39 |
!! such damages. |
40 |
!! |
41 |
|
42 |
|
43 |
!! Module neighborLists |
44 |
!! Implements verlet neighbor lists for force modules. |
45 |
!! Automagically expands neighbor list if size too small |
46 |
!! up to maxAllocations times. If after maxAllocations we try to |
47 |
!! expand the neighbor list, we get an error message and quit. |
48 |
!! @author Charles F. Vardeman II |
49 |
!! @author Matthew Meineke |
50 |
!! @author J. Daniel Gezelter |
51 |
!! @version $Id: neighborLists.F90,v 1.6 2006-12-05 00:17:24 chuckv Exp $, |
52 |
|
53 |
module neighborLists |
54 |
|
55 |
use definitions |
56 |
#ifdef IS_MPI |
57 |
use mpiSimulation |
58 |
#endif |
59 |
|
60 |
implicit none |
61 |
PRIVATE |
62 |
|
63 |
!--------------MODULE VARIABLES----------------------> |
64 |
!! Default Parameter for size > # of long range particles neighbor list |
65 |
!! should be. It is listMultiplier x nGroups for size of list |
66 |
integer, save :: listMultiplier = 80 |
67 |
!! Maximum number of times we should reallocate neighbor list. |
68 |
integer, parameter :: maxAllocations = 10 |
69 |
!! Number of times we have allocated the neighbor list. |
70 |
integer, save :: nAllocations = 0 |
71 |
!! Pointer array to location in list for atom i. |
72 |
integer, dimension(:),public, pointer :: point => null() |
73 |
!! Neighbor list for atom i. |
74 |
integer, dimension(:),public, pointer :: list => null() |
75 |
!! Position array of previous positions for check. Allocated first time |
76 |
!! into saveNeighborList. |
77 |
real( kind = dp ), dimension(:,:), allocatable, save :: q0 |
78 |
!! Current list size |
79 |
integer, save :: listSize |
80 |
|
81 |
!--------------MODULE ACCESS--------------------------> |
82 |
public :: expandNeighborList |
83 |
public :: checkNeighborList |
84 |
public :: saveNeighborList |
85 |
public :: getNeighborListSize |
86 |
public :: setNeighbors |
87 |
|
88 |
contains |
89 |
|
90 |
|
91 |
subroutine expandNeighborList(nGroups, error) |
92 |
integer, intent(out) :: error |
93 |
integer :: nGroups |
94 |
integer :: alloc_error |
95 |
integer :: oldSize = 0 |
96 |
integer :: newSize = 0 |
97 |
integer :: i |
98 |
integer, dimension(:), pointer :: newList => null() |
99 |
error = 0 |
100 |
|
101 |
!! First time through we should allocate point and list. |
102 |
!! If one is associated and one is not, something is wrong |
103 |
!! and return a error. |
104 |
|
105 |
#ifdef IS_MPI !! // MPI |
106 |
if (.not. associated(point) .and. & |
107 |
.not. associated(list) ) then |
108 |
allocate(point(getNgroupsInRow(plan_group_row)+1),stat=alloc_error) |
109 |
if (alloc_error /= 0) then |
110 |
write(default_error,*) & |
111 |
"ExpandNeighborLists: Error in allocating MPI point" |
112 |
error = -1 |
113 |
return |
114 |
end if |
115 |
allocate(list(listMultiplier * getNgroupsInCol(plan_group_col)),stat=alloc_error) |
116 |
if (alloc_error /= 0) then |
117 |
write(default_error,*) & |
118 |
"ExpandNeighborLists: Error in allocating MPI list" |
119 |
error = -1 |
120 |
return |
121 |
end if |
122 |
listSize = size(list) |
123 |
nAllocations = nAllocations + 1 |
124 |
return |
125 |
end if |
126 |
#else !! // NONMPI |
127 |
if (.not. associated(point) .and. & |
128 |
.not. associated(list) ) then |
129 |
allocate(point(nGroups),stat=alloc_error) |
130 |
if (alloc_error /= 0) then |
131 |
write(default_error,*) & |
132 |
"ExpandNeighborLists: Error in allocating point" |
133 |
error = -1 |
134 |
return |
135 |
end if |
136 |
allocate(list(listMultiplier * nGroups),stat=alloc_error) |
137 |
if (alloc_error /= 0) then |
138 |
write(default_error,*) & |
139 |
"ExpandNeighborLists: Error in allocating list" |
140 |
error = -1 |
141 |
return |
142 |
end if |
143 |
listSize = size(list) |
144 |
nAllocations = nAllocations + 1 |
145 |
return |
146 |
end if |
147 |
#endif !! //MPI |
148 |
|
149 |
! Expand the neighbor list |
150 |
|
151 |
! Check to see if we have exceeded the maximum number of allocations. |
152 |
if (nAllocations > maxAllocations) then |
153 |
write(default_error,*) & |
154 |
"ExpandNeighborList: exceeded maximum number of re-allocations" |
155 |
error = -1 |
156 |
return |
157 |
else !! Expand the list. |
158 |
oldSize = size(list) |
159 |
|
160 |
|
161 |
#ifdef IS_MPI !! MPI |
162 |
newSize = listMultiplier * getNgroupsInCol(plan_group_col) + oldSize |
163 |
#else |
164 |
newSize = listMultiplier * nGroups + oldSize |
165 |
#endif !! MPI |
166 |
|
167 |
|
168 |
|
169 |
allocate(newList(newSize), stat=alloc_error) |
170 |
if (alloc_error /= 0) then |
171 |
write(*,*) "Error allocating new neighborlist" |
172 |
error = -1 |
173 |
return |
174 |
end if |
175 |
|
176 |
!! Copy old list to new list |
177 |
do i = 1, oldSize |
178 |
newList(i) = list(i) |
179 |
end do |
180 |
!! Free old list |
181 |
if(associated(list)) deallocate(list,stat=alloc_error) |
182 |
if (alloc_error /= 0) then |
183 |
error = -1 |
184 |
return |
185 |
end if |
186 |
|
187 |
!! Point list at new list |
188 |
|
189 |
list => newList |
190 |
end if |
191 |
|
192 |
nAllocations = nAllocations + 1 |
193 |
listSize = size(list) |
194 |
end subroutine expandNeighborList |
195 |
|
196 |
!! checks to see if any long range particle has moved |
197 |
!! through the neighbor list skin thickness. |
198 |
subroutine checkNeighborList(nGroups, q, listSkin, update_nlist) |
199 |
integer :: nGroups |
200 |
real(kind = dp), intent(in) :: listSkin |
201 |
real( kind = dp ), dimension(:,:) :: q |
202 |
integer :: i |
203 |
real( kind = DP ) :: dispmx |
204 |
logical, intent(out) :: update_nlist |
205 |
real( kind = DP ) :: dispmx_tmp |
206 |
|
207 |
dispmx = 0.0E0_DP |
208 |
!! calculate the largest displacement of any atom in any direction |
209 |
|
210 |
!! If we have changed the particle idents, then we need to update |
211 |
if (.not. allocated(q0)) then |
212 |
update_nlist = .true. |
213 |
return |
214 |
end if |
215 |
|
216 |
if (size(q0,2) /= nGroups) then |
217 |
update_nlist = .true. |
218 |
return |
219 |
end if |
220 |
|
221 |
|
222 |
#ifdef MPI |
223 |
|
224 |
dispmx_tmp = 0.0E0_DP |
225 |
do i = 1, nGroups |
226 |
dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx_tmp ) |
227 |
dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx_tmp ) |
228 |
dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx_tmp ) |
229 |
end do |
230 |
#ifdef SINGLE_PRECISION |
231 |
call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_real, & |
232 |
mpi_max,mpi_comm_world,mpi_err) |
233 |
#else |
234 |
call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, & |
235 |
mpi_max,mpi_comm_world,mpi_err) |
236 |
#endif |
237 |
|
238 |
#else |
239 |
|
240 |
dispmx = 0.0_DP |
241 |
do i = 1, nGroups |
242 |
dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx ) |
243 |
dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx ) |
244 |
dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx ) |
245 |
end do |
246 |
|
247 |
#endif |
248 |
|
249 |
!! a conservative test of list skin crossings |
250 |
dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx) |
251 |
|
252 |
update_nlist = (dispmx.gt.listSkin) |
253 |
|
254 |
end subroutine checkNeighborList |
255 |
|
256 |
|
257 |
!! Saves neighbor list for comparison in check. |
258 |
!! Save_neighborList will work even if the number of |
259 |
!! local atoms has changed. |
260 |
subroutine saveNeighborList(nGroups, q) |
261 |
|
262 |
integer :: nGroups |
263 |
real(kind = dp ), dimension(3,nGroups), intent(in) :: q |
264 |
integer :: list_size |
265 |
|
266 |
!! get size of list |
267 |
list_size = nGroups |
268 |
|
269 |
if (.not. allocated(q0)) then |
270 |
allocate(q0(3,list_size)) |
271 |
else if( list_size > size(q0,2)) then |
272 |
deallocate(q0) |
273 |
allocate(q0(3,list_size)) |
274 |
endif |
275 |
q0 = q |
276 |
end subroutine saveNeighborList |
277 |
|
278 |
|
279 |
function getNeighborListSize() result(returnListSize) |
280 |
integer :: returnListSize |
281 |
returnListSize = listSize |
282 |
end function getNeighborListSize |
283 |
|
284 |
|
285 |
subroutine setNeighbors(nNeighbors) |
286 |
integer, intent(in) :: nNeighbors |
287 |
listMultiplier = nNeighbors |
288 |
end subroutine setNeighbors |
289 |
|
290 |
end module neighborLists |