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 388 by chuckv, Fri Mar 21 22:11:50 2003 UTC vs.
Revision 569 by mmeineke, Tue Jul 1 21:33:45 2003 UTC

# Line 6 | Line 6 | module simulation
6    use force_globals
7    use vector_class
8    use atype_module
9  use lj
9   #ifdef IS_MPI
10    use mpiSimulation
11   #endif
# Line 21 | Line 20 | module simulation
20  
21    logical, save :: simulation_setup_complete = .false.
22  
23 <  integer, public, save :: natoms
23 >  integer, public, save :: nLocal, nGlobal
24    integer, public, save :: nExcludes_Global = 0
25    integer, public, save :: nExcludes_Local = 0
26    integer, allocatable, dimension(:,:), public :: excludesLocal
27    integer, allocatable, dimension(:), public :: excludesGlobal
28 +  integer, allocatable, dimension(:), public :: molMembershipList
29  
30    real(kind=dp), save :: rcut2 = 0.0_DP
31    real(kind=dp), save :: rcut6 = 0.0_DP
32    real(kind=dp), save :: rlist2 = 0.0_DP
33 <  real(kind=dp), public, dimension(3), save :: box
34 <
33 >  real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
34 >  logical, save :: boxIsOrthorhombic
35    
36    public :: SimulationSetup
37    public :: getNlocal
38    public :: setBox
39  public :: setBox_3d
40  public :: getBox
39    public :: setRcut
40    public :: getRcut
41    public :: getRlist
# Line 54 | Line 52 | module simulation
52    public :: SimRequiresPrepairCalc
53    public :: SimRequiresPostpairCalc
54    public :: SimUsesDirectionalAtoms
57
58  interface getBox
59     module procedure getBox_3d
60     module procedure getBox_1d
61  end interface
55    
63  interface setBox
64     module procedure setBox_3d
65     module procedure setBox_1d
66  end interface
67  
56   contains
57    
58 <  subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
58 >  subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
59         CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
60 +       CmolMembership, &
61         status)    
62  
63      type (simtype) :: setThisSim
64 <    integer, intent(inout) :: nComponents
65 <    integer, dimension(nComponents),intent(inout) :: c_idents
64 >    integer, intent(inout) :: CnGlobal, CnLocal
65 >    integer, dimension(CnLocal),intent(inout) :: c_idents
66  
67      integer :: CnLocalExcludes
68      integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
69      integer :: CnGlobalExcludes
70      integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
71 +    integer, dimension(CnGlobal),intent(in) :: CmolMembership
72      !!  Result status, success = 0, status = -1
73      integer, intent(out) :: status
74      integer :: i, me, thisStat, alloc_stat, myNode
# Line 93 | Line 83 | contains
83      status = 0
84  
85      ! copy C struct into fortran type
86 +
87 +    nLocal = CnLocal
88 +    nGlobal = CnGlobal
89 +
90      thisSim = setThisSim
97    natoms = nComponents
91      rcut2 = thisSim%rcut * thisSim%rcut
92      rcut6 = rcut2 * rcut2 * rcut2
93      rlist2 = thisSim%rlist * thisSim%rlist
# Line 103 | Line 96 | contains
96      nExcludes_Global = CnGlobalExcludes
97      nExcludes_Local = CnLocalExcludes
98  
99 <    call InitializeForceGlobals(natoms, thisStat)
99 >    call InitializeForceGlobals(nLocal, thisStat)
100      if (thisStat /= 0) then
101 +       write(default_error,*) "SimSetup: InitializeForceGlobals error"
102         status = -1
103         return
104      endif
105  
106      call InitializeSimGlobals(thisStat)
107      if (thisStat /= 0) then
108 +       write(default_error,*) "SimSetup: InitializeSimGlobals error"
109         status = -1
110         return
111      endif
# Line 160 | Line 155 | contains
155      endif
156      
157   #else
158 <    do i = 1, nComponents
158 >    do i = 1, nLocal
159        
160         me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
161         atid(i) = me
# Line 168 | Line 163 | contains
163      enddo
164   #endif
165  
171    !! Create neighbor lists
172    call expandNeighborList(nComponents, thisStat)
173    if (thisStat /= 0) then
174       status = -1
175       return
176    endif
166  
167  
168      do i = 1, nExcludes_Local
# Line 184 | Line 173 | contains
173      do i = 1, nExcludes_Global
174         excludesGlobal(i) = CexcludesGlobal(i)
175      enddo
176 <    
176 >
177 >    do i = 1, nGlobal
178 >       molMemberShipList(i) = CmolMembership(i)
179 >     enddo
180 >
181      if (status == 0) simulation_setup_complete = .true.
182      
183    end subroutine SimulationSetup
184    
185 <  subroutine setBox_3d(new_box_size)
186 <    real(kind=dp), dimension(3) :: new_box_size
185 >  subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
186 >    real(kind=dp), dimension(3,3) :: cHmat, cHamtInv
187 >    integer :: cBoxIsOrthorhombic
188      integer :: smallest, status, i
189  
190 <    thisSim%box = new_box_size
191 <    box = thisSim%box
190 >    Hmat = cHmat
191 >    HmatInv = cHmatInv
192 >    if(cBoxisOrthorhombic .eq. 0 )then
193 >       boxIsOrthorhombic = .false.
194 >    else boxIsOrthorhombic = .true.
195 >    endif
196 >    
197  
199    smallest = 1
200    do i = 2, 3
201       if (new_box_size(i) .lt. new_box_size(smallest)) smallest = i
202    end do
203    if (thisSim%rcut .gt. 0.5_dp * new_box_size(smallest)) &
204         call setRcut(0.5_dp * new_box_size(smallest), status)
198      return    
199 <  end subroutine setBox_3d
199 >  end subroutine setBox
200  
208  subroutine setBox_1d(dim, new_box_size)
209    integer :: dim, status
210    real(kind=dp) :: new_box_size
211    thisSim%box(dim) = new_box_size
212    box(dim) = thisSim%box(dim)
213    if (thisSim%rcut .gt. 0.5_dp * new_box_size) &
214         call setRcut(0.5_dp * new_box_size, status)
215  end subroutine setBox_1d
216
201    subroutine setRcut(new_rcut, status)
202      real(kind = dp) :: new_rcut
203      integer :: myStatus, status
204      thisSim%rcut = new_rcut
205      rcut2 = thisSim%rcut * thisSim%rcut
206      rcut6 = rcut2 * rcut2 * rcut2
223    myStatus = 0
224    call LJ_new_rcut(new_rcut, myStatus)
225    if (myStatus .ne. 0) then
226       write(default_error, *) 'LJ module refused our rcut!'
227       status = -1
228       return
229    endif
207      status = 0
208      return
209    end subroutine setRcut
233
234  function getBox_3d() result(thisBox)
235    real( kind = dp ), dimension(3) :: thisBox
236    thisBox = thisSim%box
237  end function getBox_3d
238  
239  function getBox_1d(dim) result(thisBox)
240    integer, intent(in) :: dim
241    real( kind = dp ) :: thisBox
242    
243    thisBox = thisSim%box(dim)
244  end function getBox_1d
210      
211    subroutine getRcut(thisrcut,rc2,rc6,status)
212      real( kind = dp ), intent(out) :: thisrcut
# Line 281 | Line 246 | contains
246    function getRrf() result(rrf)
247      real( kind = dp ) :: rrf
248      rrf = thisSim%rrf
249 +    write(*,*) 'getRrf = ', rrf, thisSim%rrf
250    end function getRrf
251    
252    function getRt() result(rt)
# Line 363 | Line 329 | contains
329         thisStat = -1
330         return
331      endif
332 +
333 +    allocate(molMembershipList(nGlobal), stat=alloc_stat)
334 +    if (alloc_stat /= 0 ) then
335 +       thisStat = -1
336 +       return
337 +    endif
338      
339    end subroutine InitializeSimGlobals
340    
341    subroutine FreeSimGlobals()
342      
343      !We free in the opposite order in which we allocate in.
344 <    
344 >
345 >    if (allocated(molMembershipList)) deallocate(molMembershipList)    
346      if (allocated(excludesGlobal)) deallocate(excludesGlobal)
347      if (allocated(excludesLocal)) deallocate(excludesLocal)
348  
349    end subroutine FreeSimGlobals
350  
351 <  pure function getNlocal() result(nlocal)
352 <    integer :: nlocal
353 <    nlocal = natoms
351 >  pure function getNlocal() result(n)
352 >    integer :: n
353 >    n = nLocal
354    end function getNlocal
355  
356    

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines