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 1144 by tim, Sat May 1 18:52:38 2004 UTC vs.
Revision 1150 by gezelter, Fri May 7 21:35:05 2004 UTC

# Line 6 | Line 6 | module simulation
6    use force_globals
7    use vector_class
8    use atype_module
9 +  use switcheroo
10   #ifdef IS_MPI
11    use mpiSimulation
12   #endif
# Line 15 | Line 16 | module simulation
16  
17   #define __FORTRAN90
18   #include "fSimulation.h"
19 + #include "fSwitchingFunction.h"
20  
21    type (simtype), public, save :: thisSim
22  
23    logical, save :: simulation_setup_complete = .false.
24  
25    integer, public, save :: nLocal, nGlobal
26 +  integer, public, save :: nGroup, nGroupGlobal
27    integer, public, save :: nExcludes_Global = 0
28    integer, public, save :: nExcludes_Local = 0
29    integer, allocatable, dimension(:,:), public :: excludesLocal
30    integer, allocatable, dimension(:), public :: excludesGlobal
31    integer, allocatable, dimension(:), public :: molMembershipList
32 +  integer, allocatable, dimension(:), public :: groupList
33 +  integer, allocatable, dimension(:), public :: groupStart
34 +  real(kind=dp), allocatable, dimension(:), public :: mfact
35  
36    real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
37    logical, public, save :: boxIsOrthorhombic
# Line 50 | Line 56 | contains
56    
57    subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
58         CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
59 <       CmolMembership, mfact, ngroup, groupList, groupStart, &
59 >       CmolMembership, Cmfact, CnGroup, CgroupList, CgroupStart, &
60         status)    
61  
62      type (simtype) :: setThisSim
# Line 65 | Line 71 | contains
71      !!  Result status, success = 0, status = -1
72      integer, intent(out) :: status
73      integer :: i, me, thisStat, alloc_stat, myNode
74 +    integer :: gStart, gEnd, groupSize, biggestGroupSize, ia
75  
76      !! mass factors used for molecular cutoffs
77 <    real ( kind = dp ), dimension(3,nLocal) :: mfact
78 <    integer, intent(in):: ngroup
79 <    integer, dimension(nLocal),intent(in) :: groupList
80 <    integer, dimension(ngroup),intent(in) :: groupStart
77 >    real ( kind = dp ), dimension(CnLocal) :: Cmfact
78 >    integer, intent(in):: CnGroup
79 >    integer, dimension(CnLocal),intent(in) :: CgroupList
80 >    integer, dimension(CnGroup),intent(in) :: CgroupStart
81  
82   #ifdef IS_MPI
83      integer, allocatable, dimension(:) :: c_idents_Row
# Line 86 | Line 93 | contains
93  
94      nLocal = CnLocal
95      nGlobal = CnGlobal
96 +    nGroup = CnGroup
97  
98      thisSim = setThisSim
99  
# Line 160 | Line 168 | contains
168    
169      enddo
170  
163
164
165
171      do i = 1, nExcludes_Local
172         excludesLocal(1,i) = CexcludesLocal(1,i)
173         excludesLocal(2,i) = CexcludesLocal(2,i)
# Line 174 | Line 179 | contains
179  
180      do i = 1, nGlobal
181         molMemberShipList(i) = CmolMembership(i)
182 <     enddo
182 >    enddo
183 >    
184 >    biggestGroupSize = 0
185 >    do i = 1, nGroup
186 >       groupStart(i) = CgroupStart(i)
187 >       groupSize = 0
188 >       gStart = groupStart(i)      
189 >       if (i .eq. nGroup) then
190 >          gEnd = nLocal
191 >       else
192 >          gEnd = CgroupStart(i+1) - 1
193 >       endif
194 >       do ia = gStart, gEnd
195 >          groupList(ia) = CgroupList(ia)
196 >          groupSize = groupSize + 1
197 >       enddo
198 >       if (groupSize .gt. biggestGroupSize) biggestGroupSize = groupSize      
199 >    enddo
200 >    groupStart(nGroup+1) = nLocal+1
201  
202 +    do i = 1, nLocal
203 +       mfact(i) = Cmfact(i)
204 +    enddo    
205 +    
206      if (status == 0) simulation_setup_complete = .true.
207      
208    end subroutine SimulationSetup
# Line 282 | Line 309 | contains
309         thisStat = -1
310         return
311      endif
312 +
313 +    allocate(groupStart(nGroup+1), stat=alloc_stat)
314 +    if (alloc_stat /= 0 ) then
315 +       thisStat = -1
316 +       return
317 +    endif
318 +
319 +    allocate(groupList(nLocal), stat=alloc_stat)
320 +    if (alloc_stat /= 0 ) then
321 +       thisStat = -1
322 +       return
323 +    endif
324 +
325 +    allocate(mfact(nLocal), stat=alloc_stat)
326 +    if (alloc_stat /= 0 ) then
327 +       thisStat = -1
328 +       return
329 +    endif
330      
331    end subroutine InitializeSimGlobals
332    
# Line 289 | Line 334 | contains
334      
335      !We free in the opposite order in which we allocate in.
336  
337 +    if (allocated(mfact)) deallocate(mfact)
338 +    if (allocated(groupList)) deallocate(groupList)    
339 +    if (allocated(groupStart)) deallocate(groupStart)    
340      if (allocated(molMembershipList)) deallocate(molMembershipList)    
341      if (allocated(excludesGlobal)) deallocate(excludesGlobal)
342      if (allocated(excludesLocal)) deallocate(excludesLocal)
343 <
343 >    
344    end subroutine FreeSimGlobals
345 <
345 >  
346    pure function getNlocal() result(n)
347      integer :: n
348      n = nLocal
349    end function getNlocal
302
350    
351 +  
352   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines