ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 620
Committed: Tue Jul 15 22:29:50 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 9123 byte(s)
Log Message:
removed some debugging print statements.

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     #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 mmeineke 491 integer, public, save :: nLocal, nGlobal
24 mmeineke 377 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 chuckv 482 integer, allocatable, dimension(:), public :: molMembershipList
29 mmeineke 377
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 mmeineke 569 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
34 gezelter 570 logical, public, save :: boxIsOrthorhombic
35 mmeineke 377
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 mmeineke 491 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
59 gezelter 483 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
60     CmolMembership, &
61 mmeineke 377 status)
62    
63     type (simtype) :: setThisSim
64 mmeineke 491 integer, intent(inout) :: CnGlobal, CnLocal
65     integer, dimension(CnLocal),intent(inout) :: c_idents
66 mmeineke 377
67     integer :: CnLocalExcludes
68     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
69     integer :: CnGlobalExcludes
70     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
71 mmeineke 491 integer, dimension(CnGlobal),intent(in) :: CmolMembership
72 mmeineke 377 !! 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 mmeineke 491
87     nLocal = CnLocal
88     nGlobal = CnGlobal
89    
90 mmeineke 377 thisSim = setThisSim
91 mmeineke 620
92 mmeineke 377 rcut2 = thisSim%rcut * thisSim%rcut
93     rcut6 = rcut2 * rcut2 * rcut2
94     rlist2 = thisSim%rlist * thisSim%rlist
95    
96     nExcludes_Global = CnGlobalExcludes
97     nExcludes_Local = CnLocalExcludes
98    
99 mmeineke 491 call InitializeForceGlobals(nLocal, thisStat)
100 mmeineke 377 if (thisStat /= 0) then
101 chuckv 480 write(default_error,*) "SimSetup: InitializeForceGlobals error"
102 mmeineke 377 status = -1
103     return
104     endif
105    
106     call InitializeSimGlobals(thisStat)
107     if (thisStat /= 0) then
108 chuckv 480 write(default_error,*) "SimSetup: InitializeSimGlobals error"
109 mmeineke 377 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 gezelter 490 do i = 1, nLocal
159 mmeineke 377
160     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
161     atid(i) = me
162    
163     enddo
164     #endif
165    
166    
167 chuckv 388
168 mmeineke 377 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 mmeineke 435
177 gezelter 490 do i = 1, nGlobal
178 gezelter 483 molMemberShipList(i) = CmolMembership(i)
179 mmeineke 491 enddo
180 chuckv 482
181 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
182    
183     end subroutine SimulationSetup
184    
185 mmeineke 569 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
186 gezelter 570 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
187 mmeineke 569 integer :: cBoxIsOrthorhombic
188 mmeineke 377 integer :: smallest, status, i
189 gezelter 570
190 mmeineke 569 Hmat = cHmat
191     HmatInv = cHmatInv
192 gezelter 570 if (cBoxIsOrthorhombic .eq. 0 ) then
193 mmeineke 569 boxIsOrthorhombic = .false.
194 gezelter 570 else
195     boxIsOrthorhombic = .true.
196 mmeineke 569 endif
197    
198 mmeineke 377 return
199 mmeineke 569 end subroutine setBox
200 mmeineke 377
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 mmeineke 619 end function getRrf
250 mmeineke 377
251     function getRt() result(rt)
252     real( kind = dp ) :: rt
253     rt = thisSim%rt
254 mmeineke 619 end function getRt
255 mmeineke 377
256     function getDielect() result(dielect)
257     real( kind = dp ) :: dielect
258     dielect = thisSim%dielect
259     end function getDielect
260    
261     function SimUsesPBC() result(doesit)
262     logical :: doesit
263     doesit = thisSim%SIM_uses_PBC
264     end function SimUsesPBC
265    
266     function SimUsesLJ() result(doesit)
267     logical :: doesit
268     doesit = thisSim%SIM_uses_LJ
269     end function SimUsesLJ
270    
271     function SimUsesSticky() result(doesit)
272     logical :: doesit
273     doesit = thisSim%SIM_uses_sticky
274     end function SimUsesSticky
275    
276     function SimUsesDipoles() result(doesit)
277     logical :: doesit
278     doesit = thisSim%SIM_uses_dipoles
279     end function SimUsesDipoles
280    
281     function SimUsesRF() result(doesit)
282     logical :: doesit
283     doesit = thisSim%SIM_uses_RF
284     end function SimUsesRF
285    
286     function SimUsesGB() result(doesit)
287     logical :: doesit
288     doesit = thisSim%SIM_uses_GB
289     end function SimUsesGB
290    
291     function SimUsesEAM() result(doesit)
292     logical :: doesit
293     doesit = thisSim%SIM_uses_EAM
294     end function SimUsesEAM
295    
296     function SimUsesDirectionalAtoms() result(doesit)
297     logical :: doesit
298     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
299     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
300     end function SimUsesDirectionalAtoms
301    
302     function SimRequiresPrepairCalc() result(doesit)
303     logical :: doesit
304     doesit = thisSim%SIM_uses_EAM
305     end function SimRequiresPrepairCalc
306    
307     function SimRequiresPostpairCalc() result(doesit)
308     logical :: doesit
309     doesit = thisSim%SIM_uses_RF
310     end function SimRequiresPostpairCalc
311    
312     subroutine InitializeSimGlobals(thisStat)
313     integer, intent(out) :: thisStat
314     integer :: alloc_stat
315    
316     thisStat = 0
317    
318     call FreeSimGlobals()
319    
320     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
321     if (alloc_stat /= 0 ) then
322     thisStat = -1
323     return
324     endif
325    
326     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
327     if (alloc_stat /= 0 ) then
328     thisStat = -1
329     return
330     endif
331 chuckv 482
332 mmeineke 491 allocate(molMembershipList(nGlobal), stat=alloc_stat)
333 chuckv 482 if (alloc_stat /= 0 ) then
334     thisStat = -1
335     return
336     endif
337 mmeineke 377
338     end subroutine InitializeSimGlobals
339    
340     subroutine FreeSimGlobals()
341    
342     !We free in the opposite order in which we allocate in.
343 gezelter 483
344     if (allocated(molMembershipList)) deallocate(molMembershipList)
345 mmeineke 377 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
346     if (allocated(excludesLocal)) deallocate(excludesLocal)
347 gezelter 483
348 mmeineke 377 end subroutine FreeSimGlobals
349    
350 mmeineke 491 pure function getNlocal() result(n)
351     integer :: n
352     n = nLocal
353 mmeineke 377 end function getNlocal
354    
355    
356     end module simulation