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

# Content
1 !! 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 use switcheroo
10 #ifdef IS_MPI
11 use mpiSimulation
12 #endif
13
14 implicit none
15 PRIVATE
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
38
39 public :: SimulationSetup
40 public :: getNlocal
41 public :: setBox
42 public :: getDielect
43 public :: SimUsesPBC
44 public :: SimUsesLJ
45 public :: SimUsesCharges
46 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 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
58 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
59 CmolMembership, Cmfact, CnGroup, CgroupList, CgroupStart, &
60 status)
61
62 type (simtype) :: setThisSim
63 integer, intent(inout) :: CnGlobal, CnLocal
64 integer, dimension(CnLocal),intent(inout) :: c_idents
65
66 integer :: CnLocalExcludes
67 integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
68 integer :: CnGlobalExcludes
69 integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
70 integer, dimension(CnGlobal),intent(in) :: CmolMembership
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(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
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
94 nLocal = CnLocal
95 nGlobal = CnGlobal
96 nGroup = CnGroup
97
98 thisSim = setThisSim
99
100 nExcludes_Global = CnGlobalExcludes
101 nExcludes_Local = CnLocalExcludes
102
103 call InitializeForceGlobals(nLocal, thisStat)
104 if (thisStat /= 0) then
105 write(default_error,*) "SimSetup: InitializeForceGlobals error"
106 status = -1
107 return
108 endif
109
110 call InitializeSimGlobals(thisStat)
111 if (thisStat /= 0) then
112 write(default_error,*) "SimSetup: InitializeSimGlobals error"
113 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 #endif
162
163 ! We build the local atid's for both mpi and nonmpi
164 do i = 1, nLocal
165
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
180 do i = 1, nGlobal
181 molMemberShipList(i) = CmolMembership(i)
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
209
210 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
211 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
212 integer :: cBoxIsOrthorhombic
213 integer :: smallest, status, i
214
215 Hmat = cHmat
216 HmatInv = cHmatInv
217 if (cBoxIsOrthorhombic .eq. 0 ) then
218 boxIsOrthorhombic = .false.
219 else
220 boxIsOrthorhombic = .true.
221 endif
222
223 return
224 end subroutine setBox
225
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 function SimUsesCharges() result(doesit)
247 logical :: doesit
248 doesit = thisSim%SIM_uses_charges
249 end function SimUsesCharges
250
251 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
307 allocate(molMembershipList(nGlobal), stat=alloc_stat)
308 if (alloc_stat /= 0 ) then
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
333 subroutine FreeSimGlobals()
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
344 end subroutine FreeSimGlobals
345
346 pure function getNlocal() result(n)
347 integer :: n
348 n = nLocal
349 end function getNlocal
350
351
352 end module simulation