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 480 by chuckv, Tue Apr 8 17:16:22 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
31 <  real(kind=dp), save :: rlist2 = 0.0_DP
32 <  real(kind=dp), public, dimension(3), save :: box
33 <
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
38  public :: setBox_3d
39  public :: getBox
40  public :: setRcut
41  public :: getRcut
42  public :: getRlist
43  public :: getRrf
44  public :: getRt
36    public :: getDielect
37    public :: SimUsesPBC
38    public :: SimUsesLJ
# Line 53 | Line 44 | module simulation
44    public :: SimRequiresPrepairCalc
45    public :: SimRequiresPostpairCalc
46    public :: SimUsesDirectionalAtoms
56
57  interface getBox
58     module procedure getBox_3d
59     module procedure getBox_1d
60  end interface
47    
62  interface setBox
63     module procedure setBox_3d
64     module procedure setBox_1d
65  end interface
66  
48   contains
49    
50 <  subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
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(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 92 | 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
96    natoms = nComponents
97    rcut2 = thisSim%rcut * thisSim%rcut
98    rcut6 = rcut2 * rcut2 * rcut2
99    rlist2 = thisSim%rlist * thisSim%rlist
100    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 160 | 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
170 #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 180 | Line 164 | contains
164         excludesGlobal(i) = CexcludesGlobal(i)
165      enddo
166  
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
196 <
197 <  subroutine setBox_1d(dim, new_box_size)
198 <    integer :: dim, status
199 <    real(kind=dp) :: new_box_size
200 <    thisSim%box(dim) = new_box_size
201 <    box(dim) = thisSim%box(dim)
202 <  end subroutine setBox_1d
189 >  end subroutine setBox
190  
204  subroutine setRcut(new_rcut, status)
205    real(kind = dp) :: new_rcut
206    integer :: myStatus, status
207    thisSim%rcut = new_rcut
208    rcut2 = thisSim%rcut * thisSim%rcut
209    rcut6 = rcut2 * rcut2 * rcut2
210    status = 0
211    return
212  end subroutine setRcut
213
214  function getBox_3d() result(thisBox)
215    real( kind = dp ), dimension(3) :: thisBox
216    thisBox = thisSim%box
217  end function getBox_3d
218  
219  function getBox_1d(dim) result(thisBox)
220    integer, intent(in) :: dim
221    real( kind = dp ) :: thisBox
222    
223    thisBox = thisSim%box(dim)
224  end function getBox_1d
225    
226  subroutine getRcut(thisrcut,rc2,rc6,status)
227    real( kind = dp ), intent(out) :: thisrcut
228    real( kind = dp ), intent(out), optional :: rc2
229    real( kind = dp ), intent(out), optional :: rc6
230    integer, optional :: status
231
232    if (present(status)) status = 0
233    
234    if (.not.simulation_setup_complete ) then
235       if (present(status)) status = -1
236       return
237    end if
238    
239    thisrcut = thisSim%rcut
240    if(present(rc2)) rc2 = rcut2
241    if(present(rc6)) rc6 = rcut6
242  end subroutine getRcut
243  
244  subroutine getRlist(thisrlist,rl2,status)
245    real( kind = dp ), intent(out) :: thisrlist
246    real( kind = dp ), intent(out), optional :: rl2
247
248    integer, optional :: status
249
250    if (present(status)) status = 0
251
252    if (.not.simulation_setup_complete ) then
253       if (present(status)) status = -1
254       return
255    end if
256    
257    thisrlist = thisSim%rlist
258    if(present(rl2)) rl2 = rlist2
259  end subroutine getRlist
260
261  function getRrf() result(rrf)
262    real( kind = dp ) :: rrf
263    rrf = thisSim%rrf
264    write(*,*) 'getRrf = ', rrf, thisSim%rrf
265  end function getRrf
266  
267  function getRt() result(rt)
268    real( kind = dp ) :: rt
269    rt = thisSim%rt
270  end function getRt
271
191    function getDielect() result(dielect)
192      real( kind = dp ) :: dielect
193      dielect = thisSim%dielect
# Line 344 | Line 263 | contains
263         thisStat = -1
264         return
265      endif
266 +
267 +    allocate(molMembershipList(nGlobal), stat=alloc_stat)
268 +    if (alloc_stat /= 0 ) then
269 +       thisStat = -1
270 +       return
271 +    endif
272      
273    end subroutine InitializeSimGlobals
274    
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  
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