ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/neighborLists.F90
Revision: 288
Committed: Thu Feb 27 18:42:52 2003 UTC (21 years, 6 months ago) by chuckv
File size: 6084 byte(s)
Log Message:
Changed lj_FF to use new neighbor list module.

File Contents

# User Rev Content
1 mmeineke 285 !! 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 chuckv 288 !! @version $Id: neighborLists.F90,v 1.2 2003-02-27 18:42:52 chuckv Exp $,
10 mmeineke 285
11     module neighborLists
12     use simulation
13     #ifdef IS_MPI
14     use mpiSimulation
15     #endif
16    
17     implicit none
18     PRIVATE
19    
20     !--------------MODULE VARIABLES---------------------->
21     !! Parameter for size > # of long range particles neighbor list
22     !! should be.
23     integer, parameter, :: listMultiplier = 80
24     !! Maximum number of times we should reallocate neighbor list.
25     integer, parameter, :: maxAllocations = 5
26     !! Number of times we have allocated the neighbor list.
27     integer, save :: nAllocations = 0
28     !! Pointer array to location in list for atom i.
29     integer, dimension(:),public, pointer :: point => null()
30     !! Neighbor list for atom i.
31     integer, dimension(:),public, pointer :: list => null()
32     !! Position array of previous positions for check. Allocated first time
33     !! into save_neighborList.
34     real( kind = dp ),dimension(:,:), allocatable , save :: q0 = 0.0_dp
35     !! Current list size
36     integer, save :: listSize
37     !--------------MODULE ACCESS-------------------------->
38     public :: expandList
39     public :: check
40     public :: save_neighborList
41 chuckv 288 public :: getNeighborListSize
42 mmeineke 285
43     contains
44    
45    
46     subroutine expandList(error)
47     integer, intent(out) :: error
48    
49     integer :: alloc_error
50     integer :: oldSize = 0
51     integer :: newSize = 0
52     integer,dimension(:), pointer :: new_list => null()
53     error = 0
54    
55    
56     !! First time through we should allocate point and list.
57     !! If one is associated and one is not, something is wrong
58     !! and return a error.
59     #ifndef IS_MPI !!/Non MPI
60     if (.not. associated(point) .and. &
61     .not. associated(list) ) then
62     allocate(point(getNlocal()),stat=alloc_error)
63     if (alloc_error /= 0) then
64     error = -1
65     return
66     end if
67     allocate(list(listMultiplier * getNlocal()),stat=alloc_error)
68     if (alloc_error /= 0) then
69     error = -1
70     return
71     end if
72     nAllocations = nAllocations + 1
73     else
74     error = -1
75     return
76     end if
77     #else !!// MPI
78     if (.not. associated(point) .and. &
79     .not. associated(list) ) then
80     allocate(point(getNRow()),stat=alloc_error)
81     if (alloc_error /= 0) then
82     error = -1
83     return
84     end if
85     allocate(list(listMultiplier * getNCol()),stat=alloc_error)
86     if (alloc_error /= 0) then
87     error = -1
88     return
89     end if
90     nAllocations = nAllocations + 1
91     return
92     else
93     error = -1
94     return
95     end if
96     #endif !! //MPI
97    
98     ! Expand the neighbor list
99    
100     ! Check to see if we have exceeded the maximum number of allocations.
101     if (nAllocations > maxAllocations) then
102     error = -1
103     return
104     else !! Expand the list.
105     oldSize = size(list)
106    
107     #ifndef IS_MPI !!Not MPI
108     newSize = listMultiplier * getNlocal() + oldSize
109     allocate(newList(newSize), stat=alloc_error)
110     if (alloc_error /= 0) then
111     error = -1
112     return
113     end if
114     #else !! IS_MPI
115     newSize = listMultiplier * getNCol() + oldSize
116     allocate(newList(newSize), stat = alloc_error)
117     if (alloc_error /= 0) then
118     error = -1
119     return
120     end if
121     #endif !! IS_MPI
122     !! Copy old list to new list
123     do i = 1, oldSize
124     newList(i) = list(i)
125     end do
126     !! Free old list
127     deallocate(list,stat=alloc_error)
128     if (alloc_error /= 0) then
129     error = -1
130     return
131     end if
132    
133     !! Point list at new list
134     list => newList
135     end if
136    
137     listSize = size(list)
138     end subroutine expandList
139    
140     !! checks to see if any long range particle has moved
141     !! through the neighbor list skin thickness.
142     subroutine check(q,update_nlist)
143     real( kind = dp ), dimension(:,:) :: q
144     integer :: i
145     real( kind = DP ) :: dispmx
146     logical, intent(out) :: update_nlist
147     real( kind = DP ) :: dispmx_tmp
148     real( kind = dp ) :: skin_thickness
149     integer :: nlocal
150    
151     nlocal = getNlocal()
152     skin_thickness = getRcut() - getRlist()
153     dispmx = 0.0E0_DP
154     !! calculate the largest displacement of any atom in any direction
155    
156    
157    
158     #ifdef MPI
159    
160     !! If we have changed the particle idents, then we need to update
161     if (.not. allocated(q0) .or. &
162     size(q0) /= nlocal) then
163     update_nlist = .true.
164     return
165     end if
166    
167     dispmx_tmp = 0.0E0_DP
168     do i = 1, nlocal
169     dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx )
170     dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx )
171     dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx )
172     end do
173     call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, &
174     mpi_max,mpi_comm_world,mpi_err)
175     #else
176    
177     do i = 1, nlocal
178     dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
179     dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
180     dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
181     end do
182     #endif
183    
184     !! a conservative test of list skin crossings
185     dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
186    
187     update_nlist = (dispmx.gt.(skin_thickness))
188    
189     end subroutine check
190    
191    
192     !! Saves neighbor list for comparison in check.
193     !! Save_neighborList will work even if the number of
194     !! local atoms has changed.
195     subroutine save_neighborList(q)
196     real(kind = dp ), dimension(:,:), intent(in) :: q
197     integer :: list_size
198    
199     !! get size of list
200     list_size = size(q)
201    
202     if (.not. allocated(q0)) then
203     allocate(q0(3,list_size))
204     else if( list_size > size(q0)) then
205     deallocate(q0)
206     allocate(q0(3,list_size))
207     endif
208     q0 = q
209     end subroutine save_neighborList
210    
211    
212 chuckv 288 function getNeighborListSize() result(returnListSize)
213     integer :: returnListSize
214     returnListSize = listSize
215     end function getNeighborListSize
216    
217    
218 mmeineke 285 end module neighborLists