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

Comparing trunk/OOPSE/libmdtools/simulation_module.F90 (file contents):
Revision 1206 by tim, Thu May 27 19:51:18 2004 UTC vs.
Revision 1217 by gezelter, Tue Jun 1 21:45:22 2004 UTC

# Line 64 | Line 64 | contains
64    
65    subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
66         CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
67 <       CmolMembership, Cmfact, CnGroups, CgroupList, CgroupStart, &
67 >       CmolMembership, Cmfact, CnGroups, CglobalGroupMembership, &
68         status)    
69  
70      type (simtype) :: setThisSim
# Line 84 | Line 84 | contains
84      !! mass factors used for molecular cutoffs
85      real ( kind = dp ), dimension(CnLocal) :: Cmfact
86      integer, intent(in):: CnGroups
87 <    integer, dimension(CnLocal),intent(in) :: CgroupList
88 <    integer, dimension(CnGroups),intent(in) :: CgroupStart
89 <    integer :: maxSkipsForAtom
87 >    integer, dimension(CnGlobal), intent(in):: CglobalGroupMembership
88 >    integer :: maxSkipsForAtom, glPointer
89  
90   #ifdef IS_MPI
91      integer, allocatable, dimension(:) :: c_idents_Row
92      integer, allocatable, dimension(:) :: c_idents_Col
93 <    integer :: nAtomsInRow, nGroupsInRow
94 <    integer :: nAtomsInCol, nGroupsInCol
93 >    integer :: nAtomsInRow, nGroupsInRow, aid
94 >    integer :: nAtomsInCol, nGroupsInCol, gid
95   #endif  
96  
97      simulation_setup_complete = .false.
# Line 202 | Line 201 | contains
201         status = -1
202         return
203      endif
204 <    allocate(groupStartLocal(nGroups), stat=alloc_stat)
204 >    allocate(mfactLocal(nLocal),stat=alloc_stat)
205      if (alloc_stat /= 0 ) then
206         status = -1
207         return
208      endif
209 <    allocate(groupListLocal(nLocal), stat=alloc_stat)
210 <    if (alloc_stat /= 0 ) then
212 <       status = -1
213 <       return
214 <    endif    
215 <    allocate(mfactLocal(nLocal), stat=alloc_stat)
216 <    if (alloc_stat /= 0 ) then
217 <       status = -1
218 <       return
219 <    endif    
220 <            
221 <    groupStartLocal = CgroupStart
222 <    groupListLocal = CgroupList
223 <    mfactLocal = Cmfact        
224 <            
225 <    call gather(groupStartLocal, groupStartRow, plan_group_row)
226 <    call gather(groupStartLocal, groupStartCol, plan_group_col)  
227 <    groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1
228 <    groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1
229 <    call gather(groupListLocal,  groupListRow,  plan_atom_row)
230 <    call gather(groupListLocal,  groupListCol,  plan_atom_col)
209 >    
210 >    glPointer = 1
211  
212 <    ! C passes us groupList as globalID numbers for the atoms
233 <    ! we need to remap these onto the row and column ids on this
234 <    ! processor.  This is a linear search, but is only done once
235 <    ! (we hope)
212 >    do i = 1, nGroupsInRow
213  
214 <    do i = 1, nAtomsInRow
214 >       gid = GroupRowToGlobal(i)
215 >       groupStartRow(i) = glPointer      
216 >
217         do j = 1, nAtomsInRow
218 <          if (AtomRowToGlobal(j) .eq. groupListRow(i)) then
219 <             groupListRow(i) = j
218 >          aid = AtomRowToGlobal(j)
219 >          if (CglobalGroupMembership(aid) .eq. gid) then
220 >             groupListRow(glPointer) = j
221 >             glPointer = glPointer + 1
222            endif
223         enddo
224      enddo
225 <    do i = 1, nAtomsInCol
225 >    groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1
226 >
227 >    glPointer = 1
228 >
229 >    do i = 1, nGroupsInCol
230 >
231 >       gid = GroupColToGlobal(i)
232 >       groupStartCol(i) = glPointer      
233 >
234         do j = 1, nAtomsInCol
235 <          if (AtomColToGlobal(j) .eq. groupListCol(i)) then
236 <             groupListCol(i) = j
235 >          aid = AtomColToGlobal(j)
236 >          if (CglobalGroupMembership(aid) .eq. gid) then
237 >             groupListCol(glPointer) = j
238 >             glPointer = glPointer + 1
239            endif
240         enddo
241      enddo
242 +    groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1
243 +
244 +    mfactLocal = Cmfact        
245 +
246      call gather(mfactLocal,      mfactRow,      plan_atom_row)
247      call gather(mfactLocal,      mfactCol,      plan_atom_col)
248      
249      if (allocated(mfactLocal)) then
250         deallocate(mfactLocal)
251      end if
257    if (allocated(groupListLocal)) then
258       deallocate(groupListLocal)
259    endif
260    if (allocated(groupStartLocal)) then
261       deallocate(groupStartLocal)
262    endif
252   #else
253      allocate(groupStartRow(nGroups+1),stat=alloc_stat)
254      if (alloc_stat /= 0 ) then
# Line 291 | Line 280 | contains
280         status = -1
281         return
282      endif
283 +    allocate(mfactLocal(nLocal),stat=alloc_stat)
284 +    if (alloc_stat /= 0 ) then
285 +       status = -1
286 +       return
287 +    endif
288 +
289 +    glPointer = 1
290      do i = 1, nGroups
291 <       groupStartRow(i) = CgroupStart(i)
292 <       groupStartCol(i) = CgroupStart(i)
293 <    end do
291 >       groupStartRow(i) = glPointer      
292 >       groupStartCol(i) = glPointer
293 >       do j = 1, nLocal
294 >          if (CglobalGroupMembership(j) .eq. i) then
295 >             groupListRow(glPointer) = j
296 >             groupListCol(glPointer) = j
297 >             glPointer = glPointer + 1
298 >          endif
299 >       enddo
300 >    enddo
301      groupStartRow(nGroups+1) = nLocal + 1
302      groupStartCol(nGroups+1) = nLocal + 1
303 +
304      do i = 1, nLocal
301       groupListRow(i) = CgroupList(i)
302       groupListCol(i) = CgroupList(i)
305         mfactRow(i) = Cmfact(i)
306         mfactCol(i) = Cmfact(i)
307      end do
# Line 320 | Line 322 | contains
322         excludesLocal(2,i) = CexcludesLocal(2,i)
323      enddo
324  
325 + #ifdef IS_MPI
326 +    allocate(nSkipsForAtom(nAtomsInRow), stat=alloc_stat)
327 + #else
328 +    allocate(nSkipsForAtom(nLocal), stat=alloc_stat)
329 + #endif
330 +    if (alloc_stat /= 0 ) then
331 +       thisStat = -1
332 +       write(*,*) 'Could not allocate nSkipsForAtom array'
333 +       return
334 +    endif
335 +
336      maxSkipsForAtom = 0
337   #ifdef IS_MPI
338      do j = 1, nAtomsInRow
# Line 379 | Line 392 | contains
392               id2 = excludesLocal(2,i)
393               skipsForAtom(j, nSkipsForAtom(j)) = id2
394            endif
395 <          if (excludesLocal(2, i) .eq. id2 ) then
395 >          if (excludesLocal(2, i) .eq. id1 ) then
396               nSkipsForAtom(j) = nSkipsForAtom(j) + 1
397               ! exclude lists have global ID's so this line is
398               ! the same in MPI and non-MPI
# Line 485 | Line 498 | contains
498      thisStat = 0
499      
500      call FreeSimGlobals()    
488
489    allocate(nSkipsForAtom(nLocal), stat=alloc_stat)
490    if (alloc_stat /= 0 ) then
491       thisStat = -1
492       return
493    endif
501      
502      allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
503      if (alloc_stat /= 0 ) then
# Line 517 | Line 524 | contains
524      !We free in the opposite order in which we allocate in.
525  
526      if (allocated(skipsForAtom)) deallocate(skipsForAtom)
527 +    if (allocated(nSkipsForAtom)) deallocate(nSkipsForAtom)
528 +    if (allocated(mfactLocal)) deallocate(mfactLocal)
529      if (allocated(mfactCol)) deallocate(mfactCol)
530      if (allocated(mfactRow)) deallocate(mfactRow)
531      if (allocated(groupListCol)) deallocate(groupListCol)    

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines