ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 435
Committed: Fri Mar 28 19:33:37 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 9795 byte(s)
Log Message:
fixed a bug where the Excludes were not being created properly

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 lj
10 #ifdef IS_MPI
11 use mpiSimulation
12 #endif
13
14 implicit none
15 PRIVATE
16
17 #define __FORTRAN90
18 #include "fSimulation.h"
19
20 type (simtype), public :: thisSim
21
22 logical, save :: simulation_setup_complete = .false.
23
24 integer, public, save :: natoms
25 integer, public, save :: nExcludes_Global = 0
26 integer, public, save :: nExcludes_Local = 0
27 integer, allocatable, dimension(:,:), public :: excludesLocal
28 integer, allocatable, dimension(:), public :: excludesGlobal
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
35
36 public :: SimulationSetup
37 public :: getNlocal
38 public :: setBox
39 public :: setBox_3d
40 public :: getBox
41 public :: setRcut
42 public :: getRcut
43 public :: getRlist
44 public :: getRrf
45 public :: getRt
46 public :: getDielect
47 public :: SimUsesPBC
48 public :: SimUsesLJ
49 public :: SimUsesDipoles
50 public :: SimUsesSticky
51 public :: SimUsesRF
52 public :: SimUsesGB
53 public :: SimUsesEAM
54 public :: SimRequiresPrepairCalc
55 public :: SimRequiresPostpairCalc
56 public :: SimUsesDirectionalAtoms
57
58 interface getBox
59 module procedure getBox_3d
60 module procedure getBox_1d
61 end interface
62
63 interface setBox
64 module procedure setBox_3d
65 module procedure setBox_1d
66 end interface
67
68 contains
69
70 subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
71 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
72 status)
73
74 type (simtype) :: setThisSim
75 integer, intent(inout) :: nComponents
76 integer, dimension(nComponents),intent(inout) :: c_idents
77
78 integer :: CnLocalExcludes
79 integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
80 integer :: CnGlobalExcludes
81 integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
82 !! Result status, success = 0, status = -1
83 integer, intent(out) :: status
84 integer :: i, me, thisStat, alloc_stat, myNode
85 #ifdef IS_MPI
86 integer, allocatable, dimension(:) :: c_idents_Row
87 integer, allocatable, dimension(:) :: c_idents_Col
88 integer :: nrow
89 integer :: ncol
90 #endif
91
92 simulation_setup_complete = .false.
93 status = 0
94
95 ! copy C struct into fortran type
96 thisSim = setThisSim
97 natoms = nComponents
98 rcut2 = thisSim%rcut * thisSim%rcut
99 rcut6 = rcut2 * rcut2 * rcut2
100 rlist2 = thisSim%rlist * thisSim%rlist
101 box = thisSim%box
102
103 nExcludes_Global = CnGlobalExcludes
104 nExcludes_Local = CnLocalExcludes
105
106 call InitializeForceGlobals(natoms, thisStat)
107 if (thisStat /= 0) then
108 status = -1
109 return
110 endif
111
112 call InitializeSimGlobals(thisStat)
113 if (thisStat /= 0) then
114 status = -1
115 return
116 endif
117
118 #ifdef IS_MPI
119 ! We can only set up forces if mpiSimulation has been setup.
120 if (.not. isMPISimSet()) then
121 write(default_error,*) "MPI is not set"
122 status = -1
123 return
124 endif
125 nrow = getNrow(plan_row)
126 ncol = getNcol(plan_col)
127 mynode = getMyNode()
128
129 allocate(c_idents_Row(nrow),stat=alloc_stat)
130 if (alloc_stat /= 0 ) then
131 status = -1
132 return
133 endif
134
135 allocate(c_idents_Col(ncol),stat=alloc_stat)
136 if (alloc_stat /= 0 ) then
137 status = -1
138 return
139 endif
140
141 call gather(c_idents, c_idents_Row, plan_row)
142 call gather(c_idents, c_idents_Col, plan_col)
143
144 do i = 1, nrow
145 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
146 atid_Row(i) = me
147 enddo
148
149 do i = 1, ncol
150 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
151 atid_Col(i) = me
152 enddo
153
154 !! free temporary ident arrays
155 if (allocated(c_idents_Col)) then
156 deallocate(c_idents_Col)
157 end if
158 if (allocated(c_idents_Row)) then
159 deallocate(c_idents_Row)
160 endif
161
162 #else
163 do i = 1, nComponents
164
165 me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
166 atid(i) = me
167
168 enddo
169 #endif
170
171 !! Create neighbor lists
172 call expandNeighborList(nComponents, thisStat)
173 if (thisStat /= 0) then
174 status = -1
175 return
176 endif
177
178
179 do i = 1, nExcludes_Local
180 excludesLocal(1,i) = CexcludesLocal(1,i)
181 excludesLocal(2,i) = CexcludesLocal(2,i)
182 enddo
183
184 do i = 1, nExcludes_Global
185 excludesGlobal(i) = CexcludesGlobal(i)
186 enddo
187
188 if (status == 0) simulation_setup_complete = .true.
189
190 end subroutine SimulationSetup
191
192 subroutine setBox_3d(new_box_size)
193 real(kind=dp), dimension(3) :: new_box_size
194 integer :: smallest, status, i
195
196 thisSim%box = new_box_size
197 box = thisSim%box
198
199 smallest = 1
200 do i = 2, 3
201 if (new_box_size(i) .lt. new_box_size(smallest)) smallest = i
202 end do
203 if (thisSim%rcut .gt. 0.5_dp * new_box_size(smallest)) &
204 call setRcut(0.5_dp * new_box_size(smallest), status)
205 return
206 end subroutine setBox_3d
207
208 subroutine setBox_1d(dim, new_box_size)
209 integer :: dim, status
210 real(kind=dp) :: new_box_size
211 thisSim%box(dim) = new_box_size
212 box(dim) = thisSim%box(dim)
213 if (thisSim%rcut .gt. 0.5_dp * new_box_size) &
214 call setRcut(0.5_dp * new_box_size, status)
215 end subroutine setBox_1d
216
217 subroutine setRcut(new_rcut, status)
218 real(kind = dp) :: new_rcut
219 integer :: myStatus, status
220 thisSim%rcut = new_rcut
221 rcut2 = thisSim%rcut * thisSim%rcut
222 rcut6 = rcut2 * rcut2 * rcut2
223 myStatus = 0
224 call LJ_new_rcut(new_rcut, myStatus)
225 if (myStatus .ne. 0) then
226 write(default_error, *) 'LJ module refused our rcut!'
227 status = -1
228 return
229 endif
230 status = 0
231 return
232 end subroutine setRcut
233
234 function getBox_3d() result(thisBox)
235 real( kind = dp ), dimension(3) :: thisBox
236 thisBox = thisSim%box
237 end function getBox_3d
238
239 function getBox_1d(dim) result(thisBox)
240 integer, intent(in) :: dim
241 real( kind = dp ) :: thisBox
242
243 thisBox = thisSim%box(dim)
244 end function getBox_1d
245
246 subroutine getRcut(thisrcut,rc2,rc6,status)
247 real( kind = dp ), intent(out) :: thisrcut
248 real( kind = dp ), intent(out), optional :: rc2
249 real( kind = dp ), intent(out), optional :: rc6
250 integer, optional :: status
251
252 if (present(status)) status = 0
253
254 if (.not.simulation_setup_complete ) then
255 if (present(status)) status = -1
256 return
257 end if
258
259 thisrcut = thisSim%rcut
260 if(present(rc2)) rc2 = rcut2
261 if(present(rc6)) rc6 = rcut6
262 end subroutine getRcut
263
264 subroutine getRlist(thisrlist,rl2,status)
265 real( kind = dp ), intent(out) :: thisrlist
266 real( kind = dp ), intent(out), optional :: rl2
267
268 integer, optional :: status
269
270 if (present(status)) status = 0
271
272 if (.not.simulation_setup_complete ) then
273 if (present(status)) status = -1
274 return
275 end if
276
277 thisrlist = thisSim%rlist
278 if(present(rl2)) rl2 = rlist2
279 end subroutine getRlist
280
281 function getRrf() result(rrf)
282 real( kind = dp ) :: rrf
283 rrf = thisSim%rrf
284 write(*,*) 'getRrf = ', rrf, thisSim%rrf
285 end function getRrf
286
287 function getRt() result(rt)
288 real( kind = dp ) :: rt
289 rt = thisSim%rt
290 end function getRt
291
292 function getDielect() result(dielect)
293 real( kind = dp ) :: dielect
294 dielect = thisSim%dielect
295 end function getDielect
296
297 function SimUsesPBC() result(doesit)
298 logical :: doesit
299 doesit = thisSim%SIM_uses_PBC
300 end function SimUsesPBC
301
302 function SimUsesLJ() result(doesit)
303 logical :: doesit
304 doesit = thisSim%SIM_uses_LJ
305 end function SimUsesLJ
306
307 function SimUsesSticky() result(doesit)
308 logical :: doesit
309 doesit = thisSim%SIM_uses_sticky
310 end function SimUsesSticky
311
312 function SimUsesDipoles() result(doesit)
313 logical :: doesit
314 doesit = thisSim%SIM_uses_dipoles
315 end function SimUsesDipoles
316
317 function SimUsesRF() result(doesit)
318 logical :: doesit
319 doesit = thisSim%SIM_uses_RF
320 end function SimUsesRF
321
322 function SimUsesGB() result(doesit)
323 logical :: doesit
324 doesit = thisSim%SIM_uses_GB
325 end function SimUsesGB
326
327 function SimUsesEAM() result(doesit)
328 logical :: doesit
329 doesit = thisSim%SIM_uses_EAM
330 end function SimUsesEAM
331
332 function SimUsesDirectionalAtoms() result(doesit)
333 logical :: doesit
334 doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
335 thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
336 end function SimUsesDirectionalAtoms
337
338 function SimRequiresPrepairCalc() result(doesit)
339 logical :: doesit
340 doesit = thisSim%SIM_uses_EAM
341 end function SimRequiresPrepairCalc
342
343 function SimRequiresPostpairCalc() result(doesit)
344 logical :: doesit
345 doesit = thisSim%SIM_uses_RF
346 end function SimRequiresPostpairCalc
347
348 subroutine InitializeSimGlobals(thisStat)
349 integer, intent(out) :: thisStat
350 integer :: alloc_stat
351
352 thisStat = 0
353
354 call FreeSimGlobals()
355
356 allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
357 if (alloc_stat /= 0 ) then
358 thisStat = -1
359 return
360 endif
361
362 allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
363 if (alloc_stat /= 0 ) then
364 thisStat = -1
365 return
366 endif
367
368 end subroutine InitializeSimGlobals
369
370 subroutine FreeSimGlobals()
371
372 !We free in the opposite order in which we allocate in.
373
374 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
375 if (allocated(excludesLocal)) deallocate(excludesLocal)
376
377 end subroutine FreeSimGlobals
378
379 pure function getNlocal() result(nlocal)
380 integer :: nlocal
381 nlocal = natoms
382 end function getNlocal
383
384
385 end module simulation