ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/simulation_module.F90
Revision: 941
Committed: Tue Jan 13 23:01:43 2004 UTC (20 years, 5 months ago) by gezelter
File size: 7597 byte(s)
Log Message:
Changes for adding direct charge-charge interactions (with switching function)

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