ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/simulation.F90
Revision: 1346
Committed: Tue May 19 20:21:59 2009 UTC (15 years, 11 months ago) by gezelter
File size: 29048 byte(s)
Log Message:
Fixing a parallel bug in the self-self interaction.

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