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 747 by gezelter, Fri Sep 5 21:28:52 2003 UTC

# Line 16 | Line 16 | module simulation
16   #define __FORTRAN90
17   #include "fSimulation.h"
18  
19 <  type (simtype), public :: thisSim
19 >  type (simtype), public, save :: thisSim
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 <
30 >  real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
31 >  logical, public, save :: boxIsOrthorhombic
32    
33    public :: SimulationSetup
34    public :: getNlocal
35    public :: setBox
39  public :: setBox_3d
40  public :: getBox
41  public :: setRcut
42  public :: getRcut
43  public :: getRlist
44  public :: getRrf
45  public :: getRt
36    public :: getDielect
37    public :: SimUsesPBC
38    public :: SimUsesLJ
# Line 54 | Line 44 | module simulation
44    public :: SimRequiresPrepairCalc
45    public :: SimRequiresPostpairCalc
46    public :: SimUsesDirectionalAtoms
57
58  interface getBox
59     module procedure getBox_3d
60     module procedure getBox_1d
61  end interface
47    
63  interface setBox
64     module procedure setBox_3d
65     module procedure setBox_1d
66  end interface
67  
48   contains
49    
50 <  subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
51 <       CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, CmolMembership, &
50 >  subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
51 >       CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
52 >       CmolMembership, &
53         status)    
54  
55      type (simtype) :: setThisSim
56 <    integer, intent(inout) :: nComponents
57 <    integer, dimension(nComponents),intent(inout) :: c_idents
56 >    integer, intent(inout) :: CnGlobal, CnLocal
57 >    integer, dimension(CnLocal),intent(inout) :: c_idents
58  
59      integer :: CnLocalExcludes
60      integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
61      integer :: CnGlobalExcludes
62      integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
63 <    integer, dimension(nComponents),intent(in) :: CmolMembership
63 >    integer, dimension(CnGlobal),intent(in) :: CmolMembership
64      !!  Result status, success = 0, status = -1
65      integer, intent(out) :: status
66      integer :: i, me, thisStat, alloc_stat, myNode
# Line 94 | Line 75 | contains
75      status = 0
76  
77      ! copy C struct into fortran type
78 +
79 +    nLocal = CnLocal
80 +    nGlobal = CnGlobal
81 +
82      thisSim = setThisSim
98    natoms = nComponents
99    rcut2 = thisSim%rcut * thisSim%rcut
100    rcut6 = rcut2 * rcut2 * rcut2
101    rlist2 = thisSim%rlist * thisSim%rlist
102    box = thisSim%box
83  
84      nExcludes_Global = CnGlobalExcludes
85      nExcludes_Local = CnLocalExcludes
86  
87 <    call InitializeForceGlobals(natoms, thisStat)
87 >    call InitializeForceGlobals(nLocal, thisStat)
88      if (thisStat /= 0) then
89         write(default_error,*) "SimSetup: InitializeForceGlobals error"
90         status = -1
# Line 162 | Line 142 | contains
142         deallocate(c_idents_Row)
143      endif
144      
145 < #else
146 <    do i = 1, nComponents
145 > #endif
146 >
147 > ! We build the local atid's for both mpi and nonmpi
148 >    do i = 1, nLocal
149        
150         me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
151         atid(i) = me
152    
153      enddo
172 #endif
154  
155  
156  
157 +
158      do i = 1, nExcludes_Local
159         excludesLocal(1,i) = CexcludesLocal(1,i)
160         excludesLocal(2,i) = CexcludesLocal(2,i)
# Line 182 | Line 164 | contains
164         excludesGlobal(i) = CexcludesGlobal(i)
165      enddo
166  
167 <    molMemberShipList = CmolMembership
167 >    do i = 1, nGlobal
168 >       molMemberShipList(i) = CmolMembership(i)
169 >     enddo
170  
171      if (status == 0) simulation_setup_complete = .true.
172      
173    end subroutine SimulationSetup
174    
175 <  subroutine setBox_3d(new_box_size)
176 <    real(kind=dp), dimension(3) :: new_box_size
175 >  subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
176 >    real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
177 >    integer :: cBoxIsOrthorhombic
178      integer :: smallest, status, i
179 <
180 <    thisSim%box = new_box_size
181 <    box = thisSim%box
182 <
179 >    
180 >    Hmat = cHmat
181 >    HmatInv = cHmatInv
182 >    if (cBoxIsOrthorhombic .eq. 0 ) then
183 >       boxIsOrthorhombic = .false.
184 >    else
185 >       boxIsOrthorhombic = .true.
186 >    endif
187 >    
188      return    
189 <  end subroutine setBox_3d
200 <
201 <  subroutine setBox_1d(dim, new_box_size)
202 <    integer :: dim, status
203 <    real(kind=dp) :: new_box_size
204 <    thisSim%box(dim) = new_box_size
205 <    box(dim) = thisSim%box(dim)
206 <  end subroutine setBox_1d
189 >  end subroutine setBox
190  
208  subroutine setRcut(new_rcut, status)
209    real(kind = dp) :: new_rcut
210    integer :: myStatus, status
211    thisSim%rcut = new_rcut
212    rcut2 = thisSim%rcut * thisSim%rcut
213    rcut6 = rcut2 * rcut2 * rcut2
214    status = 0
215    return
216  end subroutine setRcut
217
218  function getBox_3d() result(thisBox)
219    real( kind = dp ), dimension(3) :: thisBox
220    thisBox = thisSim%box
221  end function getBox_3d
222  
223  function getBox_1d(dim) result(thisBox)
224    integer, intent(in) :: dim
225    real( kind = dp ) :: thisBox
226    
227    thisBox = thisSim%box(dim)
228  end function getBox_1d
229    
230  subroutine getRcut(thisrcut,rc2,rc6,status)
231    real( kind = dp ), intent(out) :: thisrcut
232    real( kind = dp ), intent(out), optional :: rc2
233    real( kind = dp ), intent(out), optional :: rc6
234    integer, optional :: status
235
236    if (present(status)) status = 0
237    
238    if (.not.simulation_setup_complete ) then
239       if (present(status)) status = -1
240       return
241    end if
242    
243    thisrcut = thisSim%rcut
244    if(present(rc2)) rc2 = rcut2
245    if(present(rc6)) rc6 = rcut6
246  end subroutine getRcut
247  
248  subroutine getRlist(thisrlist,rl2,status)
249    real( kind = dp ), intent(out) :: thisrlist
250    real( kind = dp ), intent(out), optional :: rl2
251
252    integer, optional :: status
253
254    if (present(status)) status = 0
255
256    if (.not.simulation_setup_complete ) then
257       if (present(status)) status = -1
258       return
259    end if
260    
261    thisrlist = thisSim%rlist
262    if(present(rl2)) rl2 = rlist2
263  end subroutine getRlist
264
265  function getRrf() result(rrf)
266    real( kind = dp ) :: rrf
267    rrf = thisSim%rrf
268    write(*,*) 'getRrf = ', rrf, thisSim%rrf
269  end function getRrf
270  
271  function getRt() result(rt)
272    real( kind = dp ) :: rt
273    rt = thisSim%rt
274  end function getRt
275
191    function getDielect() result(dielect)
192      real( kind = dp ) :: dielect
193      dielect = thisSim%dielect
# Line 349 | Line 264 | contains
264         return
265      endif
266  
267 <    allocate(molMembershipList(getNlocal()), stat=alloc_stat)
267 >    allocate(molMembershipList(nGlobal), stat=alloc_stat)
268      if (alloc_stat /= 0 ) then
269         thisStat = -1
270         return
# Line 360 | Line 275 | contains
275    subroutine FreeSimGlobals()
276      
277      !We free in the opposite order in which we allocate in.
278 <    
278 >
279 >    if (allocated(molMembershipList)) deallocate(molMembershipList)    
280      if (allocated(excludesGlobal)) deallocate(excludesGlobal)
281      if (allocated(excludesLocal)) deallocate(excludesLocal)
282 <    if (allocated(molMembershipList)) deallocate(molMembershipList)
282 >
283    end subroutine FreeSimGlobals
284  
285 <  pure function getNlocal() result(nlocal)
286 <    integer :: nlocal
287 <    nlocal = natoms
285 >  pure function getNlocal() result(n)
286 >    integer :: n
287 >    n = nLocal
288    end function getNlocal
289  
290    

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines