1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
!! publication of scientific results based in part on use of the |
11 |
!! program. An acceptable form of acknowledgement is citation of |
12 |
!! the article in which the program was described (Matthew |
13 |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
!! |
18 |
!! 2. Redistributions of source code must retain the above copyright |
19 |
!! notice, this list of conditions and the following disclaimer. |
20 |
!! |
21 |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
!! notice, this list of conditions and the following disclaimer in the |
23 |
!! documentation and/or other materials provided with the |
24 |
!! distribution. |
25 |
!! |
26 |
!! This software is provided "AS IS," without a warranty of any |
27 |
!! kind. All express or implied conditions, representations and |
28 |
!! warranties, including any implied warranty of merchantability, |
29 |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
!! be liable for any damages suffered by licensee as a result of |
32 |
!! using, modifying or distributing the software or its |
33 |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
!! damages, however caused and regardless of the theory of liability, |
37 |
!! arising out of the use of or inability to use software, even if the |
38 |
!! University of Notre Dame has been advised of the possibility of |
39 |
!! such damages. |
40 |
!! |
41 |
|
42 |
!! Fortran interface to C entry plug. |
43 |
|
44 |
module simulation |
45 |
use definitions |
46 |
use status |
47 |
use linearAlgebra |
48 |
use neighborLists |
49 |
use force_globals |
50 |
use vector_class |
51 |
use atype_module |
52 |
use switcheroo |
53 |
#ifdef IS_MPI |
54 |
use mpiSimulation |
55 |
#endif |
56 |
|
57 |
implicit none |
58 |
PRIVATE |
59 |
|
60 |
#define __FORTRAN90 |
61 |
#include "brains/fSimulation.h" |
62 |
#include "UseTheForce/fSwitchingFunction.h" |
63 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
64 |
|
65 |
type (simtype), public, save :: thisSim |
66 |
|
67 |
logical, save :: simulation_setup_complete = .false. |
68 |
|
69 |
integer, public, save :: nLocal, nGlobal |
70 |
integer, public, save :: nGroups, nGroupGlobal |
71 |
integer, public, save :: nExcludes_Global = 0 |
72 |
integer, public, save :: nExcludes_Local = 0 |
73 |
integer, allocatable, dimension(:,:), public :: excludesLocal |
74 |
integer, allocatable, dimension(:), public :: excludesGlobal |
75 |
integer, allocatable, dimension(:), public :: molMembershipList |
76 |
integer, allocatable, dimension(:), public :: groupListRow |
77 |
integer, allocatable, dimension(:), public :: groupStartRow |
78 |
integer, allocatable, dimension(:), public :: groupListCol |
79 |
integer, allocatable, dimension(:), public :: groupStartCol |
80 |
integer, allocatable, dimension(:), public :: groupListLocal |
81 |
integer, allocatable, dimension(:), public :: groupStartLocal |
82 |
integer, allocatable, dimension(:), public :: nSkipsForAtom |
83 |
integer, allocatable, dimension(:,:), public :: skipsForAtom |
84 |
real(kind=dp), allocatable, dimension(:), public :: mfactRow |
85 |
real(kind=dp), allocatable, dimension(:), public :: mfactCol |
86 |
real(kind=dp), allocatable, dimension(:), public :: mfactLocal |
87 |
|
88 |
logical, allocatable, dimension(:) :: simHasAtypeMap |
89 |
#ifdef IS_MPI |
90 |
logical, allocatable, dimension(:) :: simHasAtypeMapTemp |
91 |
#endif |
92 |
|
93 |
real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv |
94 |
real(kind=dp), save :: DangerRcut |
95 |
logical, public, save :: boxIsOrthorhombic |
96 |
|
97 |
public :: SimulationSetup |
98 |
public :: getNlocal |
99 |
public :: setBox |
100 |
public :: checkBox |
101 |
public :: getDielect |
102 |
public :: SimUsesPBC |
103 |
|
104 |
public :: SimUsesDirectionalAtoms |
105 |
public :: SimUsesLennardJones |
106 |
public :: SimUsesElectrostatics |
107 |
public :: SimUsesCharges |
108 |
public :: SimUsesDipoles |
109 |
public :: SimUsesSticky |
110 |
public :: SimUsesStickyPower |
111 |
public :: SimUsesGayBerne |
112 |
public :: SimUsesEAM |
113 |
public :: SimUsesShapes |
114 |
public :: SimUsesFLARB |
115 |
public :: SimUsesRF |
116 |
public :: SimUsesSF |
117 |
public :: SimRequiresPrepairCalc |
118 |
public :: SimRequiresPostpairCalc |
119 |
public :: SimHasAtype |
120 |
public :: SimUsesSC |
121 |
public :: SimUsesMEAM |
122 |
public :: setHmatDangerousRcutValue |
123 |
|
124 |
contains |
125 |
|
126 |
subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, & |
127 |
CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, & |
128 |
CmolMembership, Cmfact, CnGroups, CglobalGroupMembership, & |
129 |
status) |
130 |
|
131 |
type (simtype) :: setThisSim |
132 |
integer, intent(inout) :: CnGlobal, CnLocal |
133 |
integer, dimension(CnLocal),intent(inout) :: c_idents |
134 |
|
135 |
integer :: CnLocalExcludes |
136 |
integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal |
137 |
integer :: CnGlobalExcludes |
138 |
integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal |
139 |
integer, dimension(CnGlobal),intent(in) :: CmolMembership |
140 |
!! Result status, success = 0, status = -1 |
141 |
integer, intent(out) :: status |
142 |
integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2 |
143 |
integer :: ia |
144 |
|
145 |
!! mass factors used for molecular cutoffs |
146 |
real ( kind = dp ), dimension(CnLocal) :: Cmfact |
147 |
integer, intent(in):: CnGroups |
148 |
integer, dimension(CnGlobal), intent(in):: CglobalGroupMembership |
149 |
integer :: maxSkipsForAtom, glPointer |
150 |
|
151 |
#ifdef IS_MPI |
152 |
integer, allocatable, dimension(:) :: c_idents_Row |
153 |
integer, allocatable, dimension(:) :: c_idents_Col |
154 |
integer :: nAtomsInRow, nGroupsInRow, aid |
155 |
integer :: nAtomsInCol, nGroupsInCol, gid |
156 |
#endif |
157 |
|
158 |
simulation_setup_complete = .false. |
159 |
status = 0 |
160 |
|
161 |
! copy C struct into fortran type |
162 |
|
163 |
nLocal = CnLocal |
164 |
nGlobal = CnGlobal |
165 |
nGroups = CnGroups |
166 |
|
167 |
thisSim = setThisSim |
168 |
|
169 |
nExcludes_Global = CnGlobalExcludes |
170 |
nExcludes_Local = CnLocalExcludes |
171 |
|
172 |
call InitializeForceGlobals(nLocal, thisStat) |
173 |
if (thisStat /= 0) then |
174 |
write(default_error,*) "SimSetup: InitializeForceGlobals error" |
175 |
status = -1 |
176 |
return |
177 |
endif |
178 |
|
179 |
call InitializeSimGlobals(thisStat) |
180 |
if (thisStat /= 0) then |
181 |
write(default_error,*) "SimSetup: InitializeSimGlobals error" |
182 |
status = -1 |
183 |
return |
184 |
endif |
185 |
|
186 |
#ifdef IS_MPI |
187 |
! We can only set up forces if mpiSimulation has been setup. |
188 |
if (.not. isMPISimSet()) then |
189 |
write(default_error,*) "MPI is not set" |
190 |
status = -1 |
191 |
return |
192 |
endif |
193 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
194 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
195 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
196 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
197 |
mynode = getMyNode() |
198 |
|
199 |
allocate(c_idents_Row(nAtomsInRow),stat=alloc_stat) |
200 |
if (alloc_stat /= 0 ) then |
201 |
status = -1 |
202 |
return |
203 |
endif |
204 |
|
205 |
allocate(c_idents_Col(nAtomsInCol),stat=alloc_stat) |
206 |
if (alloc_stat /= 0 ) then |
207 |
status = -1 |
208 |
return |
209 |
endif |
210 |
|
211 |
call gather(c_idents, c_idents_Row, plan_atom_row) |
212 |
call gather(c_idents, c_idents_Col, plan_atom_col) |
213 |
|
214 |
do i = 1, nAtomsInRow |
215 |
me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i)) |
216 |
atid_Row(i) = me |
217 |
enddo |
218 |
|
219 |
do i = 1, nAtomsInCol |
220 |
me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i)) |
221 |
atid_Col(i) = me |
222 |
enddo |
223 |
|
224 |
!! free temporary ident arrays |
225 |
if (allocated(c_idents_Col)) then |
226 |
deallocate(c_idents_Col) |
227 |
end if |
228 |
if (allocated(c_idents_Row)) then |
229 |
deallocate(c_idents_Row) |
230 |
endif |
231 |
|
232 |
#endif |
233 |
|
234 |
#ifdef IS_MPI |
235 |
allocate(groupStartRow(nGroupsInRow+1),stat=alloc_stat) |
236 |
if (alloc_stat /= 0 ) then |
237 |
status = -1 |
238 |
return |
239 |
endif |
240 |
allocate(groupStartCol(nGroupsInCol+1),stat=alloc_stat) |
241 |
if (alloc_stat /= 0 ) then |
242 |
status = -1 |
243 |
return |
244 |
endif |
245 |
allocate(groupListRow(nAtomsInRow),stat=alloc_stat) |
246 |
if (alloc_stat /= 0 ) then |
247 |
status = -1 |
248 |
return |
249 |
endif |
250 |
allocate(groupListCol(nAtomsInCol),stat=alloc_stat) |
251 |
if (alloc_stat /= 0 ) then |
252 |
status = -1 |
253 |
return |
254 |
endif |
255 |
allocate(mfactRow(nAtomsInRow),stat=alloc_stat) |
256 |
if (alloc_stat /= 0 ) then |
257 |
status = -1 |
258 |
return |
259 |
endif |
260 |
allocate(mfactCol(nAtomsInCol),stat=alloc_stat) |
261 |
if (alloc_stat /= 0 ) then |
262 |
status = -1 |
263 |
return |
264 |
endif |
265 |
allocate(mfactLocal(nLocal),stat=alloc_stat) |
266 |
if (alloc_stat /= 0 ) then |
267 |
status = -1 |
268 |
return |
269 |
endif |
270 |
|
271 |
glPointer = 1 |
272 |
|
273 |
do i = 1, nGroupsInRow |
274 |
|
275 |
gid = GroupRowToGlobal(i) |
276 |
groupStartRow(i) = glPointer |
277 |
|
278 |
do j = 1, nAtomsInRow |
279 |
aid = AtomRowToGlobal(j) |
280 |
if (CglobalGroupMembership(aid) .eq. gid) then |
281 |
groupListRow(glPointer) = j |
282 |
glPointer = glPointer + 1 |
283 |
endif |
284 |
enddo |
285 |
enddo |
286 |
groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1 |
287 |
|
288 |
glPointer = 1 |
289 |
|
290 |
do i = 1, nGroupsInCol |
291 |
|
292 |
gid = GroupColToGlobal(i) |
293 |
groupStartCol(i) = glPointer |
294 |
|
295 |
do j = 1, nAtomsInCol |
296 |
aid = AtomColToGlobal(j) |
297 |
if (CglobalGroupMembership(aid) .eq. gid) then |
298 |
groupListCol(glPointer) = j |
299 |
glPointer = glPointer + 1 |
300 |
endif |
301 |
enddo |
302 |
enddo |
303 |
groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1 |
304 |
|
305 |
mfactLocal = Cmfact |
306 |
|
307 |
call gather(mfactLocal, mfactRow, plan_atom_row) |
308 |
call gather(mfactLocal, mfactCol, plan_atom_col) |
309 |
|
310 |
if (allocated(mfactLocal)) then |
311 |
deallocate(mfactLocal) |
312 |
end if |
313 |
#else |
314 |
allocate(groupStartRow(nGroups+1),stat=alloc_stat) |
315 |
if (alloc_stat /= 0 ) then |
316 |
status = -1 |
317 |
return |
318 |
endif |
319 |
allocate(groupStartCol(nGroups+1),stat=alloc_stat) |
320 |
if (alloc_stat /= 0 ) then |
321 |
status = -1 |
322 |
return |
323 |
endif |
324 |
allocate(groupListRow(nLocal),stat=alloc_stat) |
325 |
if (alloc_stat /= 0 ) then |
326 |
status = -1 |
327 |
return |
328 |
endif |
329 |
allocate(groupListCol(nLocal),stat=alloc_stat) |
330 |
if (alloc_stat /= 0 ) then |
331 |
status = -1 |
332 |
return |
333 |
endif |
334 |
allocate(mfactRow(nLocal),stat=alloc_stat) |
335 |
if (alloc_stat /= 0 ) then |
336 |
status = -1 |
337 |
return |
338 |
endif |
339 |
allocate(mfactCol(nLocal),stat=alloc_stat) |
340 |
if (alloc_stat /= 0 ) then |
341 |
status = -1 |
342 |
return |
343 |
endif |
344 |
allocate(mfactLocal(nLocal),stat=alloc_stat) |
345 |
if (alloc_stat /= 0 ) then |
346 |
status = -1 |
347 |
return |
348 |
endif |
349 |
|
350 |
glPointer = 1 |
351 |
do i = 1, nGroups |
352 |
groupStartRow(i) = glPointer |
353 |
groupStartCol(i) = glPointer |
354 |
do j = 1, nLocal |
355 |
if (CglobalGroupMembership(j) .eq. i) then |
356 |
groupListRow(glPointer) = j |
357 |
groupListCol(glPointer) = j |
358 |
glPointer = glPointer + 1 |
359 |
endif |
360 |
enddo |
361 |
enddo |
362 |
groupStartRow(nGroups+1) = nLocal + 1 |
363 |
groupStartCol(nGroups+1) = nLocal + 1 |
364 |
|
365 |
do i = 1, nLocal |
366 |
mfactRow(i) = Cmfact(i) |
367 |
mfactCol(i) = Cmfact(i) |
368 |
end do |
369 |
|
370 |
#endif |
371 |
|
372 |
|
373 |
! We build the local atid's for both mpi and nonmpi |
374 |
do i = 1, nLocal |
375 |
|
376 |
me = getFirstMatchingElement(atypes, "c_ident", c_idents(i)) |
377 |
atid(i) = me |
378 |
|
379 |
enddo |
380 |
|
381 |
do i = 1, nExcludes_Local |
382 |
excludesLocal(1,i) = CexcludesLocal(1,i) |
383 |
excludesLocal(2,i) = CexcludesLocal(2,i) |
384 |
enddo |
385 |
|
386 |
#ifdef IS_MPI |
387 |
allocate(nSkipsForAtom(nAtomsInRow), stat=alloc_stat) |
388 |
#else |
389 |
allocate(nSkipsForAtom(nLocal), stat=alloc_stat) |
390 |
#endif |
391 |
if (alloc_stat /= 0 ) then |
392 |
thisStat = -1 |
393 |
write(*,*) 'Could not allocate nSkipsForAtom array' |
394 |
return |
395 |
endif |
396 |
|
397 |
maxSkipsForAtom = 0 |
398 |
#ifdef IS_MPI |
399 |
do j = 1, nAtomsInRow |
400 |
#else |
401 |
do j = 1, nLocal |
402 |
#endif |
403 |
nSkipsForAtom(j) = 0 |
404 |
#ifdef IS_MPI |
405 |
id1 = AtomRowToGlobal(j) |
406 |
#else |
407 |
id1 = j |
408 |
#endif |
409 |
do i = 1, nExcludes_Local |
410 |
if (excludesLocal(1,i) .eq. id1 ) then |
411 |
nSkipsForAtom(j) = nSkipsForAtom(j) + 1 |
412 |
|
413 |
if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then |
414 |
maxSkipsForAtom = nSkipsForAtom(j) |
415 |
endif |
416 |
endif |
417 |
if (excludesLocal(2,i) .eq. id1 ) then |
418 |
nSkipsForAtom(j) = nSkipsForAtom(j) + 1 |
419 |
|
420 |
if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then |
421 |
maxSkipsForAtom = nSkipsForAtom(j) |
422 |
endif |
423 |
endif |
424 |
end do |
425 |
enddo |
426 |
|
427 |
#ifdef IS_MPI |
428 |
allocate(skipsForAtom(nAtomsInRow, maxSkipsForAtom), stat=alloc_stat) |
429 |
#else |
430 |
allocate(skipsForAtom(nLocal, maxSkipsForAtom), stat=alloc_stat) |
431 |
#endif |
432 |
if (alloc_stat /= 0 ) then |
433 |
write(*,*) 'Could not allocate skipsForAtom array' |
434 |
return |
435 |
endif |
436 |
|
437 |
#ifdef IS_MPI |
438 |
do j = 1, nAtomsInRow |
439 |
#else |
440 |
do j = 1, nLocal |
441 |
#endif |
442 |
nSkipsForAtom(j) = 0 |
443 |
#ifdef IS_MPI |
444 |
id1 = AtomRowToGlobal(j) |
445 |
#else |
446 |
id1 = j |
447 |
#endif |
448 |
do i = 1, nExcludes_Local |
449 |
if (excludesLocal(1,i) .eq. id1 ) then |
450 |
nSkipsForAtom(j) = nSkipsForAtom(j) + 1 |
451 |
! exclude lists have global ID's so this line is |
452 |
! the same in MPI and non-MPI |
453 |
id2 = excludesLocal(2,i) |
454 |
skipsForAtom(j, nSkipsForAtom(j)) = id2 |
455 |
endif |
456 |
if (excludesLocal(2, i) .eq. id1 ) then |
457 |
nSkipsForAtom(j) = nSkipsForAtom(j) + 1 |
458 |
! exclude lists have global ID's so this line is |
459 |
! the same in MPI and non-MPI |
460 |
id2 = excludesLocal(1,i) |
461 |
skipsForAtom(j, nSkipsForAtom(j)) = id2 |
462 |
endif |
463 |
end do |
464 |
enddo |
465 |
|
466 |
do i = 1, nExcludes_Global |
467 |
excludesGlobal(i) = CexcludesGlobal(i) |
468 |
enddo |
469 |
|
470 |
do i = 1, nGlobal |
471 |
molMemberShipList(i) = CmolMembership(i) |
472 |
enddo |
473 |
|
474 |
call createSimHasAtype(alloc_stat) |
475 |
if (alloc_stat /= 0) then |
476 |
status = -1 |
477 |
end if |
478 |
|
479 |
if (status == 0) simulation_setup_complete = .true. |
480 |
|
481 |
end subroutine SimulationSetup |
482 |
|
483 |
subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic) |
484 |
real(kind=dp), dimension(3,3) :: cHmat, cHmatInv |
485 |
integer :: cBoxIsOrthorhombic |
486 |
integer :: smallest, status |
487 |
|
488 |
Hmat = cHmat |
489 |
HmatInv = cHmatInv |
490 |
if (cBoxIsOrthorhombic .eq. 0 ) then |
491 |
boxIsOrthorhombic = .false. |
492 |
else |
493 |
boxIsOrthorhombic = .true. |
494 |
endif |
495 |
|
496 |
call checkBox() |
497 |
return |
498 |
end subroutine setBox |
499 |
|
500 |
subroutine checkBox() |
501 |
|
502 |
integer :: i |
503 |
real(kind=dp), dimension(3) :: hx, hy, hz, ax, ay, az, piped |
504 |
character(len = statusMsgSize) :: errMsg |
505 |
|
506 |
hx = Hmat(1,:) |
507 |
hy = Hmat(2,:) |
508 |
hz = Hmat(3,:) |
509 |
|
510 |
ax = cross_product(hy, hz) |
511 |
ay = cross_product(hx, hz) |
512 |
az = cross_product(hx, hy) |
513 |
|
514 |
ax = ax / length(ax) |
515 |
ay = ay / length(ay) |
516 |
az = az / length(az) |
517 |
|
518 |
piped(1) = abs(dot_product(ax, hx)) |
519 |
piped(2) = abs(dot_product(ay, hy)) |
520 |
piped(3) = abs(dot_product(az, hz)) |
521 |
|
522 |
do i = 1, 3 |
523 |
if ((0.5_dp * piped(i)).lt.DangerRcut) then |
524 |
write(errMsg, '(a94,f9.4,a1)') 'One of the dimensions of the Periodic ' & |
525 |
// 'Box is smaller than ' // newline // tab // & |
526 |
'the largest cutoff radius' // & |
527 |
' (rCut = ', DangerRcut, ')' |
528 |
call handleError("checkBox", errMsg) |
529 |
end if |
530 |
enddo |
531 |
return |
532 |
end subroutine checkBox |
533 |
|
534 |
function getDielect() result(dielect) |
535 |
real( kind = dp ) :: dielect |
536 |
dielect = thisSim%dielect |
537 |
end function getDielect |
538 |
|
539 |
function SimUsesPBC() result(doesit) |
540 |
logical :: doesit |
541 |
doesit = thisSim%SIM_uses_PBC |
542 |
end function SimUsesPBC |
543 |
|
544 |
function SimUsesDirectionalAtoms() result(doesit) |
545 |
logical :: doesit |
546 |
doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_Sticky .or. & |
547 |
thisSim%SIM_uses_StickyPower .or. & |
548 |
thisSim%SIM_uses_GayBerne .or. thisSim%SIM_uses_Shapes |
549 |
end function SimUsesDirectionalAtoms |
550 |
|
551 |
function SimUsesLennardJones() result(doesit) |
552 |
logical :: doesit |
553 |
doesit = thisSim%SIM_uses_LennardJones |
554 |
end function SimUsesLennardJones |
555 |
|
556 |
function SimUsesElectrostatics() result(doesit) |
557 |
logical :: doesit |
558 |
doesit = thisSim%SIM_uses_Electrostatics |
559 |
end function SimUsesElectrostatics |
560 |
|
561 |
function SimUsesCharges() result(doesit) |
562 |
logical :: doesit |
563 |
doesit = thisSim%SIM_uses_Charges |
564 |
end function SimUsesCharges |
565 |
|
566 |
function SimUsesDipoles() result(doesit) |
567 |
logical :: doesit |
568 |
doesit = thisSim%SIM_uses_Dipoles |
569 |
end function SimUsesDipoles |
570 |
|
571 |
function SimUsesSticky() result(doesit) |
572 |
logical :: doesit |
573 |
doesit = thisSim%SIM_uses_Sticky |
574 |
end function SimUsesSticky |
575 |
|
576 |
function SimUsesStickyPower() result(doesit) |
577 |
logical :: doesit |
578 |
doesit = thisSim%SIM_uses_StickyPower |
579 |
end function SimUsesStickyPower |
580 |
|
581 |
function SimUsesGayBerne() result(doesit) |
582 |
logical :: doesit |
583 |
doesit = thisSim%SIM_uses_GayBerne |
584 |
end function SimUsesGayBerne |
585 |
|
586 |
function SimUsesEAM() result(doesit) |
587 |
logical :: doesit |
588 |
doesit = thisSim%SIM_uses_EAM |
589 |
end function SimUsesEAM |
590 |
|
591 |
|
592 |
function SimUsesSC() result(doesit) |
593 |
logical :: doesit |
594 |
doesit = thisSim%SIM_uses_SC |
595 |
end function SimUsesSC |
596 |
|
597 |
function SimUsesMEAM() result(doesit) |
598 |
logical :: doesit |
599 |
doesit = thisSim%SIM_uses_MEAM |
600 |
end function SimUsesMEAM |
601 |
|
602 |
|
603 |
function SimUsesShapes() result(doesit) |
604 |
logical :: doesit |
605 |
doesit = thisSim%SIM_uses_Shapes |
606 |
end function SimUsesShapes |
607 |
|
608 |
function SimUsesFLARB() result(doesit) |
609 |
logical :: doesit |
610 |
doesit = thisSim%SIM_uses_FLARB |
611 |
end function SimUsesFLARB |
612 |
|
613 |
function SimUsesRF() result(doesit) |
614 |
logical :: doesit |
615 |
doesit = thisSim%SIM_uses_RF |
616 |
end function SimUsesRF |
617 |
|
618 |
function SimUsesSF() result(doesit) |
619 |
logical :: doesit |
620 |
doesit = thisSim%SIM_uses_SF |
621 |
end function SimUsesSF |
622 |
|
623 |
function SimRequiresPrepairCalc() result(doesit) |
624 |
logical :: doesit |
625 |
doesit = thisSim%SIM_uses_EAM .or. thisSim%SIM_uses_SC & |
626 |
.or. thisSim%SIM_uses_MEAM |
627 |
end function SimRequiresPrepairCalc |
628 |
|
629 |
function SimRequiresPostpairCalc() result(doesit) |
630 |
logical :: doesit |
631 |
doesit = thisSim%SIM_uses_RF .or. thisSim%SIM_uses_SF |
632 |
end function SimRequiresPostpairCalc |
633 |
|
634 |
! Function returns true if the simulation has this atype |
635 |
function SimHasAtype(thisAtype) result(doesit) |
636 |
logical :: doesit |
637 |
integer :: thisAtype |
638 |
doesit = .false. |
639 |
if(.not.allocated(SimHasAtypeMap)) return |
640 |
|
641 |
doesit = SimHasAtypeMap(thisAtype) |
642 |
|
643 |
end function SimHasAtype |
644 |
|
645 |
subroutine createSimHasAtype(status) |
646 |
integer, intent(out) :: status |
647 |
integer :: alloc_stat |
648 |
integer :: me_i |
649 |
integer :: mpiErrors |
650 |
integer :: nAtypes |
651 |
status = 0 |
652 |
|
653 |
nAtypes = getSize(atypes) |
654 |
! Setup logical map for atypes in simulation |
655 |
if (.not.allocated(SimHasAtypeMap)) then |
656 |
allocate(SimHasAtypeMap(nAtypes),stat=alloc_stat) |
657 |
if (alloc_stat /= 0 ) then |
658 |
status = -1 |
659 |
return |
660 |
end if |
661 |
SimHasAtypeMap = .false. |
662 |
end if |
663 |
|
664 |
! Loop through the local atoms and grab the atypes present |
665 |
do me_i = 1,nLocal |
666 |
SimHasAtypeMap(atid(me_i)) = .true. |
667 |
end do |
668 |
! For MPI, we need to know all possible atypes present in |
669 |
! simulation on all processors. Use LOR operation to set map. |
670 |
#ifdef IS_MPI |
671 |
if (.not.allocated(SimHasAtypeMapTemp)) then |
672 |
allocate(SimHasAtypeMapTemp(nAtypes),stat=alloc_stat) |
673 |
if (alloc_stat /= 0 ) then |
674 |
status = -1 |
675 |
return |
676 |
end if |
677 |
end if |
678 |
call mpi_allreduce(SimHasAtypeMap, SimHasAtypeMaptemp, nAtypes, & |
679 |
mpi_logical, MPI_LOR, mpi_comm_world, mpiErrors) |
680 |
simHasAtypeMap = simHasAtypeMapTemp |
681 |
deallocate(simHasAtypeMapTemp) |
682 |
#endif |
683 |
end subroutine createSimHasAtype |
684 |
|
685 |
subroutine InitializeSimGlobals(thisStat) |
686 |
integer, intent(out) :: thisStat |
687 |
integer :: alloc_stat |
688 |
|
689 |
thisStat = 0 |
690 |
|
691 |
call FreeSimGlobals() |
692 |
|
693 |
allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat) |
694 |
if (alloc_stat /= 0 ) then |
695 |
thisStat = -1 |
696 |
return |
697 |
endif |
698 |
|
699 |
allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat) |
700 |
if (alloc_stat /= 0 ) then |
701 |
thisStat = -1 |
702 |
return |
703 |
endif |
704 |
|
705 |
allocate(molMembershipList(nGlobal), stat=alloc_stat) |
706 |
if (alloc_stat /= 0 ) then |
707 |
thisStat = -1 |
708 |
return |
709 |
endif |
710 |
|
711 |
end subroutine InitializeSimGlobals |
712 |
|
713 |
subroutine FreeSimGlobals() |
714 |
|
715 |
!We free in the opposite order in which we allocate in. |
716 |
|
717 |
if (allocated(skipsForAtom)) deallocate(skipsForAtom) |
718 |
if (allocated(nSkipsForAtom)) deallocate(nSkipsForAtom) |
719 |
if (allocated(mfactLocal)) deallocate(mfactLocal) |
720 |
if (allocated(mfactCol)) deallocate(mfactCol) |
721 |
if (allocated(mfactRow)) deallocate(mfactRow) |
722 |
if (allocated(groupListCol)) deallocate(groupListCol) |
723 |
if (allocated(groupListRow)) deallocate(groupListRow) |
724 |
if (allocated(groupStartCol)) deallocate(groupStartCol) |
725 |
if (allocated(groupStartRow)) deallocate(groupStartRow) |
726 |
if (allocated(molMembershipList)) deallocate(molMembershipList) |
727 |
if (allocated(excludesGlobal)) deallocate(excludesGlobal) |
728 |
if (allocated(excludesLocal)) deallocate(excludesLocal) |
729 |
|
730 |
end subroutine FreeSimGlobals |
731 |
|
732 |
pure function getNlocal() result(n) |
733 |
integer :: n |
734 |
n = nLocal |
735 |
end function getNlocal |
736 |
|
737 |
subroutine setHmatDangerousRcutValue(dangerWillRobinson) |
738 |
real(kind=dp), intent(in) :: dangerWillRobinson |
739 |
DangerRcut = dangerWillRobinson |
740 |
|
741 |
call checkBox() |
742 |
|
743 |
return |
744 |
end subroutine setHmatDangerousRcutValue |
745 |
|
746 |
|
747 |
|
748 |
end module simulation |