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