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