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