ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 491
Committed: Fri Apr 11 18:46:37 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 9700 byte(s)
Log Message:
fixed a memory bug in Fortran, where molMembershipArray was declared nLocal instead of nGlobal.

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