ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/simulation.F90
Revision: 998
Committed: Mon Jul 3 13:18:43 2006 UTC (18 years, 10 months ago) by chrisfen
File size: 23628 byte(s)
Log Message:
Added simulation box dipole moment accumulation for the purposes of calculating dielectric constants

File Contents

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