ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 1183
Committed: Fri May 21 15:58:48 2004 UTC (20 years, 1 month ago) by gezelter
File size: 11358 byte(s)
Log Message:
Major changes to skipThisPair for efficiency

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 gezelter 1183 integer, allocatable, dimension(:), public :: excludesGlobal
31     integer, allocatable, dimension(:), public :: molMembershipList
32     integer, allocatable, dimension(:), public :: groupList
33     integer, allocatable, dimension(:), public :: groupStart
34     integer, allocatable, dimension(:), public :: nSkipsForAtom
35     integer, allocatable, dimension(:,:), public :: skipsForAtom
36 gezelter 1150 real(kind=dp), allocatable, dimension(:), public :: mfact
37 mmeineke 377
38 mmeineke 569 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
39 gezelter 570 logical, public, save :: boxIsOrthorhombic
40 mmeineke 377
41     public :: SimulationSetup
42     public :: getNlocal
43     public :: setBox
44     public :: getDielect
45     public :: SimUsesPBC
46     public :: SimUsesLJ
47 gezelter 941 public :: SimUsesCharges
48 mmeineke 377 public :: SimUsesDipoles
49     public :: SimUsesSticky
50     public :: SimUsesRF
51     public :: SimUsesGB
52     public :: SimUsesEAM
53     public :: SimRequiresPrepairCalc
54     public :: SimRequiresPostpairCalc
55     public :: SimUsesDirectionalAtoms
56    
57     contains
58    
59 mmeineke 491 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
60 gezelter 483 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
61 gezelter 1150 CmolMembership, Cmfact, CnGroup, CgroupList, CgroupStart, &
62 mmeineke 377 status)
63    
64     type (simtype) :: setThisSim
65 mmeineke 491 integer, intent(inout) :: CnGlobal, CnLocal
66     integer, dimension(CnLocal),intent(inout) :: c_idents
67 mmeineke 377
68     integer :: CnLocalExcludes
69     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
70     integer :: CnGlobalExcludes
71     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
72 mmeineke 491 integer, dimension(CnGlobal),intent(in) :: CmolMembership
73 mmeineke 377 !! Result status, success = 0, status = -1
74     integer, intent(out) :: status
75 gezelter 1183 integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2
76 gezelter 1150 integer :: gStart, gEnd, groupSize, biggestGroupSize, ia
77 tim 1144
78     !! mass factors used for molecular cutoffs
79 gezelter 1150 real ( kind = dp ), dimension(CnLocal) :: Cmfact
80     integer, intent(in):: CnGroup
81     integer, dimension(CnLocal),intent(in) :: CgroupList
82     integer, dimension(CnGroup),intent(in) :: CgroupStart
83 gezelter 1183 integer :: maxSkipsForAtom
84 tim 1144
85 mmeineke 377 #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 mmeineke 491
97     nLocal = CnLocal
98     nGlobal = CnGlobal
99 gezelter 1150 nGroup = CnGroup
100 mmeineke 491
101 mmeineke 377 thisSim = setThisSim
102 mmeineke 620
103 mmeineke 377 nExcludes_Global = CnGlobalExcludes
104     nExcludes_Local = CnLocalExcludes
105    
106 mmeineke 491 call InitializeForceGlobals(nLocal, thisStat)
107 mmeineke 377 if (thisStat /= 0) then
108 chuckv 480 write(default_error,*) "SimSetup: InitializeForceGlobals error"
109 mmeineke 377 status = -1
110     return
111     endif
112    
113     call InitializeSimGlobals(thisStat)
114     if (thisStat /= 0) then
115 chuckv 480 write(default_error,*) "SimSetup: InitializeSimGlobals error"
116 mmeineke 377 status = -1
117     return
118     endif
119    
120     #ifdef IS_MPI
121     ! We can only set up forces if mpiSimulation has been setup.
122     if (.not. isMPISimSet()) then
123     write(default_error,*) "MPI is not set"
124     status = -1
125     return
126     endif
127     nrow = getNrow(plan_row)
128     ncol = getNcol(plan_col)
129     mynode = getMyNode()
130    
131     allocate(c_idents_Row(nrow),stat=alloc_stat)
132     if (alloc_stat /= 0 ) then
133     status = -1
134     return
135     endif
136    
137     allocate(c_idents_Col(ncol),stat=alloc_stat)
138     if (alloc_stat /= 0 ) then
139     status = -1
140     return
141     endif
142    
143     call gather(c_idents, c_idents_Row, plan_row)
144     call gather(c_idents, c_idents_Col, plan_col)
145    
146     do i = 1, nrow
147     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
148     atid_Row(i) = me
149     enddo
150    
151     do i = 1, ncol
152     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
153     atid_Col(i) = me
154     enddo
155    
156     !! free temporary ident arrays
157     if (allocated(c_idents_Col)) then
158     deallocate(c_idents_Col)
159     end if
160     if (allocated(c_idents_Row)) then
161     deallocate(c_idents_Row)
162     endif
163    
164 chuckv 648 #endif
165    
166     ! We build the local atid's for both mpi and nonmpi
167 gezelter 490 do i = 1, nLocal
168 mmeineke 377
169     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
170     atid(i) = me
171    
172     enddo
173    
174     do i = 1, nExcludes_Local
175     excludesLocal(1,i) = CexcludesLocal(1,i)
176     excludesLocal(2,i) = CexcludesLocal(2,i)
177     enddo
178 gezelter 1183
179     maxSkipsForAtom = 0
180     do j = 1, nLocal
181     nSkipsForAtom(j) = 0
182     #ifdef IS_MPI
183     id1 = tagRow(j)
184     #else
185     id1 = j
186     #endif
187     do i = 1, nExcludes_Local
188     if (excludesLocal(1,i) .eq. id1 ) then
189     nSkipsForAtom(j) = nSkipsForAtom(j) + 1
190    
191     if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then
192     maxSkipsForAtom = nSkipsForAtom(j)
193     endif
194     endif
195     if (excludesLocal(2,i) .eq. id1 ) then
196     nSkipsForAtom(j) = nSkipsForAtom(j) + 1
197    
198     if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then
199     maxSkipsForAtom = nSkipsForAtom(j)
200     endif
201     endif
202     end do
203     enddo
204    
205     allocate(skipsForAtom(nLocal, maxSkipsForAtom), stat=alloc_stat)
206     if (alloc_stat /= 0 ) then
207     write(*,*) 'Could not allocate skipsForAtom array'
208     return
209     endif
210    
211     do j = 1, nLocal
212     nSkipsForAtom(j) = 0
213     #ifdef IS_MPI
214     id1 = tagRow(j)
215     #else
216     id1 = j
217     #endif
218     do i = 1, nExcludes_Local
219     if (excludesLocal(1,i) .eq. id1 ) then
220     nSkipsForAtom(j) = nSkipsForAtom(j) + 1
221     #ifdef IS_MPI
222     id2 = tagColumn(excludesLocal(2,i))
223     #else
224     id2 = excludesLocal(2,i)
225     #endif
226     skipsForAtom(j, nSkipsForAtom(j)) = id2
227     endif
228     if (excludesLocal(2, i) .eq. id2 ) then
229     nSkipsForAtom(j) = nSkipsForAtom(j) + 1
230     #ifdef IS_MPI
231     id2 = tagColumn(excludesLocal(1,i))
232     #else
233     id2 = excludesLocal(1,i)
234     #endif
235     skipsForAtom(j, nSkipsForAtom(j)) = id2
236     endif
237     end do
238     enddo
239 mmeineke 377
240     do i = 1, nExcludes_Global
241     excludesGlobal(i) = CexcludesGlobal(i)
242     enddo
243 mmeineke 435
244 gezelter 490 do i = 1, nGlobal
245 gezelter 483 molMemberShipList(i) = CmolMembership(i)
246 gezelter 1150 enddo
247    
248     biggestGroupSize = 0
249     do i = 1, nGroup
250     groupStart(i) = CgroupStart(i)
251     groupSize = 0
252     gStart = groupStart(i)
253     if (i .eq. nGroup) then
254     gEnd = nLocal
255     else
256     gEnd = CgroupStart(i+1) - 1
257     endif
258     do ia = gStart, gEnd
259     groupList(ia) = CgroupList(ia)
260     groupSize = groupSize + 1
261     enddo
262     if (groupSize .gt. biggestGroupSize) biggestGroupSize = groupSize
263     enddo
264     groupStart(nGroup+1) = nLocal+1
265 chuckv 482
266 gezelter 1150 do i = 1, nLocal
267     mfact(i) = Cmfact(i)
268     enddo
269    
270 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
271    
272     end subroutine SimulationSetup
273    
274 mmeineke 569 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
275 gezelter 570 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
276 mmeineke 569 integer :: cBoxIsOrthorhombic
277 mmeineke 377 integer :: smallest, status, i
278 gezelter 570
279 mmeineke 569 Hmat = cHmat
280     HmatInv = cHmatInv
281 gezelter 570 if (cBoxIsOrthorhombic .eq. 0 ) then
282 mmeineke 569 boxIsOrthorhombic = .false.
283 gezelter 570 else
284     boxIsOrthorhombic = .true.
285 mmeineke 569 endif
286    
287 mmeineke 377 return
288 mmeineke 569 end subroutine setBox
289 mmeineke 377
290     function getDielect() result(dielect)
291     real( kind = dp ) :: dielect
292     dielect = thisSim%dielect
293     end function getDielect
294    
295     function SimUsesPBC() result(doesit)
296     logical :: doesit
297     doesit = thisSim%SIM_uses_PBC
298     end function SimUsesPBC
299    
300     function SimUsesLJ() result(doesit)
301     logical :: doesit
302     doesit = thisSim%SIM_uses_LJ
303     end function SimUsesLJ
304    
305     function SimUsesSticky() result(doesit)
306     logical :: doesit
307     doesit = thisSim%SIM_uses_sticky
308     end function SimUsesSticky
309    
310 gezelter 941 function SimUsesCharges() result(doesit)
311     logical :: doesit
312     doesit = thisSim%SIM_uses_charges
313     end function SimUsesCharges
314    
315 mmeineke 377 function SimUsesDipoles() result(doesit)
316     logical :: doesit
317     doesit = thisSim%SIM_uses_dipoles
318     end function SimUsesDipoles
319    
320     function SimUsesRF() result(doesit)
321     logical :: doesit
322     doesit = thisSim%SIM_uses_RF
323     end function SimUsesRF
324    
325     function SimUsesGB() result(doesit)
326     logical :: doesit
327     doesit = thisSim%SIM_uses_GB
328     end function SimUsesGB
329    
330     function SimUsesEAM() result(doesit)
331     logical :: doesit
332     doesit = thisSim%SIM_uses_EAM
333     end function SimUsesEAM
334    
335     function SimUsesDirectionalAtoms() result(doesit)
336     logical :: doesit
337     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
338     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
339     end function SimUsesDirectionalAtoms
340    
341     function SimRequiresPrepairCalc() result(doesit)
342     logical :: doesit
343     doesit = thisSim%SIM_uses_EAM
344     end function SimRequiresPrepairCalc
345    
346     function SimRequiresPostpairCalc() result(doesit)
347     logical :: doesit
348     doesit = thisSim%SIM_uses_RF
349     end function SimRequiresPostpairCalc
350    
351     subroutine InitializeSimGlobals(thisStat)
352     integer, intent(out) :: thisStat
353     integer :: alloc_stat
354    
355     thisStat = 0
356    
357     call FreeSimGlobals()
358 gezelter 1183
359     allocate(nSkipsForAtom(nLocal), stat=alloc_stat)
360     if (alloc_stat /= 0 ) then
361     thisStat = -1
362     return
363     endif
364 mmeineke 377
365     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
366     if (alloc_stat /= 0 ) then
367     thisStat = -1
368     return
369     endif
370    
371     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
372     if (alloc_stat /= 0 ) then
373     thisStat = -1
374     return
375     endif
376 chuckv 482
377 mmeineke 491 allocate(molMembershipList(nGlobal), stat=alloc_stat)
378 chuckv 482 if (alloc_stat /= 0 ) then
379     thisStat = -1
380     return
381     endif
382 gezelter 1150
383     allocate(groupStart(nGroup+1), stat=alloc_stat)
384     if (alloc_stat /= 0 ) then
385     thisStat = -1
386     return
387     endif
388    
389     allocate(groupList(nLocal), stat=alloc_stat)
390     if (alloc_stat /= 0 ) then
391     thisStat = -1
392     return
393     endif
394    
395     allocate(mfact(nLocal), stat=alloc_stat)
396     if (alloc_stat /= 0 ) then
397     thisStat = -1
398     return
399     endif
400 mmeineke 377
401     end subroutine InitializeSimGlobals
402    
403     subroutine FreeSimGlobals()
404    
405     !We free in the opposite order in which we allocate in.
406 gezelter 483
407 gezelter 1183 if (allocated(skipsForAtom)) deallocate(skipsForAtom)
408 gezelter 1150 if (allocated(mfact)) deallocate(mfact)
409     if (allocated(groupList)) deallocate(groupList)
410     if (allocated(groupStart)) deallocate(groupStart)
411 gezelter 483 if (allocated(molMembershipList)) deallocate(molMembershipList)
412 mmeineke 377 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
413     if (allocated(excludesLocal)) deallocate(excludesLocal)
414 gezelter 1150
415 mmeineke 377 end subroutine FreeSimGlobals
416 gezelter 1150
417 mmeineke 491 pure function getNlocal() result(n)
418     integer :: n
419     n = nLocal
420 mmeineke 377 end function getNlocal
421    
422 gezelter 1150
423 mmeineke 377 end module simulation