ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/neighborLists.F90
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 3 months ago) by gezelter
File size: 8943 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

# Content
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.3 2005-04-15 22:03:48 gezelter 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 !! Parameter for size > # of long range particles neighbor list
65 !! should be.
66 integer, parameter :: listMultiplier = 80
67 !! Maximum number of times we should reallocate neighbor list.
68 integer, parameter :: maxAllocations = 5
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 !--------------MODULE ACCESS-------------------------->
81 public :: expandNeighborList
82 public :: checkNeighborList
83 public :: saveNeighborList
84 public :: getNeighborListSize
85
86 contains
87
88
89 subroutine expandNeighborList(nGroups, error)
90 integer, intent(out) :: error
91 integer :: nGroups
92 integer :: alloc_error
93 integer :: oldSize = 0
94 integer :: newSize = 0
95 integer :: i
96 integer, dimension(:), pointer :: newList => null()
97 error = 0
98
99 !! First time through we should allocate point and list.
100 !! If one is associated and one is not, something is wrong
101 !! and return a error.
102
103 #ifdef IS_MPI !! // MPI
104 if (.not. associated(point) .and. &
105 .not. associated(list) ) then
106 allocate(point(getNgroupsInRow(plan_group_row)+1),stat=alloc_error)
107 if (alloc_error /= 0) then
108 write(default_error,*) &
109 "ExpandNeighborLists: Error in allocating MPI point"
110 error = -1
111 return
112 end if
113 allocate(list(listMultiplier * getNgroupsInCol(plan_group_col)),stat=alloc_error)
114 if (alloc_error /= 0) then
115 write(default_error,*) &
116 "ExpandNeighborLists: Error in allocating MPI list"
117 error = -1
118 return
119 end if
120 listSize = size(list)
121 nAllocations = nAllocations + 1
122 return
123 end if
124 #else !! // NONMPI
125 if (.not. associated(point) .and. &
126 .not. associated(list) ) then
127 allocate(point(nGroups),stat=alloc_error)
128 if (alloc_error /= 0) then
129 write(default_error,*) &
130 "ExpandNeighborLists: Error in allocating point"
131 error = -1
132 return
133 end if
134 allocate(list(listMultiplier * nGroups),stat=alloc_error)
135 if (alloc_error /= 0) then
136 write(default_error,*) &
137 "ExpandNeighborLists: Error in allocating list"
138 error = -1
139 return
140 end if
141 listSize = size(list)
142 nAllocations = nAllocations + 1
143 return
144 end if
145 #endif !! //MPI
146
147 ! Expand the neighbor list
148
149 ! Check to see if we have exceeded the maximum number of allocations.
150 if (nAllocations > maxAllocations) then
151 write(default_error,*) &
152 "ExpandNeighborList: exceeded maximum number of re-allocations"
153 error = -1
154 return
155 else !! Expand the list.
156 oldSize = size(list)
157
158
159 #ifdef IS_MPI !! MPI
160 newSize = listMultiplier * getNgroupsInCol(plan_group_col) + oldSize
161 #else
162 newSize = listMultiplier * nGroups + oldSize
163 #endif !! MPI
164
165
166
167 allocate(newList(newSize), stat=alloc_error)
168 if (alloc_error /= 0) then
169 write(*,*) "Error allocating new neighborlist"
170 error = -1
171 return
172 end if
173
174 !! Copy old list to new list
175 do i = 1, oldSize
176 newList(i) = list(i)
177 end do
178 !! Free old list
179 if(associated(list)) deallocate(list,stat=alloc_error)
180 if (alloc_error /= 0) then
181 error = -1
182 return
183 end if
184
185 !! Point list at new list
186
187 list => newList
188 end if
189
190 nAllocations = nAllocations + 1
191 listSize = size(list)
192 end subroutine expandNeighborList
193
194 !! checks to see if any long range particle has moved
195 !! through the neighbor list skin thickness.
196 subroutine checkNeighborList(nGroups, q, listSkin, update_nlist)
197 integer :: nGroups
198 real(kind = dp), intent(in) :: listSkin
199 real( kind = dp ), dimension(:,:) :: q
200 integer :: i
201 real( kind = DP ) :: dispmx
202 logical, intent(out) :: update_nlist
203 real( kind = DP ) :: dispmx_tmp
204
205 dispmx = 0.0E0_DP
206 !! calculate the largest displacement of any atom in any direction
207
208 !! If we have changed the particle idents, then we need to update
209 if (.not. allocated(q0)) then
210 update_nlist = .true.
211 return
212 end if
213
214 if (size(q0,2) /= nGroups) then
215 update_nlist = .true.
216 return
217 end if
218
219
220 #ifdef MPI
221
222 dispmx_tmp = 0.0E0_DP
223 do i = 1, nGroups
224 dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx_tmp )
225 dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx_tmp )
226 dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx_tmp )
227 end do
228 call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, &
229 mpi_max,mpi_comm_world,mpi_err)
230
231 #else
232
233 dispmx = 0.0_DP
234 do i = 1, nGroups
235 dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
236 dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
237 dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
238 end do
239
240 #endif
241
242 !! a conservative test of list skin crossings
243 dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
244
245 update_nlist = (dispmx.gt.listSkin)
246
247 end subroutine checkNeighborList
248
249
250 !! Saves neighbor list for comparison in check.
251 !! Save_neighborList will work even if the number of
252 !! local atoms has changed.
253 subroutine saveNeighborList(nGroups, q)
254
255 integer :: nGroups
256 real(kind = dp ), dimension(3,nGroups), intent(in) :: q
257 integer :: list_size
258
259 !! get size of list
260 list_size = nGroups
261
262 if (.not. allocated(q0)) then
263 allocate(q0(3,list_size))
264 else if( list_size > size(q0,2)) then
265 deallocate(q0)
266 allocate(q0(3,list_size))
267 endif
268 q0 = q
269 end subroutine saveNeighborList
270
271
272 function getNeighborListSize() result(returnListSize)
273 integer :: returnListSize
274 returnListSize = listSize
275 end function getNeighborListSize
276
277 end module neighborLists