ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 626
Committed: Wed Jul 16 21:30:56 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 7382 byte(s)
Log Message:
Changed how cutoffs were handled from C. Now notifyCutoffs in Fortran notifies those who need the information of any changes to cutoffs.

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 mmeineke 569 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
31 gezelter 570 logical, public, save :: boxIsOrthorhombic
32 mmeineke 377
33     public :: SimulationSetup
34     public :: getNlocal
35     public :: setBox
36     public :: getDielect
37     public :: SimUsesPBC
38     public :: SimUsesLJ
39     public :: SimUsesDipoles
40     public :: SimUsesSticky
41     public :: SimUsesRF
42     public :: SimUsesGB
43     public :: SimUsesEAM
44     public :: SimRequiresPrepairCalc
45     public :: SimRequiresPostpairCalc
46     public :: SimUsesDirectionalAtoms
47    
48     contains
49    
50 mmeineke 491 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
51 gezelter 483 CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
52     CmolMembership, &
53 mmeineke 377 status)
54    
55     type (simtype) :: setThisSim
56 mmeineke 491 integer, intent(inout) :: CnGlobal, CnLocal
57     integer, dimension(CnLocal),intent(inout) :: c_idents
58 mmeineke 377
59     integer :: CnLocalExcludes
60     integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
61     integer :: CnGlobalExcludes
62     integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
63 mmeineke 491 integer, dimension(CnGlobal),intent(in) :: CmolMembership
64 mmeineke 377 !! Result status, success = 0, status = -1
65     integer, intent(out) :: status
66     integer :: i, me, thisStat, alloc_stat, myNode
67     #ifdef IS_MPI
68     integer, allocatable, dimension(:) :: c_idents_Row
69     integer, allocatable, dimension(:) :: c_idents_Col
70     integer :: nrow
71     integer :: ncol
72     #endif
73    
74     simulation_setup_complete = .false.
75     status = 0
76    
77     ! copy C struct into fortran type
78 mmeineke 491
79     nLocal = CnLocal
80     nGlobal = CnGlobal
81    
82 mmeineke 377 thisSim = setThisSim
83 mmeineke 620
84 mmeineke 377 nExcludes_Global = CnGlobalExcludes
85     nExcludes_Local = CnLocalExcludes
86    
87 mmeineke 491 call InitializeForceGlobals(nLocal, thisStat)
88 mmeineke 377 if (thisStat /= 0) then
89 chuckv 480 write(default_error,*) "SimSetup: InitializeForceGlobals error"
90 mmeineke 377 status = -1
91     return
92     endif
93    
94     call InitializeSimGlobals(thisStat)
95     if (thisStat /= 0) then
96 chuckv 480 write(default_error,*) "SimSetup: InitializeSimGlobals error"
97 mmeineke 377 status = -1
98     return
99     endif
100    
101     #ifdef IS_MPI
102     ! We can only set up forces if mpiSimulation has been setup.
103     if (.not. isMPISimSet()) then
104     write(default_error,*) "MPI is not set"
105     status = -1
106     return
107     endif
108     nrow = getNrow(plan_row)
109     ncol = getNcol(plan_col)
110     mynode = getMyNode()
111    
112     allocate(c_idents_Row(nrow),stat=alloc_stat)
113     if (alloc_stat /= 0 ) then
114     status = -1
115     return
116     endif
117    
118     allocate(c_idents_Col(ncol),stat=alloc_stat)
119     if (alloc_stat /= 0 ) then
120     status = -1
121     return
122     endif
123    
124     call gather(c_idents, c_idents_Row, plan_row)
125     call gather(c_idents, c_idents_Col, plan_col)
126    
127     do i = 1, nrow
128     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
129     atid_Row(i) = me
130     enddo
131    
132     do i = 1, ncol
133     me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
134     atid_Col(i) = me
135     enddo
136    
137     !! free temporary ident arrays
138     if (allocated(c_idents_Col)) then
139     deallocate(c_idents_Col)
140     end if
141     if (allocated(c_idents_Row)) then
142     deallocate(c_idents_Row)
143     endif
144    
145     #else
146 gezelter 490 do i = 1, nLocal
147 mmeineke 377
148     me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
149     atid(i) = me
150    
151     enddo
152     #endif
153    
154    
155 chuckv 388
156 mmeineke 377 do i = 1, nExcludes_Local
157     excludesLocal(1,i) = CexcludesLocal(1,i)
158     excludesLocal(2,i) = CexcludesLocal(2,i)
159     enddo
160    
161     do i = 1, nExcludes_Global
162     excludesGlobal(i) = CexcludesGlobal(i)
163     enddo
164 mmeineke 435
165 gezelter 490 do i = 1, nGlobal
166 gezelter 483 molMemberShipList(i) = CmolMembership(i)
167 mmeineke 491 enddo
168 chuckv 482
169 mmeineke 377 if (status == 0) simulation_setup_complete = .true.
170    
171     end subroutine SimulationSetup
172    
173 mmeineke 569 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
174 gezelter 570 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
175 mmeineke 569 integer :: cBoxIsOrthorhombic
176 mmeineke 377 integer :: smallest, status, i
177 gezelter 570
178 mmeineke 569 Hmat = cHmat
179     HmatInv = cHmatInv
180 gezelter 570 if (cBoxIsOrthorhombic .eq. 0 ) then
181 mmeineke 569 boxIsOrthorhombic = .false.
182 gezelter 570 else
183     boxIsOrthorhombic = .true.
184 mmeineke 569 endif
185    
186 mmeineke 377 return
187 mmeineke 569 end subroutine setBox
188 mmeineke 377
189     function getDielect() result(dielect)
190     real( kind = dp ) :: dielect
191     dielect = thisSim%dielect
192     end function getDielect
193    
194     function SimUsesPBC() result(doesit)
195     logical :: doesit
196     doesit = thisSim%SIM_uses_PBC
197     end function SimUsesPBC
198    
199     function SimUsesLJ() result(doesit)
200     logical :: doesit
201     doesit = thisSim%SIM_uses_LJ
202     end function SimUsesLJ
203    
204     function SimUsesSticky() result(doesit)
205     logical :: doesit
206     doesit = thisSim%SIM_uses_sticky
207     end function SimUsesSticky
208    
209     function SimUsesDipoles() result(doesit)
210     logical :: doesit
211     doesit = thisSim%SIM_uses_dipoles
212     end function SimUsesDipoles
213    
214     function SimUsesRF() result(doesit)
215     logical :: doesit
216     doesit = thisSim%SIM_uses_RF
217     end function SimUsesRF
218    
219     function SimUsesGB() result(doesit)
220     logical :: doesit
221     doesit = thisSim%SIM_uses_GB
222     end function SimUsesGB
223    
224     function SimUsesEAM() result(doesit)
225     logical :: doesit
226     doesit = thisSim%SIM_uses_EAM
227     end function SimUsesEAM
228    
229     function SimUsesDirectionalAtoms() result(doesit)
230     logical :: doesit
231     doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. &
232     thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF
233     end function SimUsesDirectionalAtoms
234    
235     function SimRequiresPrepairCalc() result(doesit)
236     logical :: doesit
237     doesit = thisSim%SIM_uses_EAM
238     end function SimRequiresPrepairCalc
239    
240     function SimRequiresPostpairCalc() result(doesit)
241     logical :: doesit
242     doesit = thisSim%SIM_uses_RF
243     end function SimRequiresPostpairCalc
244    
245     subroutine InitializeSimGlobals(thisStat)
246     integer, intent(out) :: thisStat
247     integer :: alloc_stat
248    
249     thisStat = 0
250    
251     call FreeSimGlobals()
252    
253     allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat)
254     if (alloc_stat /= 0 ) then
255     thisStat = -1
256     return
257     endif
258    
259     allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat)
260     if (alloc_stat /= 0 ) then
261     thisStat = -1
262     return
263     endif
264 chuckv 482
265 mmeineke 491 allocate(molMembershipList(nGlobal), stat=alloc_stat)
266 chuckv 482 if (alloc_stat /= 0 ) then
267     thisStat = -1
268     return
269     endif
270 mmeineke 377
271     end subroutine InitializeSimGlobals
272    
273     subroutine FreeSimGlobals()
274    
275     !We free in the opposite order in which we allocate in.
276 gezelter 483
277     if (allocated(molMembershipList)) deallocate(molMembershipList)
278 mmeineke 377 if (allocated(excludesGlobal)) deallocate(excludesGlobal)
279     if (allocated(excludesLocal)) deallocate(excludesLocal)
280 gezelter 483
281 mmeineke 377 end subroutine FreeSimGlobals
282    
283 mmeineke 491 pure function getNlocal() result(n)
284     integer :: n
285     n = nLocal
286 mmeineke 377 end function getNlocal
287    
288    
289     end module simulation