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 :: 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 |
subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, & |
51 |
CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, & |
52 |
CmolMembership, & |
53 |
status) |
54 |
|
55 |
type (simtype) :: setThisSim |
56 |
integer, intent(inout) :: CnGlobal, CnLocal |
57 |
integer, dimension(CnLocal),intent(inout) :: c_idents |
58 |
|
59 |
integer :: CnLocalExcludes |
60 |
integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal |
61 |
integer :: CnGlobalExcludes |
62 |
integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal |
63 |
integer, dimension(CnGlobal),intent(in) :: CmolMembership |
64 |
!! 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 |
|
79 |
nLocal = CnLocal |
80 |
nGlobal = CnGlobal |
81 |
|
82 |
thisSim = setThisSim |
83 |
|
84 |
nExcludes_Global = CnGlobalExcludes |
85 |
nExcludes_Local = CnLocalExcludes |
86 |
|
87 |
call InitializeForceGlobals(nLocal, thisStat) |
88 |
if (thisStat /= 0) then |
89 |
write(default_error,*) "SimSetup: InitializeForceGlobals error" |
90 |
status = -1 |
91 |
return |
92 |
endif |
93 |
|
94 |
call InitializeSimGlobals(thisStat) |
95 |
if (thisStat /= 0) then |
96 |
write(default_error,*) "SimSetup: InitializeSimGlobals error" |
97 |
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 |
#endif |
146 |
|
147 |
! We build the local atid's for both mpi and nonmpi |
148 |
do i = 1, nLocal |
149 |
|
150 |
me = getFirstMatchingElement(atypes, "c_ident", c_idents(i)) |
151 |
atid(i) = me |
152 |
|
153 |
enddo |
154 |
|
155 |
|
156 |
|
157 |
|
158 |
do i = 1, nExcludes_Local |
159 |
excludesLocal(1,i) = CexcludesLocal(1,i) |
160 |
excludesLocal(2,i) = CexcludesLocal(2,i) |
161 |
enddo |
162 |
|
163 |
do i = 1, nExcludes_Global |
164 |
excludesGlobal(i) = CexcludesGlobal(i) |
165 |
enddo |
166 |
|
167 |
do i = 1, nGlobal |
168 |
molMemberShipList(i) = CmolMembership(i) |
169 |
enddo |
170 |
|
171 |
if (status == 0) simulation_setup_complete = .true. |
172 |
|
173 |
end subroutine SimulationSetup |
174 |
|
175 |
subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic) |
176 |
real(kind=dp), dimension(3,3) :: cHmat, cHmatInv |
177 |
integer :: cBoxIsOrthorhombic |
178 |
integer :: smallest, status, i |
179 |
|
180 |
Hmat = cHmat |
181 |
HmatInv = cHmatInv |
182 |
if (cBoxIsOrthorhombic .eq. 0 ) then |
183 |
boxIsOrthorhombic = .false. |
184 |
else |
185 |
boxIsOrthorhombic = .true. |
186 |
endif |
187 |
|
188 |
return |
189 |
end subroutine setBox |
190 |
|
191 |
function getDielect() result(dielect) |
192 |
real( kind = dp ) :: dielect |
193 |
dielect = thisSim%dielect |
194 |
end function getDielect |
195 |
|
196 |
function SimUsesPBC() result(doesit) |
197 |
logical :: doesit |
198 |
doesit = thisSim%SIM_uses_PBC |
199 |
end function SimUsesPBC |
200 |
|
201 |
function SimUsesLJ() result(doesit) |
202 |
logical :: doesit |
203 |
doesit = thisSim%SIM_uses_LJ |
204 |
end function SimUsesLJ |
205 |
|
206 |
function SimUsesSticky() result(doesit) |
207 |
logical :: doesit |
208 |
doesit = thisSim%SIM_uses_sticky |
209 |
end function SimUsesSticky |
210 |
|
211 |
function SimUsesDipoles() result(doesit) |
212 |
logical :: doesit |
213 |
doesit = thisSim%SIM_uses_dipoles |
214 |
end function SimUsesDipoles |
215 |
|
216 |
function SimUsesRF() result(doesit) |
217 |
logical :: doesit |
218 |
doesit = thisSim%SIM_uses_RF |
219 |
end function SimUsesRF |
220 |
|
221 |
function SimUsesGB() result(doesit) |
222 |
logical :: doesit |
223 |
doesit = thisSim%SIM_uses_GB |
224 |
end function SimUsesGB |
225 |
|
226 |
function SimUsesEAM() result(doesit) |
227 |
logical :: doesit |
228 |
doesit = thisSim%SIM_uses_EAM |
229 |
end function SimUsesEAM |
230 |
|
231 |
function SimUsesDirectionalAtoms() result(doesit) |
232 |
logical :: doesit |
233 |
doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. & |
234 |
thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF |
235 |
end function SimUsesDirectionalAtoms |
236 |
|
237 |
function SimRequiresPrepairCalc() result(doesit) |
238 |
logical :: doesit |
239 |
doesit = thisSim%SIM_uses_EAM |
240 |
end function SimRequiresPrepairCalc |
241 |
|
242 |
function SimRequiresPostpairCalc() result(doesit) |
243 |
logical :: doesit |
244 |
doesit = thisSim%SIM_uses_RF |
245 |
end function SimRequiresPostpairCalc |
246 |
|
247 |
subroutine InitializeSimGlobals(thisStat) |
248 |
integer, intent(out) :: thisStat |
249 |
integer :: alloc_stat |
250 |
|
251 |
thisStat = 0 |
252 |
|
253 |
call FreeSimGlobals() |
254 |
|
255 |
allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat) |
256 |
if (alloc_stat /= 0 ) then |
257 |
thisStat = -1 |
258 |
return |
259 |
endif |
260 |
|
261 |
allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat) |
262 |
if (alloc_stat /= 0 ) then |
263 |
thisStat = -1 |
264 |
return |
265 |
endif |
266 |
|
267 |
allocate(molMembershipList(nGlobal), stat=alloc_stat) |
268 |
if (alloc_stat /= 0 ) then |
269 |
thisStat = -1 |
270 |
return |
271 |
endif |
272 |
|
273 |
end subroutine InitializeSimGlobals |
274 |
|
275 |
subroutine FreeSimGlobals() |
276 |
|
277 |
!We free in the opposite order in which we allocate in. |
278 |
|
279 |
if (allocated(molMembershipList)) deallocate(molMembershipList) |
280 |
if (allocated(excludesGlobal)) deallocate(excludesGlobal) |
281 |
if (allocated(excludesLocal)) deallocate(excludesLocal) |
282 |
|
283 |
end subroutine FreeSimGlobals |
284 |
|
285 |
pure function getNlocal() result(n) |
286 |
integer :: n |
287 |
n = nLocal |
288 |
end function getNlocal |
289 |
|
290 |
|
291 |
end module simulation |