ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/simulation.F90
Revision: 3020
Committed: Fri Sep 22 22:19:59 2006 UTC (17 years, 11 months ago) by chrisfen
File size: 23451 byte(s)
Log Message:
refined the reaction field dielectric passing

File Contents

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