ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 1150
Committed: Fri May 7 21:35:05 2004 UTC (20 years, 2 months ago) by gezelter
File size: 9365 byte(s)
Log Message:
Many changes to get group-based cutoffs to work

File Contents

# User Rev Content
1 mmeineke 377 !! Fortran interface to C entry plug.
2    
3     module simulation
4     use definitions
5     use neighborLists
6     use force_globals
7     use vector_class
8     use atype_module
9 gezelter 1150 use switcheroo
10 mmeineke 377 #ifdef IS_MPI
11     use mpiSimulation
12     #endif
13    
14     implicit none
15     PRIVATE
16    
17     #define __FORTRAN90
18     #include "fSimulation.h"
19 gezelter 1150 #include "fSwitchingFunction.h"
20 mmeineke 377
21 gezelter 747 type (simtype), public, save :: thisSim
22 mmeineke 377
23     logical, save :: simulation_setup_complete = .false.
24    
25 mmeineke 491 integer, public, save :: nLocal, nGlobal
26 gezelter 1150 integer, public, save :: nGroup, nGroupGlobal
27 mmeineke 377 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 chuckv 482 integer, allocatable, dimension(:), public :: molMembershipList
32 gezelter 1150 integer, allocatable, dimension(:), public :: groupList
33     integer, allocatable, dimension(:), public :: groupStart
34     real(kind=dp), allocatable, dimension(:), public :: mfact
35 mmeineke 377
36 mmeineke 569 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
37 gezelter 570 logical, public, save :: boxIsOrthorhombic
38 mmeineke 377
39     public :: SimulationSetup
40     public :: getNlocal
41     public :: setBox
42     public :: getDielect
43     public :: SimUsesPBC
44     public :: SimUsesLJ
45 gezelter 941 public :: SimUsesCharges
46 mmeineke 377 public :: SimUsesDipoles
47     public :: SimUsesSticky
48     public :: SimUsesRF
49     public :: SimUsesGB
50     public :: SimUsesEAM
51     public :: SimRequiresPrepairCalc
52     public :: SimRequiresPostpairCalc
53     public :: SimUsesDirectionalAtoms
54    
55     contains
56    
57 mmeineke 491 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
58 gezelter 483 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
59 gezelter 1150 CmolMembership, Cmfact, CnGroup, CgroupList, CgroupStart, &
60 mmeineke 377 status)
61    
62     type (simtype) :: setThisSim
63 mmeineke 491 integer, intent(inout) :: CnGlobal, CnLocal
64     integer, dimension(CnLocal),intent(inout) :: c_idents
65 mmeineke 377
66     integer :: CnLocalExcludes
67     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
68     integer :: CnGlobalExcludes
69     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
70 mmeineke 491 integer, dimension(CnGlobal),intent(in) :: CmolMembership
71 mmeineke 377 !! Result status, success = 0, status = -1
72     integer, intent(out) :: status
73     integer :: i, me, thisStat, alloc_stat, myNode
74 gezelter 1150 integer :: gStart, gEnd, groupSize, biggestGroupSize, ia
75 tim 1144
76     !! mass factors used for molecular cutoffs
77 gezelter 1150 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 tim 1144
82 mmeineke 377 #ifdef IS_MPI
83     integer, allocatable, dimension(:) :: c_idents_Row
84     integer, allocatable, dimension(:) :: c_idents_Col
85     integer :: nrow
86     integer :: ncol
87     #endif
88    
89     simulation_setup_complete = .false.
90     status = 0
91    
92     ! copy C struct into fortran type
93 mmeineke 491
94     nLocal = CnLocal
95     nGlobal = CnGlobal
96 gezelter 1150 nGroup = CnGroup
97 mmeineke 491
98 mmeineke 377 thisSim = setThisSim
99 mmeineke 620
100 mmeineke 377 nExcludes_Global = CnGlobalExcludes
101     nExcludes_Local = CnLocalExcludes
102    
103 mmeineke 491 call InitializeForceGlobals(nLocal, thisStat)
104 mmeineke 377 if (thisStat /= 0) then
105 chuckv 480 write(default_error,*) "SimSetup: InitializeForceGlobals error"
106 mmeineke 377 status = -1
107     return
108     endif
109    
110     call InitializeSimGlobals(thisStat)
111     if (thisStat /= 0) then
112 chuckv 480 write(default_error,*) "SimSetup: InitializeSimGlobals error"
113 mmeineke 377 status = -1
114     return
115     endif
116    
117     #ifdef IS_MPI
118     ! We can only set up forces if mpiSimulation has been setup.
119     if (.not. isMPISimSet()) then
120     write(default_error,*) "MPI is not set"
121     status = -1
122     return
123     endif
124     nrow = getNrow(plan_row)
125     ncol = getNcol(plan_col)
126     mynode = getMyNode()
127    
128     allocate(c_idents_Row(nrow),stat=alloc_stat)
129     if (alloc_stat /= 0 ) then
130     status = -1
131     return
132     endif
133    
134     allocate(c_idents_Col(ncol),stat=alloc_stat)
135     if (alloc_stat /= 0 ) then
136     status = -1
137     return
138     endif
139    
140     call gather(c_idents, c_idents_Row, plan_row)
141     call gather(c_idents, c_idents_Col, plan_col)
142    
143     do i = 1, nrow
144     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
145     atid_Row(i) = me
146     enddo
147    
148     do i = 1, ncol
149     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
150     atid_Col(i) = me
151     enddo
152    
153     !! free temporary ident arrays
154     if (allocated(c_idents_Col)) then
155     deallocate(c_idents_Col)
156     end if
157     if (allocated(c_idents_Row)) then
158     deallocate(c_idents_Row)
159     endif
160    
161 chuckv 648 #endif
162    
163     ! We build the local atid's for both mpi and nonmpi
164 gezelter 490 do i = 1, nLocal
165 mmeineke 377
166     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
167     atid(i) = me
168    
169     enddo
170    
171     do i = 1, nExcludes_Local
172     excludesLocal(1,i) = CexcludesLocal(1,i)
173     excludesLocal(2,i) = CexcludesLocal(2,i)
174     enddo
175    
176     do i = 1, nExcludes_Global
177     excludesGlobal(i) = CexcludesGlobal(i)
178     enddo
179 mmeineke 435
180 gezelter 490 do i = 1, nGlobal
181 gezelter 483 molMemberShipList(i) = CmolMembership(i)
182 gezelter 1150 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 chuckv 482
202 gezelter 1150 do i = 1, nLocal
203     mfact(i) = Cmfact(i)
204     enddo
205    
206 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
207    
208     end subroutine SimulationSetup
209    
210 mmeineke 569 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
211 gezelter 570 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
212 mmeineke 569 integer :: cBoxIsOrthorhombic
213 mmeineke 377 integer :: smallest, status, i
214 gezelter 570
215 mmeineke 569 Hmat = cHmat
216     HmatInv = cHmatInv
217 gezelter 570 if (cBoxIsOrthorhombic .eq. 0 ) then
218 mmeineke 569 boxIsOrthorhombic = .false.
219 gezelter 570 else
220     boxIsOrthorhombic = .true.
221 mmeineke 569 endif
222    
223 mmeineke 377 return
224 mmeineke 569 end subroutine setBox
225 mmeineke 377
226     function getDielect() result(dielect)
227     real( kind = dp ) :: dielect
228     dielect = thisSim%dielect
229     end function getDielect
230    
231     function SimUsesPBC() result(doesit)
232     logical :: doesit
233     doesit = thisSim%SIM_uses_PBC
234     end function SimUsesPBC
235    
236     function SimUsesLJ() result(doesit)
237     logical :: doesit
238     doesit = thisSim%SIM_uses_LJ
239     end function SimUsesLJ
240    
241     function SimUsesSticky() result(doesit)
242     logical :: doesit
243     doesit = thisSim%SIM_uses_sticky
244     end function SimUsesSticky
245    
246 gezelter 941 function SimUsesCharges() result(doesit)
247     logical :: doesit
248     doesit = thisSim%SIM_uses_charges
249     end function SimUsesCharges
250    
251 mmeineke 377 function SimUsesDipoles() result(doesit)
252     logical :: doesit
253     doesit = thisSim%SIM_uses_dipoles
254     end function SimUsesDipoles
255    
256     function SimUsesRF() result(doesit)
257     logical :: doesit
258     doesit = thisSim%SIM_uses_RF
259     end function SimUsesRF
260    
261     function SimUsesGB() result(doesit)
262     logical :: doesit
263     doesit = thisSim%SIM_uses_GB
264     end function SimUsesGB
265    
266     function SimUsesEAM() result(doesit)
267     logical :: doesit
268     doesit = thisSim%SIM_uses_EAM
269     end function SimUsesEAM
270    
271     function SimUsesDirectionalAtoms() result(doesit)
272     logical :: doesit
273     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
274     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
275     end function SimUsesDirectionalAtoms
276    
277     function SimRequiresPrepairCalc() result(doesit)
278     logical :: doesit
279     doesit = thisSim%SIM_uses_EAM
280     end function SimRequiresPrepairCalc
281    
282     function SimRequiresPostpairCalc() result(doesit)
283     logical :: doesit
284     doesit = thisSim%SIM_uses_RF
285     end function SimRequiresPostpairCalc
286    
287     subroutine InitializeSimGlobals(thisStat)
288     integer, intent(out) :: thisStat
289     integer :: alloc_stat
290    
291     thisStat = 0
292    
293     call FreeSimGlobals()
294    
295     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
296     if (alloc_stat /= 0 ) then
297     thisStat = -1
298     return
299     endif
300    
301     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
302     if (alloc_stat /= 0 ) then
303     thisStat = -1
304     return
305     endif
306 chuckv 482
307 mmeineke 491 allocate(molMembershipList(nGlobal), stat=alloc_stat)
308 chuckv 482 if (alloc_stat /= 0 ) then
309     thisStat = -1
310     return
311     endif
312 gezelter 1150
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 mmeineke 377
331     end subroutine InitializeSimGlobals
332    
333     subroutine FreeSimGlobals()
334    
335     !We free in the opposite order in which we allocate in.
336 gezelter 483
337 gezelter 1150 if (allocated(mfact)) deallocate(mfact)
338     if (allocated(groupList)) deallocate(groupList)
339     if (allocated(groupStart)) deallocate(groupStart)
340 gezelter 483 if (allocated(molMembershipList)) deallocate(molMembershipList)
341 mmeineke 377 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
342     if (allocated(excludesLocal)) deallocate(excludesLocal)
343 gezelter 1150
344 mmeineke 377 end subroutine FreeSimGlobals
345 gezelter 1150
346 mmeineke 491 pure function getNlocal() result(n)
347     integer :: n
348     n = nLocal
349 mmeineke 377 end function getNlocal
350    
351 gezelter 1150
352 mmeineke 377 end module simulation