ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/neighborLists.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/neighborLists.F90 (file contents):
Revision 1197 by gezelter, Fri May 7 21:35:05 2004 UTC vs.
Revision 1198 by tim, Thu May 27 00:48:12 2004 UTC

# Line 6 | Line 6
6   !! @author Charles F. Vardeman II
7   !! @author Matthew Meineke
8   !! @author J. Daniel Gezelter
9 < !! @version $Id: neighborLists.F90,v 1.10 2004-05-07 21:35:05 gezelter Exp $,
9 > !! @version $Id: neighborLists.F90,v 1.11 2004-05-27 00:48:12 tim Exp $,
10  
11   module neighborLists
12    
# Line 44 | Line 44 | contains
44   contains
45  
46  
47 <  subroutine expandNeighborList(nGroup, error)
47 >  subroutine expandNeighborList(nGroups, error)
48      integer, intent(out) :: error
49 <    integer :: nGroup
49 >    integer :: nGroups
50      integer :: alloc_error
51      integer :: oldSize = 0
52      integer :: newSize = 0
# Line 61 | Line 61 | contains
61   #ifdef IS_MPI !! // MPI
62      if (.not. associated(point) .and. &
63           .not. associated(list) ) then
64 <       allocate(point(getNRowGroup(plan_row)),stat=alloc_error)
64 >       allocate(point(getNgroupsInRow(plan_group_row)),stat=alloc_error)
65         if (alloc_error /= 0) then
66            write(default_error,*) &
67                 "ExpandNeighborLists: Error in allocating MPI point"
68             error = -1
69            return
70         end if
71 <       allocate(list(listMultiplier * getNColGroup(plan_col)),stat=alloc_error)
71 >       allocate(list(listMultiplier * getNgroupsInCol(plan_group_col)),stat=alloc_error)
72         if (alloc_error /= 0) then
73            write(default_error,*) &
74                 "ExpandNeighborLists: Error in allocating MPI list"
# Line 82 | Line 82 | contains
82   #else !! // NONMPI
83      if (.not. associated(point) .and. &
84           .not. associated(list) ) then
85 <       allocate(point(nGroup),stat=alloc_error)
85 >       allocate(point(nGroups),stat=alloc_error)
86         if (alloc_error /= 0) then
87            write(default_error,*) &
88                 "ExpandNeighborLists: Error in allocating point"
89            error = -1
90            return
91         end if
92 <       allocate(list(listMultiplier * nGroup),stat=alloc_error)
92 >       allocate(list(listMultiplier * nGroups),stat=alloc_error)
93         if (alloc_error /= 0) then
94            write(default_error,*) &
95                 "ExpandNeighborLists: Error in allocating list"
# Line 115 | Line 115 | contains
115  
116  
117   #ifdef IS_MPI !! MPI
118 <       newSize = listMultiplier * getNCol(plan_col) + oldSize
118 >       newSize = listMultiplier * getNgroupsInCol(plan_group_col) + oldSize
119   #else
120 <       newSize = listMultiplier * nGroup + oldSize
120 >       newSize = listMultiplier * nGroups + oldSize
121   #endif !! MPI
122  
123        
# Line 151 | Line 151 | contains
151    
152    !! checks to see if any long range particle has moved
153    !! through the neighbor list skin thickness.
154 <  subroutine checkNeighborList(nGroup, q, listSkin, update_nlist)
155 <    integer :: nGroup
154 >  subroutine checkNeighborList(nGroups, q, listSkin, update_nlist)
155 >    integer :: nGroups
156      real(kind = dp), intent(in) :: listSkin
157      real( kind = dp ), dimension(:,:)  :: q
158      integer :: i
159      real( kind = DP ) :: dispmx
160      logical, intent(out) :: update_nlist
161      real( kind = DP ) :: dispmx_tmp
162    integer :: nlocal
162  
164    nlocal = nGroup
163      dispmx = 0.0E0_DP
164      !! calculate the largest displacement of any atom in any direction
165      
# Line 171 | Line 169 | contains
169         return
170      end if
171  
172 <    if (size(q0,2) /= nlocal) then
172 >    if (size(q0,2) /= nGroups) then
173         update_nlist = .true.
174         return
175      end if
# Line 180 | Line 178 | contains
178   #ifdef MPI
179      
180      dispmx_tmp = 0.0E0_DP
181 <    do i = 1, nlocal
181 >    do i = 1, nGroups
182         dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx_tmp )
183         dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx_tmp )
184         dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx_tmp )
# Line 191 | Line 189 | contains
189   #else
190      
191      dispmx = 0.0_DP
192 <    do i = 1, nlocal
192 >    do i = 1, nGroups
193         dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
194         dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
195         dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
# Line 210 | Line 208 | contains
208    !! Saves neighbor list for comparison in check.
209    !! Save_neighborList will work even if the number of
210    !! local atoms has changed.
211 <  subroutine saveNeighborList(nGroup, q)
211 >  subroutine saveNeighborList(nGroups, q)
212  
213 <    integer :: nGroup
214 <    real(kind = dp ), dimension(3,nGroup), intent(in)  :: q
213 >    integer :: nGroups
214 >    real(kind = dp ), dimension(3,nGroups), intent(in)  :: q
215      integer :: list_size
216      
217      !! get size of list
218 <    list_size = nGroup
218 >    list_size = nGroups
219      
220      if (.not. allocated(q0)) then
221         allocate(q0(3,list_size))

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines