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 482 by chuckv, Tue Apr 8 22:38:43 2003 UTC vs.
Revision 490 by gezelter, Fri Apr 11 15:16:59 2003 UTC

# Line 67 | Line 67 | contains
67    
68   contains
69    
70 <  subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
71 <       CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, CmolMembership, &
70 >  subroutine SimulationSetup(setThisSim, nGlobal, nLocal, c_idents, &
71 >       CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
72 >       CmolMembership, &
73         status)    
74  
75      type (simtype) :: setThisSim
76 <    integer, intent(inout) :: nComponents
77 <    integer, dimension(nComponents),intent(inout) :: c_idents
76 >    integer, intent(inout) :: nGlobal, nLocal
77 >    integer, dimension(nLocal),intent(inout) :: c_idents
78  
79      integer :: CnLocalExcludes
80      integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
81      integer :: CnGlobalExcludes
82      integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
83 <    integer, dimension(nComponents),intent(in) :: CmolMembership
83 >    integer, dimension(nGlobal),intent(in) :: CmolMembership
84      !!  Result status, success = 0, status = -1
85      integer, intent(out) :: status
86      integer :: i, me, thisStat, alloc_stat, myNode
# Line 95 | Line 96 | contains
96  
97      ! copy C struct into fortran type
98      thisSim = setThisSim
99 <    natoms = nComponents
99 >    natoms = nLocal
100      rcut2 = thisSim%rcut * thisSim%rcut
101      rcut6 = rcut2 * rcut2 * rcut2
102      rlist2 = thisSim%rlist * thisSim%rlist
# Line 163 | Line 164 | contains
164      endif
165      
166   #else
167 <    do i = 1, nComponents
167 >    do i = 1, nLocal
168        
169         me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
170         atid(i) = me
# Line 182 | Line 183 | contains
183         excludesGlobal(i) = CexcludesGlobal(i)
184      enddo
185  
186 <    molMemberShipList = CmolMembership
186 >    do i = 1, nGlobal
187 >       molMemberShipList(i) = CmolMembership(i)
188 >       ! write(0,*) 'molMembershipList(',i,') = ', molMemberShipList(i)
189 >    enddo
190  
191      if (status == 0) simulation_setup_complete = .true.
192      
# Line 360 | Line 364 | contains
364    subroutine FreeSimGlobals()
365      
366      !We free in the opposite order in which we allocate in.
367 <    
367 >
368 >    if (allocated(molMembershipList)) deallocate(molMembershipList)    
369      if (allocated(excludesGlobal)) deallocate(excludesGlobal)
370      if (allocated(excludesLocal)) deallocate(excludesLocal)
371 <    if (allocated(molMembershipList)) deallocate(molMembershipList)
371 >
372    end subroutine FreeSimGlobals
373  
374    pure function getNlocal() result(nlocal)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines