ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/simulation.F90
Revision: 1286
Committed: Wed Sep 10 17:57:55 2008 UTC (16 years, 8 months ago) by gezelter
File size: 27564 byte(s)
Log Message:
Added support for scaled 1-2, 1-3 and 1-4 interactions.

File Contents

# Content
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 = 0
72 integer, public, save :: nOneTwo = 0
73 integer, public, save :: nOneThree = 0
74 integer, public, save :: nOneFour = 0
75
76 integer, allocatable, dimension(:,:), public :: excludes
77 integer, allocatable, dimension(:), public :: molMembershipList
78 integer, allocatable, dimension(:), public :: groupListRow
79 integer, allocatable, dimension(:), public :: groupStartRow
80 integer, allocatable, dimension(:), public :: groupListCol
81 integer, allocatable, dimension(:), public :: groupStartCol
82 integer, allocatable, dimension(:), public :: groupListLocal
83 integer, allocatable, dimension(:), public :: groupStartLocal
84 integer, allocatable, dimension(:), public :: nSkipsForAtom
85 integer, allocatable, dimension(:,:), public :: skipsForAtom
86 integer, allocatable, dimension(:), public :: nTopoPairsForAtom
87 integer, allocatable, dimension(:,:), public :: toposForAtom
88 integer, allocatable, dimension(:,:), public :: topoDistance
89 real(kind=dp), allocatable, dimension(:), public :: mfactRow
90 real(kind=dp), allocatable, dimension(:), public :: mfactCol
91 real(kind=dp), allocatable, dimension(:), public :: mfactLocal
92
93 logical, allocatable, dimension(:) :: simHasAtypeMap
94 #ifdef IS_MPI
95 logical, allocatable, dimension(:) :: simHasAtypeMapTemp
96 #endif
97
98 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
99 real(kind=dp), save :: DangerRcut
100 logical, public, save :: boxIsOrthorhombic
101
102 public :: SimulationSetup
103 public :: getNlocal
104 public :: setBox
105 public :: checkBox
106 public :: SimUsesPBC
107 public :: SimUsesAtomicVirial
108
109 public :: SimUsesDirectionalAtoms
110 public :: SimUsesLennardJones
111 public :: SimUsesElectrostatics
112 public :: SimUsesCharges
113 public :: SimUsesDipoles
114 public :: SimUsesSticky
115 public :: SimUsesStickyPower
116 public :: SimUsesGayBerne
117 public :: SimUsesEAM
118 public :: SimUsesShapes
119 public :: SimUsesFLARB
120 public :: SimUsesRF
121 public :: SimUsesSF
122 public :: SimUsesSP
123 public :: SimUsesBoxDipole
124 public :: SimRequiresPrepairCalc
125 public :: SimRequiresPostpairCalc
126 public :: SimHasAtype
127 public :: SimUsesSC
128 public :: SimUsesMNM
129 public :: setHmatDangerousRcutValue
130
131 contains
132
133 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
134 CnExcludes, Cexcludes, CnOneTwo, ConeTwo, CnOneThree, ConeThree, &
135 CnOneFour, ConeFour, CmolMembership, Cmfact, CnGroups, &
136 CglobalGroupMembership, status)
137
138 type (simtype) :: setThisSim
139 integer, intent(inout) :: CnGlobal, CnLocal
140 integer, dimension(CnLocal),intent(inout) :: c_idents
141
142 integer :: CnExcludes
143 integer, dimension(2,CnExcludes), intent(in) :: Cexcludes
144 integer :: CnOneTwo
145 integer, dimension(2,CnOneTwo), intent(in) :: ConeTwo
146 integer :: CnOneThree
147 integer, dimension(2,CnOneThree), intent(in) :: ConeThree
148 integer :: CnOneFour
149 integer, dimension(2,CnOneFour), intent(in) :: ConeFour
150
151 integer, dimension(CnGlobal),intent(in) :: CmolMembership
152 !! Result status, success = 0, status = -1
153 integer, intent(out) :: status
154 integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2
155 integer :: ia, jend
156
157 !! mass factors used for molecular cutoffs
158 real ( kind = dp ), dimension(CnLocal) :: Cmfact
159 integer, intent(in):: CnGroups
160 integer, dimension(CnGlobal), intent(in):: CglobalGroupMembership
161 integer :: maxSkipsForAtom, maxToposForAtom, glPointer
162
163 #ifdef IS_MPI
164 integer, allocatable, dimension(:) :: c_idents_Row
165 integer, allocatable, dimension(:) :: c_idents_Col
166 integer :: nAtomsInRow, nGroupsInRow, aid
167 integer :: nAtomsInCol, nGroupsInCol, gid
168 #endif
169
170 simulation_setup_complete = .false.
171 status = 0
172
173 ! copy C struct into fortran type
174
175 nLocal = CnLocal
176 nGlobal = CnGlobal
177 nGroups = CnGroups
178
179 thisSim = setThisSim
180
181 nExcludes = CnExcludes
182
183 call InitializeForceGlobals(nLocal, thisStat)
184 if (thisStat /= 0) then
185 write(default_error,*) "SimSetup: InitializeForceGlobals error"
186 status = -1
187 return
188 endif
189
190 call InitializeSimGlobals(thisStat)
191 if (thisStat /= 0) then
192 write(default_error,*) "SimSetup: InitializeSimGlobals error"
193 status = -1
194 return
195 endif
196
197 #ifdef IS_MPI
198 ! We can only set up forces if mpiSimulation has been setup.
199 if (.not. isMPISimSet()) then
200 write(default_error,*) "MPI is not set"
201 status = -1
202 return
203 endif
204 nAtomsInRow = getNatomsInRow(plan_atom_row)
205 nAtomsInCol = getNatomsInCol(plan_atom_col)
206 nGroupsInRow = getNgroupsInRow(plan_group_row)
207 nGroupsInCol = getNgroupsInCol(plan_group_col)
208 mynode = getMyNode()
209
210 allocate(c_idents_Row(nAtomsInRow),stat=alloc_stat)
211 if (alloc_stat /= 0 ) then
212 status = -1
213 return
214 endif
215
216 allocate(c_idents_Col(nAtomsInCol),stat=alloc_stat)
217 if (alloc_stat /= 0 ) then
218 status = -1
219 return
220 endif
221
222 call gather(c_idents, c_idents_Row, plan_atom_row)
223 call gather(c_idents, c_idents_Col, plan_atom_col)
224
225 do i = 1, nAtomsInRow
226 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
227 atid_Row(i) = me
228 enddo
229
230 do i = 1, nAtomsInCol
231 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
232 atid_Col(i) = me
233 enddo
234
235 !! free temporary ident arrays
236 if (allocated(c_idents_Col)) then
237 deallocate(c_idents_Col)
238 end if
239 if (allocated(c_idents_Row)) then
240 deallocate(c_idents_Row)
241 endif
242
243 #endif
244
245 #ifdef IS_MPI
246 allocate(groupStartRow(nGroupsInRow+1),stat=alloc_stat)
247 if (alloc_stat /= 0 ) then
248 status = -1
249 return
250 endif
251 allocate(groupStartCol(nGroupsInCol+1),stat=alloc_stat)
252 if (alloc_stat /= 0 ) then
253 status = -1
254 return
255 endif
256 allocate(groupListRow(nAtomsInRow),stat=alloc_stat)
257 if (alloc_stat /= 0 ) then
258 status = -1
259 return
260 endif
261 allocate(groupListCol(nAtomsInCol),stat=alloc_stat)
262 if (alloc_stat /= 0 ) then
263 status = -1
264 return
265 endif
266 allocate(mfactRow(nAtomsInRow),stat=alloc_stat)
267 if (alloc_stat /= 0 ) then
268 status = -1
269 return
270 endif
271 allocate(mfactCol(nAtomsInCol),stat=alloc_stat)
272 if (alloc_stat /= 0 ) then
273 status = -1
274 return
275 endif
276 allocate(mfactLocal(nLocal),stat=alloc_stat)
277 if (alloc_stat /= 0 ) then
278 status = -1
279 return
280 endif
281
282 glPointer = 1
283
284 do i = 1, nGroupsInRow
285
286 gid = GroupRowToGlobal(i)
287 groupStartRow(i) = glPointer
288
289 do j = 1, nAtomsInRow
290 aid = AtomRowToGlobal(j)
291 if (CglobalGroupMembership(aid) .eq. gid) then
292 groupListRow(glPointer) = j
293 glPointer = glPointer + 1
294 endif
295 enddo
296 enddo
297 groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1
298
299 glPointer = 1
300
301 do i = 1, nGroupsInCol
302
303 gid = GroupColToGlobal(i)
304 groupStartCol(i) = glPointer
305
306 do j = 1, nAtomsInCol
307 aid = AtomColToGlobal(j)
308 if (CglobalGroupMembership(aid) .eq. gid) then
309 groupListCol(glPointer) = j
310 glPointer = glPointer + 1
311 endif
312 enddo
313 enddo
314 groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1
315
316 mfactLocal = Cmfact
317
318 call gather(mfactLocal, mfactRow, plan_atom_row)
319 call gather(mfactLocal, mfactCol, plan_atom_col)
320
321 if (allocated(mfactLocal)) then
322 deallocate(mfactLocal)
323 end if
324 #else
325 allocate(groupStartRow(nGroups+1),stat=alloc_stat)
326 if (alloc_stat /= 0 ) then
327 status = -1
328 return
329 endif
330 allocate(groupStartCol(nGroups+1),stat=alloc_stat)
331 if (alloc_stat /= 0 ) then
332 status = -1
333 return
334 endif
335 allocate(groupListRow(nLocal),stat=alloc_stat)
336 if (alloc_stat /= 0 ) then
337 status = -1
338 return
339 endif
340 allocate(groupListCol(nLocal),stat=alloc_stat)
341 if (alloc_stat /= 0 ) then
342 status = -1
343 return
344 endif
345 allocate(mfactRow(nLocal),stat=alloc_stat)
346 if (alloc_stat /= 0 ) then
347 status = -1
348 return
349 endif
350 allocate(mfactCol(nLocal),stat=alloc_stat)
351 if (alloc_stat /= 0 ) then
352 status = -1
353 return
354 endif
355 allocate(mfactLocal(nLocal),stat=alloc_stat)
356 if (alloc_stat /= 0 ) then
357 status = -1
358 return
359 endif
360
361 glPointer = 1
362 do i = 1, nGroups
363 groupStartRow(i) = glPointer
364 groupStartCol(i) = glPointer
365 do j = 1, nLocal
366 if (CglobalGroupMembership(j) .eq. i) then
367 groupListRow(glPointer) = j
368 groupListCol(glPointer) = j
369 glPointer = glPointer + 1
370 endif
371 enddo
372 enddo
373 groupStartRow(nGroups+1) = nLocal + 1
374 groupStartCol(nGroups+1) = nLocal + 1
375
376 do i = 1, nLocal
377 mfactRow(i) = Cmfact(i)
378 mfactCol(i) = Cmfact(i)
379 end do
380
381 #endif
382
383 ! We build the local atid's for both mpi and nonmpi
384 do i = 1, nLocal
385 me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
386 atid(i) = me
387 enddo
388
389 do i = 1, nExcludes
390 excludes(1,i) = Cexcludes(1,i)
391 excludes(2,i) = Cexcludes(2,i)
392 enddo
393
394 #ifdef IS_MPI
395 allocate(nSkipsForAtom(nAtomsInRow), stat=alloc_stat)
396 #else
397 allocate(nSkipsForAtom(nLocal), stat=alloc_stat)
398 #endif
399 if (alloc_stat /= 0 ) then
400 thisStat = -1
401 write(*,*) 'Could not allocate nSkipsForAtom array'
402 return
403 endif
404
405 maxSkipsForAtom = 0
406
407 #ifdef IS_MPI
408 jend = nAtomsInRow
409 #else
410 jend = nLocal
411 #endif
412
413 do j = 1, jend
414 nSkipsForAtom(j) = 0
415 #ifdef IS_MPI
416 id1 = AtomRowToGlobal(j)
417 #else
418 id1 = j
419 #endif
420 do i = 1, nExcludes
421 if (excludes(1,i) .eq. id1 ) then
422 nSkipsForAtom(j) = nSkipsForAtom(j) + 1
423
424 if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then
425 maxSkipsForAtom = nSkipsForAtom(j)
426 endif
427 endif
428 if (excludes(2,i) .eq. id1 ) then
429 nSkipsForAtom(j) = nSkipsForAtom(j) + 1
430
431 if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then
432 maxSkipsForAtom = nSkipsForAtom(j)
433 endif
434 endif
435 end do
436 enddo
437
438 #ifdef IS_MPI
439 allocate(skipsForAtom(nAtomsInRow, maxSkipsForAtom), stat=alloc_stat)
440 #else
441 allocate(skipsForAtom(nLocal, maxSkipsForAtom), stat=alloc_stat)
442 #endif
443 if (alloc_stat /= 0 ) then
444 write(*,*) 'Could not allocate skipsForAtom array'
445 return
446 endif
447
448 #ifdef IS_MPI
449 jend = nAtomsInRow
450 #else
451 jend = nLocal
452 #endif
453 do j = 1, jend
454 nSkipsForAtom(j) = 0
455 #ifdef IS_MPI
456 id1 = AtomRowToGlobal(j)
457 #else
458 id1 = j
459 #endif
460 do i = 1, nExcludes
461 if (excludes(1,i) .eq. id1 ) then
462 nSkipsForAtom(j) = nSkipsForAtom(j) + 1
463 ! exclude lists have global ID's so this line is
464 ! the same in MPI and non-MPI
465 id2 = excludes(2,i)
466 skipsForAtom(j, nSkipsForAtom(j)) = id2
467 endif
468 if (excludes(2, i) .eq. id1 ) then
469 nSkipsForAtom(j) = nSkipsForAtom(j) + 1
470 ! exclude lists have global ID's so this line is
471 ! the same in MPI and non-MPI
472 id2 = excludes(1,i)
473 skipsForAtom(j, nSkipsForAtom(j)) = id2
474 endif
475 end do
476 enddo
477
478 do i = 1, nGlobal
479 molMemberShipList(i) = CmolMembership(i)
480 enddo
481
482 #ifdef IS_MPI
483 allocate(nTopoPairsForAtom(nAtomsInRow), stat=alloc_stat)
484 #else
485 allocate(nTopoPairsForAtom(nLocal), stat=alloc_stat)
486 #endif
487 if (alloc_stat /= 0 ) then
488 thisStat = -1
489 write(*,*) 'Could not allocate nTopoPairsForAtom array'
490 return
491 endif
492
493 #ifdef IS_MPI
494 jend = nAtomsInRow
495 #else
496 jend = nLocal
497 #endif
498
499 do j = 1, jend
500 nTopoPairsForAtom(j) = 0
501 #ifdef IS_MPI
502 id1 = AtomRowToGlobal(j)
503 #else
504 id1 = j
505 #endif
506 do i = 1, CnOneTwo
507 if (ConeTwo(1,i) .eq. id1 ) then
508 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
509 endif
510 if (ConeTwo(2,i) .eq. id1 ) then
511 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
512 endif
513 end do
514
515 do i = 1, CnOneThree
516 if (ConeThree(1,i) .eq. id1 ) then
517 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
518 endif
519 if (ConeThree(2,i) .eq. id1 ) then
520 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
521 endif
522 end do
523
524 do i = 1, CnOneFour
525 if (ConeFour(1,i) .eq. id1 ) then
526 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
527 endif
528 if (ConeFour(2,i) .eq. id1 ) then
529 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
530 endif
531 end do
532 enddo
533
534 maxToposForAtom = maxval(nTopoPairsForAtom)
535 #ifdef IS_MPI
536 allocate(toposForAtom(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
537 allocate(topoDistance(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
538 #else
539 allocate(toposForAtom(nLocal, maxToposForAtom), stat=alloc_stat)
540 allocate(topoDistance(nLocal, maxToposForAtom), stat=alloc_stat)
541 #endif
542 if (alloc_stat /= 0 ) then
543 write(*,*) 'Could not allocate topoDistance array'
544 return
545 endif
546
547 #ifdef IS_MPI
548 jend = nAtomsInRow
549 #else
550 jend = nLocal
551 #endif
552 do j = 1, jend
553 nTopoPairsForAtom(j) = 0
554 #ifdef IS_MPI
555 id1 = AtomRowToGlobal(j)
556 #else
557 id1 = j
558 #endif
559 do i = 1, CnOneTwo
560 if (ConeTwo(1,i) .eq. id1 ) then
561 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
562 id2 = ConeTwo(2,i)
563 toposForAtom(j, nTopoPairsForAtom(j)) = id2
564 topoDistance(j, nTopoPairsForAtom(j)) = 1
565 endif
566 if (ConeTwo(2, i) .eq. id1 ) then
567 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
568 id2 = ConeTwo(1,i)
569 toposForAtom(j, nTopoPairsForAtom(j)) = id2
570 topoDistance(j, nTopoPairsForAtom(j)) = 1
571 endif
572 end do
573
574 do i = 1, CnOneThree
575 if (ConeThree(1,i) .eq. id1 ) then
576 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
577 id2 = ConeThree(2,i)
578 toposForAtom(j, nTopoPairsForAtom(j)) = id2
579 topoDistance(j, nTopoPairsForAtom(j)) = 2
580 endif
581 if (ConeThree(2, i) .eq. id1 ) then
582 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
583 id2 = ConeThree(1,i)
584 toposForAtom(j, nTopoPairsForAtom(j)) = id2
585 topoDistance(j, nTopoPairsForAtom(j)) = 2
586 endif
587 end do
588
589 do i = 1, CnOneFour
590 if (ConeFour(1,i) .eq. id1 ) then
591 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
592 id2 = ConeFour(2,i)
593 toposForAtom(j, nTopoPairsForAtom(j)) = id2
594 topoDistance(j, nTopoPairsForAtom(j)) = 3
595 endif
596 if (ConeFour(2, i) .eq. id1 ) then
597 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
598 id2 = ConeFour(1,i)
599 toposForAtom(j, nTopoPairsForAtom(j)) = id2
600 topoDistance(j, nTopoPairsForAtom(j)) = 3
601 endif
602 end do
603 enddo
604
605
606 do i = 1, nLocal
607 do j = 1, nTopoPairsForAtom(i)
608
609 write(*,*) 'pair: ', i, ', ', toposForAtom(i,j), ' = ', topoDistance(i,j)
610 enddo
611 enddo
612
613
614 call createSimHasAtype(alloc_stat)
615 if (alloc_stat /= 0) then
616 status = -1
617 end if
618
619 if (status == 0) simulation_setup_complete = .true.
620
621 end subroutine SimulationSetup
622
623 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
624 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
625 integer :: cBoxIsOrthorhombic
626 integer :: smallest, status
627
628 Hmat = cHmat
629 HmatInv = cHmatInv
630 if (cBoxIsOrthorhombic .eq. 0 ) then
631 boxIsOrthorhombic = .false.
632 else
633 boxIsOrthorhombic = .true.
634 endif
635
636 call checkBox()
637 return
638 end subroutine setBox
639
640 subroutine checkBox()
641
642 integer :: i
643 real(kind=dp), dimension(3) :: hx, hy, hz, ax, ay, az, piped
644 character(len = statusMsgSize) :: errMsg
645
646 hx = Hmat(1,:)
647 hy = Hmat(2,:)
648 hz = Hmat(3,:)
649
650 ax = cross_product(hy, hz)
651 ay = cross_product(hx, hz)
652 az = cross_product(hx, hy)
653
654 ax = ax / length(ax)
655 ay = ay / length(ay)
656 az = az / length(az)
657
658 piped(1) = abs(dot_product(ax, hx))
659 piped(2) = abs(dot_product(ay, hy))
660 piped(3) = abs(dot_product(az, hz))
661
662 do i = 1, 3
663 if ((0.5_dp * piped(i)).lt.DangerRcut) then
664 write(errMsg, '(a94,f9.4,a1)') 'One of the dimensions of the Periodic ' &
665 // 'Box is smaller than ' // newline // tab // &
666 'the largest cutoff radius' // &
667 ' (rCut = ', DangerRcut, ')'
668 call handleError("checkBox", errMsg)
669
670 end if
671 enddo
672 return
673 end subroutine checkBox
674
675 function SimUsesPBC() result(doesit)
676 logical :: doesit
677 doesit = thisSim%SIM_uses_PBC
678 end function SimUsesPBC
679
680 function SimUsesAtomicVirial() result(doesit)
681 logical :: doesit
682 doesit = thisSim%SIM_uses_AtomicVirial
683 end function SimUsesAtomicVirial
684
685 function SimUsesDirectionalAtoms() result(doesit)
686 logical :: doesit
687 doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_Sticky .or. &
688 thisSim%SIM_uses_StickyPower .or. &
689 thisSim%SIM_uses_GayBerne .or. thisSim%SIM_uses_Shapes
690 end function SimUsesDirectionalAtoms
691
692 function SimUsesLennardJones() result(doesit)
693 logical :: doesit
694 doesit = thisSim%SIM_uses_LennardJones
695 end function SimUsesLennardJones
696
697 function SimUsesElectrostatics() result(doesit)
698 logical :: doesit
699 doesit = thisSim%SIM_uses_Electrostatics
700 end function SimUsesElectrostatics
701
702 function SimUsesCharges() result(doesit)
703 logical :: doesit
704 doesit = thisSim%SIM_uses_Charges
705 end function SimUsesCharges
706
707 function SimUsesDipoles() result(doesit)
708 logical :: doesit
709 doesit = thisSim%SIM_uses_Dipoles
710 end function SimUsesDipoles
711
712 function SimUsesSticky() result(doesit)
713 logical :: doesit
714 doesit = thisSim%SIM_uses_Sticky
715 end function SimUsesSticky
716
717 function SimUsesStickyPower() result(doesit)
718 logical :: doesit
719 doesit = thisSim%SIM_uses_StickyPower
720 end function SimUsesStickyPower
721
722 function SimUsesGayBerne() result(doesit)
723 logical :: doesit
724 doesit = thisSim%SIM_uses_GayBerne
725 end function SimUsesGayBerne
726
727 function SimUsesEAM() result(doesit)
728 logical :: doesit
729 doesit = thisSim%SIM_uses_EAM
730 end function SimUsesEAM
731
732
733 function SimUsesSC() result(doesit)
734 logical :: doesit
735 doesit = thisSim%SIM_uses_SC
736 end function SimUsesSC
737
738 function SimUsesMNM() result(doesit)
739 logical :: doesit
740 doesit = thisSim%SIM_uses_MNM
741 end function SimUsesMNM
742
743
744 function SimUsesShapes() result(doesit)
745 logical :: doesit
746 doesit = thisSim%SIM_uses_Shapes
747 end function SimUsesShapes
748
749 function SimUsesFLARB() result(doesit)
750 logical :: doesit
751 doesit = thisSim%SIM_uses_FLARB
752 end function SimUsesFLARB
753
754 function SimUsesRF() result(doesit)
755 logical :: doesit
756 doesit = thisSim%SIM_uses_RF
757 end function SimUsesRF
758
759 function SimUsesSF() result(doesit)
760 logical :: doesit
761 doesit = thisSim%SIM_uses_SF
762 end function SimUsesSF
763
764 function SimUsesSP() result(doesit)
765 logical :: doesit
766 doesit = thisSim%SIM_uses_SP
767 end function SimUsesSP
768
769 function SimUsesBoxDipole() result(doesit)
770 logical :: doesit
771 doesit = thisSim%SIM_uses_BoxDipole
772 end function SimUsesBoxDipole
773
774 function SimRequiresPrepairCalc() result(doesit)
775 logical :: doesit
776 doesit = thisSim%SIM_uses_EAM .or. thisSim%SIM_uses_SC
777 end function SimRequiresPrepairCalc
778
779 function SimRequiresPostpairCalc() result(doesit)
780 logical :: doesit
781 doesit = thisSim%SIM_uses_RF .or. thisSim%SIM_uses_SF &
782 .or. thisSim%SIM_uses_SP .or. thisSim%SIM_uses_BoxDipole
783 end function SimRequiresPostpairCalc
784
785 ! Function returns true if the simulation has this atype
786 function SimHasAtype(thisAtype) result(doesit)
787 logical :: doesit
788 integer :: thisAtype
789 doesit = .false.
790 if(.not.allocated(SimHasAtypeMap)) return
791
792 doesit = SimHasAtypeMap(thisAtype)
793
794 end function SimHasAtype
795
796 subroutine createSimHasAtype(status)
797 integer, intent(out) :: status
798 integer :: alloc_stat
799 integer :: me_i
800 integer :: mpiErrors
801 integer :: nAtypes
802 status = 0
803
804 nAtypes = getSize(atypes)
805 ! Setup logical map for atypes in simulation
806 if (.not.allocated(SimHasAtypeMap)) then
807 allocate(SimHasAtypeMap(nAtypes),stat=alloc_stat)
808 if (alloc_stat /= 0 ) then
809 status = -1
810 return
811 end if
812 SimHasAtypeMap = .false.
813 end if
814
815 ! Loop through the local atoms and grab the atypes present
816 do me_i = 1,nLocal
817 SimHasAtypeMap(atid(me_i)) = .true.
818 end do
819 ! For MPI, we need to know all possible atypes present in
820 ! simulation on all processors. Use LOR operation to set map.
821 #ifdef IS_MPI
822 if (.not.allocated(SimHasAtypeMapTemp)) then
823 allocate(SimHasAtypeMapTemp(nAtypes),stat=alloc_stat)
824 if (alloc_stat /= 0 ) then
825 status = -1
826 return
827 end if
828 end if
829 call mpi_allreduce(SimHasAtypeMap, SimHasAtypeMaptemp, nAtypes, &
830 mpi_logical, MPI_LOR, mpi_comm_world, mpiErrors)
831 simHasAtypeMap = simHasAtypeMapTemp
832 deallocate(simHasAtypeMapTemp)
833 #endif
834 end subroutine createSimHasAtype
835
836 subroutine InitializeSimGlobals(thisStat)
837 integer, intent(out) :: thisStat
838 integer :: alloc_stat
839
840 thisStat = 0
841
842 call FreeSimGlobals()
843
844 allocate(excludes(2,nExcludes), stat=alloc_stat)
845 if (alloc_stat /= 0 ) then
846 thisStat = -1
847 return
848 endif
849
850 allocate(molMembershipList(nGlobal), stat=alloc_stat)
851 if (alloc_stat /= 0 ) then
852 thisStat = -1
853 return
854 endif
855
856 end subroutine InitializeSimGlobals
857
858 subroutine FreeSimGlobals()
859
860 !We free in the opposite order in which we allocate in.
861 if (allocated(topoDistance)) deallocate(topoDistance)
862 if (allocated(toposForAtom)) deallocate(toposForAtom)
863 if (allocated(nTopoPairsForAtom)) deallocate(nTopoPairsForAtom)
864 if (allocated(skipsForAtom)) deallocate(skipsForAtom)
865 if (allocated(nSkipsForAtom)) deallocate(nSkipsForAtom)
866 if (allocated(skipsForAtom)) deallocate(skipsForAtom)
867 if (allocated(nSkipsForAtom)) deallocate(nSkipsForAtom)
868 if (allocated(mfactLocal)) deallocate(mfactLocal)
869 if (allocated(mfactCol)) deallocate(mfactCol)
870 if (allocated(mfactRow)) deallocate(mfactRow)
871 if (allocated(groupListCol)) deallocate(groupListCol)
872 if (allocated(groupListRow)) deallocate(groupListRow)
873 if (allocated(groupStartCol)) deallocate(groupStartCol)
874 if (allocated(groupStartRow)) deallocate(groupStartRow)
875 if (allocated(molMembershipList)) deallocate(molMembershipList)
876 if (allocated(excludes)) deallocate(excludes)
877
878 end subroutine FreeSimGlobals
879
880 pure function getNlocal() result(n)
881 integer :: n
882 n = nLocal
883 end function getNlocal
884
885 subroutine setHmatDangerousRcutValue(dangerWillRobinson)
886 real(kind=dp), intent(in) :: dangerWillRobinson
887 DangerRcut = dangerWillRobinson
888
889 call checkBox()
890
891 return
892 end subroutine setHmatDangerousRcutValue
893
894
895
896 end module simulation