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

Comparing:
branches/mmeineke/OOPSE/libmdtools/neighborLists.F90 (file contents), Revision 377 by mmeineke, Fri Mar 21 17:42:12 2003 UTC vs.
trunk/OOPSE/libmdtools/neighborLists.F90 (file contents), Revision 845 by gezelter, Thu Oct 30 18:59:20 2003 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.1.1.1 2003-03-21 17:42:12 mmeineke Exp $,
9 > !! @version $Id: neighborLists.F90,v 1.8 2003-10-30 18:59:20 gezelter Exp $,
10  
11   module neighborLists
12    
# Line 32 | Line 32 | module neighborLists
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
35 >  real( kind = dp ), dimension(:,:), allocatable, save  :: q0
36    !! Current list size
37    integer, save :: listSize
38    !--------------MODULE ACCESS-------------------------->
# Line 52 | Line 52 | contains
52      integer :: newSize = 0
53      integer :: i
54      integer, dimension(:), pointer :: newList => null()
55 <    error = 0
56 <    
55 >    error = 0  
56  
57      !! First time through we should allocate point and list.
58      !! If one is associated and one is not, something is wrong
59      !! and return a error.
60  
61 < #ifndef IS_MPI !!/Non MPI
61 > #ifdef IS_MPI !! // MPI
62      if (.not. associated(point) .and. &
63           .not. associated(list) ) then
64 <       allocate(point(natoms),stat=alloc_error)
64 >       allocate(point(getNRow(plan_row)),stat=alloc_error)
65         if (alloc_error /= 0) then
66 <          error = -1
66 >          write(default_error,*) &
67 >               "ExpandNeighborLists: Error in allocating MPI point"
68 >           error = -1
69            return
70         end if
71 <       allocate(list(listMultiplier * natoms),stat=alloc_error)
71 >       allocate(list(listMultiplier * getNCol(plan_col)),stat=alloc_error)
72         if (alloc_error /= 0) then
73 +          write(default_error,*) &
74 +               "ExpandNeighborLists: Error in allocating MPI list"
75            error = -1
76            return
77         end if
# Line 76 | Line 79 | contains
79         nAllocations = nAllocations + 1
80         return
81      end if
82 < #else !!// MPI
82 > #else !! // NONMPI
83      if (.not. associated(point) .and. &
84           .not. associated(list) ) then
85 <       allocate(point(getNRow(plan_row)),stat=alloc_error)
85 >       allocate(point(natoms),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 * getNCol(plan_col)),stat=alloc_error)
92 >       allocate(list(listMultiplier * natoms),stat=alloc_error)
93         if (alloc_error /= 0) then
94 +          write(default_error,*) &
95 +               "ExpandNeighborLists: Error in allocating list"
96            error = -1
97            return
98         end if
# Line 99 | Line 106 | contains
106      
107      ! Check to see if we have exceeded the maximum number of allocations.
108      if (nAllocations > maxAllocations) then
109 +       write(default_error,*) &
110 +            "ExpandNeighborList: exceeded maximum number of re-allocations"
111         error = -1
112         return
113      else !! Expand the list.
114         oldSize = size(list)
115  
116 +
117   #ifdef IS_MPI !! MPI
118         newSize = listMultiplier * getNCol(plan_col) + oldSize
119   #else
# Line 141 | 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(natoms, q, rcut, rlist, update_nlist)
154 >  subroutine checkNeighborList(natoms, q, listSkin, update_nlist)
155      integer :: natoms
156 <    real(kind = dp), intent(in) :: rcut, rlist
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
152    real( kind = dp ) :: skin_thickness
162      integer :: nlocal
163 <    
163 >
164      nlocal = natoms
156    skin_thickness = rlist - rcut
165      dispmx = 0.0E0_DP
166      !! calculate the largest displacement of any atom in any direction
167      
168      !! If we have changed the particle idents, then we need to update    
169 <    if (.not. allocated(q0) .or. size(q0) /= nlocal) then
169 >    if (.not. allocated(q0)) then      
170         update_nlist = .true.
171         return
172      end if
173 <    
173 >
174 >    if (size(q0,2) /= nlocal) then
175 >       update_nlist = .true.
176 >       return
177 >    end if
178 >
179 >
180   #ifdef MPI
181      
182      dispmx_tmp = 0.0E0_DP
# Line 183 | Line 197 | contains
197         dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
198      end do
199  
200 < #endif
201 <    
200 > #endif  
201 >
202      !! a conservative test of list skin crossings
203 <    dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
203 >    dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)  
204 >
205 >    update_nlist = (dispmx.gt.listSkin)
206      
207 <    update_nlist = (dispmx.gt.(skin_thickness))
192 <    
193 <  end subroutine checkNeighborList
207 >   end subroutine checkNeighborList
208    
209    
210    !! Saves neighbor list for comparison in check.
211    !! Save_neighborList will work even if the number of
212    !! local atoms has changed.
213 <  subroutine saveNeighborList(q)
214 <    real(kind = dp ), dimension(:,:), intent(in)  :: q
213 >  subroutine saveNeighborList(natoms, q)
214 >
215 >    integer :: natoms
216 >    real(kind = dp ), dimension(3,natoms), intent(in)  :: q
217      integer :: list_size
218      
219      !! get size of list
220 <    list_size = size(q)
220 >    list_size = natoms
221      
222      if (.not. allocated(q0)) then
223         allocate(q0(3,list_size))
224 <    else if( list_size > size(q0)) then
224 >    else if( list_size > size(q0,2)) then
225         deallocate(q0)
226         allocate(q0(3,list_size))
227      endif

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines