ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/neighborLists.F90
Revision: 334
Committed: Thu Mar 13 17:45:54 2003 UTC (21 years, 4 months ago) by gezelter
File size: 6433 byte(s)
Log Message:
fixes for the exclude lists

File Contents

# Content
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.4 2003-03-13 17:45:54 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 #ifndef IS_MPI !!/Non MPI
62 if (.not. associated(point) .and. &
63 .not. associated(list) ) then
64 allocate(point(natoms),stat=alloc_error)
65 if (alloc_error /= 0) then
66 error = -1
67 return
68 end if
69 allocate(list(listMultiplier * natoms),stat=alloc_error)
70 if (alloc_error /= 0) then
71 error = -1
72 return
73 end if
74 nAllocations = nAllocations + 1
75 else
76 error = -1
77 return
78 end if
79 #else !!// MPI
80 if (.not. associated(point) .and. &
81 .not. associated(list) ) then
82 allocate(point(getNRow(plan_row)),stat=alloc_error)
83 if (alloc_error /= 0) then
84 error = -1
85 return
86 end if
87 allocate(list(listMultiplier * getNCol(plan_col)),stat=alloc_error)
88 if (alloc_error /= 0) then
89 error = -1
90 return
91 end if
92 nAllocations = nAllocations + 1
93 return
94 else
95 error = -1
96 return
97 end if
98 #endif !! //MPI
99
100 ! Expand the neighbor list
101
102 ! Check to see if we have exceeded the maximum number of allocations.
103 if (nAllocations > maxAllocations) then
104 error = -1
105 return
106 else !! Expand the list.
107 oldSize = size(list)
108
109 #ifndef IS_MPI !!Not MPI
110 newSize = listMultiplier * natoms + oldSize
111 allocate(newList(newSize), stat=alloc_error)
112 if (alloc_error /= 0) then
113 error = -1
114 return
115 end if
116 #else !! IS_MPI
117 newSize = listMultiplier * getNCol(plan_col) + oldSize
118 allocate(newList(newSize), stat = alloc_error)
119 if (alloc_error /= 0) then
120 error = -1
121 return
122 end if
123 #endif !! IS_MPI
124 !! Copy old list to new list
125 do i = 1, oldSize
126 newList(i) = list(i)
127 end do
128 !! Free old list
129 deallocate(list,stat=alloc_error)
130 if (alloc_error /= 0) then
131 error = -1
132 return
133 end if
134
135 !! Point list at new list
136 list => newList
137 end if
138
139 listSize = size(list)
140 end subroutine expandNeighborList
141
142 !! checks to see if any long range particle has moved
143 !! through the neighbor list skin thickness.
144 subroutine checkNeighborList(natoms, q, rcut, rlist, update_nlist)
145 integer :: natoms
146 real(kind = dp), intent(in) :: rcut, rlist
147 real( kind = dp ), dimension(:,:) :: q
148 integer :: i
149 real( kind = DP ) :: dispmx
150 logical, intent(out) :: update_nlist
151 real( kind = DP ) :: dispmx_tmp
152 real( kind = dp ) :: skin_thickness
153 integer :: nlocal
154
155 nlocal = natoms
156 skin_thickness = rlist - rcut
157 dispmx = 0.0E0_DP
158 !! calculate the largest displacement of any atom in any direction
159
160
161
162 #ifdef MPI
163
164 !! If we have changed the particle idents, then we need to update
165 if (.not. allocated(q0) .or. &
166 size(q0) /= nlocal) then
167 update_nlist = .true.
168 return
169 end if
170
171 dispmx_tmp = 0.0E0_DP
172 do i = 1, nlocal
173 dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx )
174 dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx )
175 dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx )
176 end do
177 call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, &
178 mpi_max,mpi_comm_world,mpi_err)
179 #else
180
181 do i = 1, nlocal
182 dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
183 dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
184 dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
185 end do
186 #endif
187
188 !! a conservative test of list skin crossings
189 dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
190
191 update_nlist = (dispmx.gt.(skin_thickness))
192
193 end subroutine checkNeighborList
194
195
196 !! Saves neighbor list for comparison in check.
197 !! Save_neighborList will work even if the number of
198 !! local atoms has changed.
199 subroutine saveNeighborList(q)
200 real(kind = dp ), dimension(:,:), intent(in) :: q
201 integer :: list_size
202
203 !! get size of list
204 list_size = size(q)
205
206 if (.not. allocated(q0)) then
207 allocate(q0(3,list_size))
208 else if( list_size > size(q0)) then
209 deallocate(q0)
210 allocate(q0(3,list_size))
211 endif
212 q0 = q
213 end subroutine saveNeighborList
214
215
216 function getNeighborListSize() result(returnListSize)
217 integer :: returnListSize
218 returnListSize = listSize
219 end function getNeighborListSize
220
221 end module neighborLists