ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 569
Committed: Tue Jul 1 21:33:45 2003 UTC (21 years ago) by mmeineke
File size: 9166 byte(s)
Log Message:
working on adding the box matrix to everything.

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 #ifdef IS_MPI
10 use mpiSimulation
11 #endif
12
13 implicit none
14 PRIVATE
15
16 #define __FORTRAN90
17 #include "fSimulation.h"
18
19 type (simtype), public :: thisSim
20
21 logical, save :: simulation_setup_complete = .false.
22
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,3), save :: Hmat, HmatInv
34 logical, save :: boxIsOrthorhombic
35
36 public :: SimulationSetup
37 public :: getNlocal
38 public :: setBox
39 public :: setRcut
40 public :: getRcut
41 public :: getRlist
42 public :: getRrf
43 public :: getRt
44 public :: getDielect
45 public :: SimUsesPBC
46 public :: SimUsesLJ
47 public :: SimUsesDipoles
48 public :: SimUsesSticky
49 public :: SimUsesRF
50 public :: SimUsesGB
51 public :: SimUsesEAM
52 public :: SimRequiresPrepairCalc
53 public :: SimRequiresPostpairCalc
54 public :: SimUsesDirectionalAtoms
55
56 contains
57
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) :: 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
75 #ifdef IS_MPI
76 integer, allocatable, dimension(:) :: c_idents_Row
77 integer, allocatable, dimension(:) :: c_idents_Col
78 integer :: nrow
79 integer :: ncol
80 #endif
81
82 simulation_setup_complete = .false.
83 status = 0
84
85 ! copy C struct into fortran type
86
87 nLocal = CnLocal
88 nGlobal = CnGlobal
89
90 thisSim = setThisSim
91 rcut2 = thisSim%rcut * thisSim%rcut
92 rcut6 = rcut2 * rcut2 * rcut2
93 rlist2 = thisSim%rlist * thisSim%rlist
94 box = thisSim%box
95
96 nExcludes_Global = CnGlobalExcludes
97 nExcludes_Local = CnLocalExcludes
98
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
112
113 #ifdef IS_MPI
114 ! We can only set up forces if mpiSimulation has been setup.
115 if (.not. isMPISimSet()) then
116 write(default_error,*) "MPI is not set"
117 status = -1
118 return
119 endif
120 nrow = getNrow(plan_row)
121 ncol = getNcol(plan_col)
122 mynode = getMyNode()
123
124 allocate(c_idents_Row(nrow),stat=alloc_stat)
125 if (alloc_stat /= 0 ) then
126 status = -1
127 return
128 endif
129
130 allocate(c_idents_Col(ncol),stat=alloc_stat)
131 if (alloc_stat /= 0 ) then
132 status = -1
133 return
134 endif
135
136 call gather(c_idents, c_idents_Row, plan_row)
137 call gather(c_idents, c_idents_Col, plan_col)
138
139 do i = 1, nrow
140 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
141 atid_Row(i) = me
142 enddo
143
144 do i = 1, ncol
145 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
146 atid_Col(i) = me
147 enddo
148
149 !! free temporary ident arrays
150 if (allocated(c_idents_Col)) then
151 deallocate(c_idents_Col)
152 end if
153 if (allocated(c_idents_Row)) then
154 deallocate(c_idents_Row)
155 endif
156
157 #else
158 do i = 1, nLocal
159
160 me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
161 atid(i) = me
162
163 enddo
164 #endif
165
166
167
168 do i = 1, nExcludes_Local
169 excludesLocal(1,i) = CexcludesLocal(1,i)
170 excludesLocal(2,i) = CexcludesLocal(2,i)
171 enddo
172
173 do i = 1, nExcludes_Global
174 excludesGlobal(i) = CexcludesGlobal(i)
175 enddo
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(cHmat, cHmatInv, cBoxIsOrthorhombic)
186 real(kind=dp), dimension(3,3) :: cHmat, cHamtInv
187 integer :: cBoxIsOrthorhombic
188 integer :: smallest, status, i
189
190 Hmat = cHmat
191 HmatInv = cHmatInv
192 if(cBoxisOrthorhombic .eq. 0 )then
193 boxIsOrthorhombic = .false.
194 else boxIsOrthorhombic = .true.
195 endif
196
197
198 return
199 end subroutine setBox
200
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
207 status = 0
208 return
209 end subroutine setRcut
210
211 subroutine getRcut(thisrcut,rc2,rc6,status)
212 real( kind = dp ), intent(out) :: thisrcut
213 real( kind = dp ), intent(out), optional :: rc2
214 real( kind = dp ), intent(out), optional :: rc6
215 integer, optional :: status
216
217 if (present(status)) status = 0
218
219 if (.not.simulation_setup_complete ) then
220 if (present(status)) status = -1
221 return
222 end if
223
224 thisrcut = thisSim%rcut
225 if(present(rc2)) rc2 = rcut2
226 if(present(rc6)) rc6 = rcut6
227 end subroutine getRcut
228
229 subroutine getRlist(thisrlist,rl2,status)
230 real( kind = dp ), intent(out) :: thisrlist
231 real( kind = dp ), intent(out), optional :: rl2
232
233 integer, optional :: status
234
235 if (present(status)) status = 0
236
237 if (.not.simulation_setup_complete ) then
238 if (present(status)) status = -1
239 return
240 end if
241
242 thisrlist = thisSim%rlist
243 if(present(rl2)) rl2 = rlist2
244 end subroutine getRlist
245
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)
253 real( kind = dp ) :: rt
254 rt = thisSim%rt
255 end function getRt
256
257 function getDielect() result(dielect)
258 real( kind = dp ) :: dielect
259 dielect = thisSim%dielect
260 end function getDielect
261
262 function SimUsesPBC() result(doesit)
263 logical :: doesit
264 doesit = thisSim%SIM_uses_PBC
265 end function SimUsesPBC
266
267 function SimUsesLJ() result(doesit)
268 logical :: doesit
269 doesit = thisSim%SIM_uses_LJ
270 end function SimUsesLJ
271
272 function SimUsesSticky() result(doesit)
273 logical :: doesit
274 doesit = thisSim%SIM_uses_sticky
275 end function SimUsesSticky
276
277 function SimUsesDipoles() result(doesit)
278 logical :: doesit
279 doesit = thisSim%SIM_uses_dipoles
280 end function SimUsesDipoles
281
282 function SimUsesRF() result(doesit)
283 logical :: doesit
284 doesit = thisSim%SIM_uses_RF
285 end function SimUsesRF
286
287 function SimUsesGB() result(doesit)
288 logical :: doesit
289 doesit = thisSim%SIM_uses_GB
290 end function SimUsesGB
291
292 function SimUsesEAM() result(doesit)
293 logical :: doesit
294 doesit = thisSim%SIM_uses_EAM
295 end function SimUsesEAM
296
297 function SimUsesDirectionalAtoms() result(doesit)
298 logical :: doesit
299 doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
300 thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
301 end function SimUsesDirectionalAtoms
302
303 function SimRequiresPrepairCalc() result(doesit)
304 logical :: doesit
305 doesit = thisSim%SIM_uses_EAM
306 end function SimRequiresPrepairCalc
307
308 function SimRequiresPostpairCalc() result(doesit)
309 logical :: doesit
310 doesit = thisSim%SIM_uses_RF
311 end function SimRequiresPostpairCalc
312
313 subroutine InitializeSimGlobals(thisStat)
314 integer, intent(out) :: thisStat
315 integer :: alloc_stat
316
317 thisStat = 0
318
319 call FreeSimGlobals()
320
321 allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
322 if (alloc_stat /= 0 ) then
323 thisStat = -1
324 return
325 endif
326
327 allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
328 if (alloc_stat /= 0 ) then
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
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(n)
352 integer :: n
353 n = nLocal
354 end function getNlocal
355
356
357 end module simulation