ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2301
Committed: Thu Sep 15 22:05:21 2005 UTC (18 years, 10 months ago) by gezelter
File size: 42297 byte(s)
Log Message:
Working on adding WOLF

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 1610 !! doForces.F90
43     !! module doForces
44     !! Calculates Long Range forces.
45    
46     !! @author Charles F. Vardeman II
47     !! @author Matthew Meineke
48 gezelter 2301 !! @version $Id: doForces.F90,v 1.44 2005-09-15 22:05:17 gezelter Exp $, $Date: 2005-09-15 22:05:17 $, $Name: not supported by cvs2svn $, $Revision: 1.44 $
49 gezelter 1610
50 gezelter 1930
51 gezelter 1610 module doForces
52     use force_globals
53     use simulation
54     use definitions
55     use atype_module
56     use switcheroo
57     use neighborLists
58     use lj
59 gezelter 1930 use sticky
60 gezelter 2085 use electrostatic_module
61 chrisfen 2295 use reaction_field_module
62 gezelter 1610 use gb_pair
63 chrisfen 1636 use shapes
64 gezelter 1610 use vector_class
65     use eam
66     use status
67     #ifdef IS_MPI
68     use mpiSimulation
69     #endif
70    
71     implicit none
72     PRIVATE
73    
74     #define __FORTRAN90
75     #include "UseTheForce/fSwitchingFunction.h"
76 gezelter 2273 #include "UseTheForce/fCutoffPolicy.h"
77 gezelter 2259 #include "UseTheForce/DarkSide/fInteractionMap.h"
78 gezelter 1610
79 gezelter 2273
80 gezelter 1610 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
81     INTEGER, PARAMETER:: PAIR_LOOP = 2
82    
83     logical, save :: haveNeighborList = .false.
84     logical, save :: haveSIMvariables = .false.
85     logical, save :: haveSaneForceField = .false.
86 gezelter 2270 logical, save :: haveInteractionHash = .false.
87     logical, save :: haveGtypeCutoffMap = .false.
88 gezelter 2280 logical, save :: haveRlist = .false.
89 chrisfen 2229
90 gezelter 1634 logical, save :: FF_uses_DirectionalAtoms
91 gezelter 2085 logical, save :: FF_uses_Dipoles
92 gezelter 1634 logical, save :: FF_uses_GayBerne
93     logical, save :: FF_uses_EAM
94 gezelter 1610 logical, save :: FF_uses_RF
95 gezelter 1634
96     logical, save :: SIM_uses_DirectionalAtoms
97     logical, save :: SIM_uses_EAM
98 gezelter 1610 logical, save :: SIM_uses_RF
99     logical, save :: SIM_requires_postpair_calc
100     logical, save :: SIM_requires_prepair_calc
101     logical, save :: SIM_uses_PBC
102    
103 chrisfen 2279 integer, save :: corrMethod
104    
105 gezelter 1610 public :: init_FF
106 gezelter 2273 public :: setDefaultCutoffs
107 gezelter 1610 public :: do_force_loop
108 gezelter 2270 public :: createInteractionHash
109     public :: createGtypeCutoffMap
110 chrisfen 2277 public :: getStickyCut
111     public :: getStickyPowerCut
112     public :: getGayBerneCut
113     public :: getEAMCut
114     public :: getShapeCut
115 gezelter 1610
116     #ifdef PROFILE
117     public :: getforcetime
118     real, save :: forceTime = 0
119     real :: forceTimeInitial, forceTimeFinal
120     integer :: nLoops
121     #endif
122 chuckv 2260
123 gezelter 2270 !! Variables for cutoff mapping and interaction mapping
124     ! Bit hash to determine pair-pair interactions.
125     integer, dimension(:,:), allocatable :: InteractionHash
126     real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff
127 chuckv 2269 real(kind=dp), dimension(:), allocatable :: groupMaxCutoff
128 gezelter 2270 integer, dimension(:), allocatable :: groupToGtype
129 chuckv 2269 real(kind=dp), dimension(:), allocatable :: gtypeMaxCutoff
130 gezelter 2270 type ::gtypeCutoffs
131     real(kind=dp) :: rcut
132     real(kind=dp) :: rcutsq
133     real(kind=dp) :: rlistsq
134     end type gtypeCutoffs
135     type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap
136 gezelter 2273
137     integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
138     real(kind=dp),save :: defaultRcut, defaultRsw, defaultRlist
139 chrisfen 2295 real(kind=dp),save :: rcuti
140 gezelter 2261
141 gezelter 1610 contains
142    
143 gezelter 2270 subroutine createInteractionHash(status)
144 chuckv 2260 integer :: nAtypes
145 chuckv 2266 integer, intent(out) :: status
146 chuckv 2260 integer :: i
147     integer :: j
148 gezelter 2270 integer :: iHash
149 tim 2267 !! Test Types
150 chuckv 2260 logical :: i_is_LJ
151     logical :: i_is_Elect
152     logical :: i_is_Sticky
153     logical :: i_is_StickyP
154     logical :: i_is_GB
155     logical :: i_is_EAM
156     logical :: i_is_Shape
157     logical :: j_is_LJ
158     logical :: j_is_Elect
159     logical :: j_is_Sticky
160     logical :: j_is_StickyP
161     logical :: j_is_GB
162     logical :: j_is_EAM
163     logical :: j_is_Shape
164 gezelter 2275 real(kind=dp) :: myRcut
165    
166 gezelter 2268 status = 0
167    
168 chuckv 2260 if (.not. associated(atypes)) then
169 gezelter 2270 call handleError("atype", "atypes was not present before call of createInteractionHash!")
170 chuckv 2260 status = -1
171     return
172     endif
173    
174     nAtypes = getSize(atypes)
175    
176     if (nAtypes == 0) then
177     status = -1
178     return
179     end if
180 chrisfen 2229
181 chuckv 2269 if (.not. allocated(InteractionHash)) then
182     allocate(InteractionHash(nAtypes,nAtypes))
183 chuckv 2260 endif
184 gezelter 2270
185     if (.not. allocated(atypeMaxCutoff)) then
186     allocate(atypeMaxCutoff(nAtypes))
187     endif
188 chuckv 2260
189     do i = 1, nAtypes
190     call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
191     call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
192     call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
193     call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
194     call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
195     call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
196     call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
197 gezelter 1610
198 chuckv 2260 do j = i, nAtypes
199 chrisfen 2229
200 chuckv 2260 iHash = 0
201     myRcut = 0.0_dp
202 gezelter 1610
203 chuckv 2260 call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ)
204     call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect)
205     call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky)
206     call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP)
207     call getElementProperty(atypes, j, "is_GayBerne", j_is_GB)
208     call getElementProperty(atypes, j, "is_EAM", j_is_EAM)
209     call getElementProperty(atypes, j, "is_Shape", j_is_Shape)
210 gezelter 1610
211 chuckv 2260 if (i_is_LJ .and. j_is_LJ) then
212 gezelter 2261 iHash = ior(iHash, LJ_PAIR)
213     endif
214    
215     if (i_is_Elect .and. j_is_Elect) then
216     iHash = ior(iHash, ELECTROSTATIC_PAIR)
217     endif
218    
219     if (i_is_Sticky .and. j_is_Sticky) then
220     iHash = ior(iHash, STICKY_PAIR)
221     endif
222 chuckv 2260
223 gezelter 2261 if (i_is_StickyP .and. j_is_StickyP) then
224     iHash = ior(iHash, STICKYPOWER_PAIR)
225     endif
226 chuckv 2260
227 gezelter 2261 if (i_is_EAM .and. j_is_EAM) then
228     iHash = ior(iHash, EAM_PAIR)
229 chuckv 2260 endif
230    
231     if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR)
232     if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ)
233     if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ)
234    
235     if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR)
236     if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ)
237     if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ)
238    
239    
240 chuckv 2269 InteractionHash(i,j) = iHash
241     InteractionHash(j,i) = iHash
242 chuckv 2260
243     end do
244    
245     end do
246 tim 2267
247 gezelter 2270 haveInteractionHash = .true.
248     end subroutine createInteractionHash
249 chuckv 2260
250 gezelter 2275 subroutine createGtypeCutoffMap(stat)
251 gezelter 2268
252 gezelter 2275 integer, intent(out), optional :: stat
253 gezelter 2273 logical :: i_is_LJ
254     logical :: i_is_Elect
255     logical :: i_is_Sticky
256     logical :: i_is_StickyP
257     logical :: i_is_GB
258     logical :: i_is_EAM
259     logical :: i_is_Shape
260 gezelter 2286 logical :: GtypeFound
261 chuckv 2260
262 gezelter 2275 integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend
263 gezelter 2280 integer :: n_in_i, me_i, ia, g, atom1, nGroupTypes
264 chuckv 2288 integer :: nGroupsInRow
265 gezelter 2280 real(kind=dp):: thisSigma, bigSigma, thisRcut, tol, skin
266 gezelter 2275 real(kind=dp) :: biggestAtypeCutoff
267 gezelter 2270
268 chuckv 2266 stat = 0
269 gezelter 2270 if (.not. haveInteractionHash) then
270     call createInteractionHash(myStatus)
271 chuckv 2266 if (myStatus .ne. 0) then
272 gezelter 2270 write(default_error, *) 'createInteractionHash failed in doForces!'
273 chuckv 2266 stat = -1
274     return
275     endif
276     endif
277 chuckv 2288 #ifdef IS_MPI
278     nGroupsInRow = getNgroupsInRow(plan_group_row)
279     #endif
280 chuckv 2262 nAtypes = getSize(atypes)
281 chuckv 2298 ! Set all of the initial cutoffs to zero.
282     atypeMaxCutoff = 0.0_dp
283 gezelter 2270 do i = 1, nAtypes
284 gezelter 2281 if (SimHasAtype(i)) then
285 gezelter 2274 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
286     call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
287     call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
288     call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
289     call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
290     call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
291     call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
292    
293 chuckv 2298
294 gezelter 2274 if (i_is_LJ) then
295     thisRcut = getSigma(i) * 2.5_dp
296     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
297     endif
298     if (i_is_Elect) then
299     thisRcut = defaultRcut
300     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
301     endif
302     if (i_is_Sticky) then
303     thisRcut = getStickyCut(i)
304     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
305     endif
306 gezelter 2275 if (i_is_StickyP) then
307 gezelter 2274 thisRcut = getStickyPowerCut(i)
308     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
309     endif
310 gezelter 2275 if (i_is_GB) then
311 gezelter 2274 thisRcut = getGayBerneCut(i)
312     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
313     endif
314     if (i_is_EAM) then
315     thisRcut = getEAMCut(i)
316     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
317     endif
318     if (i_is_Shape) then
319     thisRcut = getShapeCut(i)
320     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
321     endif
322    
323     if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
324     biggestAtypeCutoff = atypeMaxCutoff(i)
325     endif
326 gezelter 2273 endif
327 gezelter 2274 enddo
328 gezelter 2280
329     nGroupTypes = 0
330    
331 gezelter 2274 istart = 1
332     #ifdef IS_MPI
333     iend = nGroupsInRow
334     #else
335     iend = nGroups
336     #endif
337 gezelter 2281
338 gezelter 2280 !! allocate the groupToGtype and gtypeMaxCutoff here.
339 chuckv 2282 if(.not.allocated(groupToGtype)) then
340     allocate(groupToGtype(iend))
341     allocate(groupMaxCutoff(iend))
342     allocate(gtypeMaxCutoff(iend))
343 chuckv 2298 groupMaxCutoff = 0.0_dp
344     gtypeMaxCutoff = 0.0_dp
345 chuckv 2282 endif
346 gezelter 2281 !! first we do a single loop over the cutoff groups to find the
347     !! largest cutoff for any atypes present in this group. We also
348     !! create gtypes at this point.
349    
350 gezelter 2280 tol = 1.0d-6
351    
352     do i = istart, iend
353 gezelter 2274 n_in_i = groupStartRow(i+1) - groupStartRow(i)
354 gezelter 2280 groupMaxCutoff(i) = 0.0_dp
355     do ia = groupStartRow(i), groupStartRow(i+1)-1
356     atom1 = groupListRow(ia)
357 gezelter 2274 #ifdef IS_MPI
358 gezelter 2280 me_i = atid_row(atom1)
359 gezelter 2274 #else
360 gezelter 2280 me_i = atid(atom1)
361     #endif
362     if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then
363     groupMaxCutoff(i)=atypeMaxCutoff(me_i)
364 gezelter 2286 endif
365 gezelter 2280 enddo
366 gezelter 2286
367 gezelter 2280 if (nGroupTypes.eq.0) then
368     nGroupTypes = nGroupTypes + 1
369     gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
370     groupToGtype(i) = nGroupTypes
371     else
372 gezelter 2286 GtypeFound = .false.
373 gezelter 2280 do g = 1, nGroupTypes
374 gezelter 2286 if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).lt.tol) then
375 gezelter 2280 groupToGtype(i) = g
376 gezelter 2286 GtypeFound = .true.
377 gezelter 2280 endif
378     enddo
379 gezelter 2286 if (.not.GtypeFound) then
380     nGroupTypes = nGroupTypes + 1
381     gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
382     groupToGtype(i) = nGroupTypes
383     endif
384 gezelter 2280 endif
385 gezelter 2286 enddo
386    
387 gezelter 2280 !! allocate the gtypeCutoffMap here.
388 chuckv 2282 allocate(gtypeCutoffMap(nGroupTypes,nGroupTypes))
389 gezelter 2280 !! then we do a double loop over all the group TYPES to find the cutoff
390     !! map between groups of two types
391    
392     do i = 1, nGroupTypes
393     do j = 1, nGroupTypes
394 gezelter 2275
395 gezelter 2280 select case(cutoffPolicy)
396 gezelter 2281 case(TRADITIONAL_CUTOFF_POLICY)
397     thisRcut = maxval(gtypeMaxCutoff)
398     case(MIX_CUTOFF_POLICY)
399     thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j))
400     case(MAX_CUTOFF_POLICY)
401     thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j))
402     case default
403     call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
404     return
405     end select
406     gtypeCutoffMap(i,j)%rcut = thisRcut
407     gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
408     skin = defaultRlist - defaultRcut
409     gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2
410 gezelter 2284
411 gezelter 2280 enddo
412     enddo
413 gezelter 2275
414 gezelter 2280 haveGtypeCutoffMap = .true.
415 chrisfen 2295 end subroutine createGtypeCutoffMap
416 chrisfen 2277
417 chrisfen 2295 subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
418     real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
419 gezelter 2273 integer, intent(in) :: cutPolicy
420 chrisfen 2295
421     defaultRcut = defRcut
422     defaultRsw = defRsw
423     defaultRlist = defRlist
424 gezelter 2273 cutoffPolicy = cutPolicy
425 chrisfen 2295 rcuti = 1.0_dp / defaultRcut
426     end subroutine setDefaultCutoffs
427    
428     subroutine setCutoffPolicy(cutPolicy)
429    
430     integer, intent(in) :: cutPolicy
431     cutoffPolicy = cutPolicy
432 gezelter 2273 call createGtypeCutoffMap()
433 gezelter 2275 end subroutine setCutoffPolicy
434 gezelter 2273
435    
436 gezelter 1610 subroutine setSimVariables()
437 gezelter 1634 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
438     SIM_uses_EAM = SimUsesEAM()
439 gezelter 1610 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
440     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
441     SIM_uses_PBC = SimUsesPBC()
442 chrisfen 2295 SIM_uses_RF = SimUsesRF()
443 gezelter 1610
444     haveSIMvariables = .true.
445    
446     return
447     end subroutine setSimVariables
448    
449     subroutine doReadyCheck(error)
450     integer, intent(out) :: error
451    
452     integer :: myStatus
453    
454     error = 0
455 chrisfen 2229
456 gezelter 2270 if (.not. haveInteractionHash) then
457 gezelter 2268 myStatus = 0
458 gezelter 2270 call createInteractionHash(myStatus)
459 gezelter 1610 if (myStatus .ne. 0) then
460 gezelter 2270 write(default_error, *) 'createInteractionHash failed in doForces!'
461 gezelter 1610 error = -1
462     return
463     endif
464     endif
465    
466 gezelter 2270 if (.not. haveGtypeCutoffMap) then
467     myStatus = 0
468     call createGtypeCutoffMap(myStatus)
469     if (myStatus .ne. 0) then
470     write(default_error, *) 'createGtypeCutoffMap failed in doForces!'
471     error = -1
472     return
473     endif
474     endif
475    
476 gezelter 1610 if (.not. haveSIMvariables) then
477     call setSimVariables()
478     endif
479    
480 chuckv 2282 ! if (.not. haveRlist) then
481     ! write(default_error, *) 'rList has not been set in doForces!'
482     ! error = -1
483     ! return
484     ! endif
485 gezelter 1610
486     if (.not. haveNeighborList) then
487     write(default_error, *) 'neighbor list has not been initialized in doForces!'
488     error = -1
489     return
490     end if
491    
492     if (.not. haveSaneForceField) then
493     write(default_error, *) 'Force Field is not sane in doForces!'
494     error = -1
495     return
496     end if
497    
498     #ifdef IS_MPI
499     if (.not. isMPISimSet()) then
500     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
501     error = -1
502     return
503     endif
504     #endif
505     return
506     end subroutine doReadyCheck
507    
508 chrisfen 2229
509 chrisfen 2295 subroutine init_FF(use_RF, correctionMethod, dampingAlpha, thisStat)
510 gezelter 1610
511 gezelter 2280 logical, intent(in) :: use_RF
512 chrisfen 2295 integer, intent(in) :: correctionMethod
513     real(kind=dp), intent(in) :: dampingAlpha
514 gezelter 1610 integer, intent(out) :: thisStat
515     integer :: my_status, nMatches
516     integer, pointer :: MatchList(:) => null()
517     real(kind=dp) :: rcut, rrf, rt, dielect
518    
519     !! assume things are copacetic, unless they aren't
520     thisStat = 0
521    
522     !! Fortran's version of a cast:
523 gezelter 2280 FF_uses_RF = use_RF
524 chrisfen 2229
525 chrisfen 2295
526 gezelter 1610 !! init_FF is called *after* all of the atom types have been
527     !! defined in atype_module using the new_atype subroutine.
528     !!
529     !! this will scan through the known atypes and figure out what
530     !! interactions are used by the force field.
531 chrisfen 2229
532 gezelter 1634 FF_uses_DirectionalAtoms = .false.
533     FF_uses_Dipoles = .false.
534     FF_uses_GayBerne = .false.
535 gezelter 1610 FF_uses_EAM = .false.
536 chrisfen 2229
537 gezelter 1634 call getMatchingElementList(atypes, "is_Directional", .true., &
538     nMatches, MatchList)
539     if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
540    
541     call getMatchingElementList(atypes, "is_Dipole", .true., &
542     nMatches, MatchList)
543 gezelter 2270 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
544 chrisfen 2220
545 gezelter 1634 call getMatchingElementList(atypes, "is_GayBerne", .true., &
546     nMatches, MatchList)
547 gezelter 2270 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
548 chrisfen 2229
549 gezelter 1610 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
550     if (nMatches .gt. 0) FF_uses_EAM = .true.
551 chrisfen 2229
552 gezelter 1634
553 gezelter 1610 haveSaneForceField = .true.
554 chrisfen 2229
555 gezelter 1610 !! check to make sure the FF_uses_RF setting makes sense
556 chrisfen 2229
557 gezelter 2270 if (FF_uses_Dipoles) then
558 gezelter 1610 if (FF_uses_RF) then
559     dielect = getDielect()
560     call initialize_rf(dielect)
561     endif
562     else
563 chrisfen 2295 if ((corrMethod == 3) .or. FF_uses_RF) then
564 gezelter 1610 write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
565     thisStat = -1
566     haveSaneForceField = .false.
567     return
568     endif
569 chrisfen 2229 endif
570 gezelter 1610
571     if (FF_uses_EAM) then
572 chrisfen 2229 call init_EAM_FF(my_status)
573 gezelter 1610 if (my_status /= 0) then
574     write(default_error, *) "init_EAM_FF returned a bad status"
575     thisStat = -1
576     haveSaneForceField = .false.
577     return
578     end if
579     endif
580    
581 gezelter 1634 if (FF_uses_GayBerne) then
582 gezelter 1610 call check_gb_pair_FF(my_status)
583     if (my_status .ne. 0) then
584     thisStat = -1
585     haveSaneForceField = .false.
586     return
587     endif
588     endif
589    
590     if (.not. haveNeighborList) then
591     !! Create neighbor lists
592     call expandNeighborList(nLocal, my_status)
593     if (my_Status /= 0) then
594     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
595     thisStat = -1
596     return
597     endif
598     haveNeighborList = .true.
599 chrisfen 2229 endif
600    
601 gezelter 1610 end subroutine init_FF
602    
603 chrisfen 2229
604 gezelter 1610 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
605     !------------------------------------------------------------->
606 gezelter 1930 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
607 gezelter 1610 do_pot_c, do_stress_c, error)
608     !! Position array provided by C, dimensioned by getNlocal
609     real ( kind = dp ), dimension(3, nLocal) :: q
610     !! molecular center-of-mass position array
611     real ( kind = dp ), dimension(3, nGroups) :: q_group
612     !! Rotation Matrix for each long range particle in simulation.
613     real( kind = dp), dimension(9, nLocal) :: A
614     !! Unit vectors for dipoles (lab frame)
615 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
616 gezelter 1610 !! Force array provided by C, dimensioned by getNlocal
617     real ( kind = dp ), dimension(3,nLocal) :: f
618     !! Torsion array provided by C, dimensioned by getNlocal
619     real( kind = dp ), dimension(3,nLocal) :: t
620    
621     !! Stress Tensor
622     real( kind = dp), dimension(9) :: tau
623     real ( kind = dp ) :: pot
624     logical ( kind = 2) :: do_pot_c, do_stress_c
625     logical :: do_pot
626     logical :: do_stress
627     logical :: in_switching_region
628     #ifdef IS_MPI
629     real( kind = DP ) :: pot_local
630     integer :: nAtomsInRow
631     integer :: nAtomsInCol
632     integer :: nprocs
633     integer :: nGroupsInRow
634     integer :: nGroupsInCol
635     #endif
636     integer :: natoms
637     logical :: update_nlist
638     integer :: i, j, jstart, jend, jnab
639     integer :: istart, iend
640     integer :: ia, jb, atom1, atom2
641     integer :: nlist
642     real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
643     real( kind = DP ) :: sw, dswdr, swderiv, mf
644     real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij
645     real(kind=dp) :: rfpot, mu_i, virial
646     integer :: me_i, me_j, n_in_i, n_in_j
647     logical :: is_dp_i
648     integer :: neighborListSize
649     integer :: listerror, error
650     integer :: localError
651     integer :: propPack_i, propPack_j
652     integer :: loopStart, loopEnd, loop
653 gezelter 2270 integer :: iHash
654 gezelter 1610 real(kind=dp) :: listSkin = 1.0
655 chrisfen 2229
656 gezelter 1610 !! initialize local variables
657 chrisfen 2229
658 gezelter 1610 #ifdef IS_MPI
659     pot_local = 0.0_dp
660     nAtomsInRow = getNatomsInRow(plan_atom_row)
661     nAtomsInCol = getNatomsInCol(plan_atom_col)
662     nGroupsInRow = getNgroupsInRow(plan_group_row)
663     nGroupsInCol = getNgroupsInCol(plan_group_col)
664     #else
665     natoms = nlocal
666     #endif
667 chrisfen 2229
668 gezelter 1610 call doReadyCheck(localError)
669     if ( localError .ne. 0 ) then
670     call handleError("do_force_loop", "Not Initialized")
671     error = -1
672     return
673     end if
674     call zero_work_arrays()
675 chrisfen 2229
676 gezelter 1610 do_pot = do_pot_c
677     do_stress = do_stress_c
678 chrisfen 2229
679 gezelter 1610 ! Gather all information needed by all force loops:
680 chrisfen 2229
681 gezelter 1610 #ifdef IS_MPI
682 chrisfen 2229
683 gezelter 1610 call gather(q, q_Row, plan_atom_row_3d)
684     call gather(q, q_Col, plan_atom_col_3d)
685    
686     call gather(q_group, q_group_Row, plan_group_row_3d)
687     call gather(q_group, q_group_Col, plan_group_col_3d)
688 chrisfen 2229
689 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
690 gezelter 1930 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
691     call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
692 chrisfen 2229
693 gezelter 1610 call gather(A, A_Row, plan_atom_row_rotation)
694     call gather(A, A_Col, plan_atom_col_rotation)
695     endif
696 chrisfen 2229
697 gezelter 1610 #endif
698 chrisfen 2229
699 gezelter 1610 !! Begin force loop timing:
700     #ifdef PROFILE
701     call cpu_time(forceTimeInitial)
702     nloops = nloops + 1
703     #endif
704 chrisfen 2229
705 gezelter 1610 loopEnd = PAIR_LOOP
706     if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
707     loopStart = PREPAIR_LOOP
708     else
709     loopStart = PAIR_LOOP
710     endif
711    
712     do loop = loopStart, loopEnd
713    
714     ! See if we need to update neighbor lists
715     ! (but only on the first time through):
716     if (loop .eq. loopStart) then
717     #ifdef IS_MPI
718     call checkNeighborList(nGroupsInRow, q_group_row, listSkin, &
719 chrisfen 2229 update_nlist)
720 gezelter 1610 #else
721     call checkNeighborList(nGroups, q_group, listSkin, &
722 chrisfen 2229 update_nlist)
723 gezelter 1610 #endif
724     endif
725 chrisfen 2229
726 gezelter 1610 if (update_nlist) then
727     !! save current configuration and construct neighbor list
728     #ifdef IS_MPI
729     call saveNeighborList(nGroupsInRow, q_group_row)
730     #else
731     call saveNeighborList(nGroups, q_group)
732     #endif
733     neighborListSize = size(list)
734     nlist = 0
735     endif
736 chrisfen 2229
737 gezelter 1610 istart = 1
738     #ifdef IS_MPI
739     iend = nGroupsInRow
740     #else
741     iend = nGroups - 1
742     #endif
743     outer: do i = istart, iend
744    
745     if (update_nlist) point(i) = nlist + 1
746 chrisfen 2229
747 gezelter 1610 n_in_i = groupStartRow(i+1) - groupStartRow(i)
748 chrisfen 2229
749 gezelter 1610 if (update_nlist) then
750     #ifdef IS_MPI
751     jstart = 1
752     jend = nGroupsInCol
753     #else
754     jstart = i+1
755     jend = nGroups
756     #endif
757     else
758     jstart = point(i)
759     jend = point(i+1) - 1
760     ! make sure group i has neighbors
761     if (jstart .gt. jend) cycle outer
762     endif
763 chrisfen 2229
764 gezelter 1610 do jnab = jstart, jend
765     if (update_nlist) then
766     j = jnab
767     else
768     j = list(jnab)
769     endif
770    
771     #ifdef IS_MPI
772 chuckv 2266 me_j = atid_col(j)
773 gezelter 1610 call get_interatomic_vector(q_group_Row(:,i), &
774     q_group_Col(:,j), d_grp, rgrpsq)
775     #else
776 chuckv 2266 me_j = atid(j)
777 gezelter 1610 call get_interatomic_vector(q_group(:,i), &
778     q_group(:,j), d_grp, rgrpsq)
779     #endif
780    
781 gezelter 2280 if (rgrpsq < gtypeCutoffMap(groupToGtype(i),groupToGtype(j))%rListsq) then
782 gezelter 1610 if (update_nlist) then
783     nlist = nlist + 1
784 chrisfen 2229
785 gezelter 1610 if (nlist > neighborListSize) then
786     #ifdef IS_MPI
787     call expandNeighborList(nGroupsInRow, listerror)
788     #else
789     call expandNeighborList(nGroups, listerror)
790     #endif
791     if (listerror /= 0) then
792     error = -1
793     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
794     return
795     end if
796     neighborListSize = size(list)
797     endif
798 chrisfen 2229
799 gezelter 1610 list(nlist) = j
800     endif
801 chrisfen 2229
802 gezelter 1610 if (loop .eq. PAIR_LOOP) then
803     vij = 0.0d0
804     fij(1:3) = 0.0d0
805     endif
806 chrisfen 2229
807 gezelter 1610 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
808     in_switching_region)
809 chrisfen 2229
810 gezelter 1610 n_in_j = groupStartCol(j+1) - groupStartCol(j)
811 chrisfen 2229
812 gezelter 1610 do ia = groupStartRow(i), groupStartRow(i+1)-1
813 chrisfen 2229
814 gezelter 1610 atom1 = groupListRow(ia)
815 chrisfen 2229
816 gezelter 1610 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
817 chrisfen 2229
818 gezelter 1610 atom2 = groupListCol(jb)
819 chrisfen 2229
820 gezelter 1610 if (skipThisPair(atom1, atom2)) cycle inner
821    
822     if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
823     d_atm(1:3) = d_grp(1:3)
824     ratmsq = rgrpsq
825     else
826     #ifdef IS_MPI
827     call get_interatomic_vector(q_Row(:,atom1), &
828     q_Col(:,atom2), d_atm, ratmsq)
829     #else
830     call get_interatomic_vector(q(:,atom1), &
831     q(:,atom2), d_atm, ratmsq)
832     #endif
833     endif
834    
835     if (loop .eq. PREPAIR_LOOP) then
836     #ifdef IS_MPI
837     call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
838     rgrpsq, d_grp, do_pot, do_stress, &
839 gezelter 1930 eFrame, A, f, t, pot_local)
840 gezelter 1610 #else
841     call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
842     rgrpsq, d_grp, do_pot, do_stress, &
843 gezelter 1930 eFrame, A, f, t, pot)
844 gezelter 1610 #endif
845     else
846     #ifdef IS_MPI
847     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
848     do_pot, &
849 gezelter 1930 eFrame, A, f, t, pot_local, vpair, fpair)
850 gezelter 1610 #else
851     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
852     do_pot, &
853 gezelter 1930 eFrame, A, f, t, pot, vpair, fpair)
854 gezelter 1610 #endif
855    
856     vij = vij + vpair
857     fij(1:3) = fij(1:3) + fpair(1:3)
858     endif
859     enddo inner
860     enddo
861 chrisfen 2229
862 gezelter 1610 if (loop .eq. PAIR_LOOP) then
863     if (in_switching_region) then
864     swderiv = vij*dswdr/rgrp
865     fij(1) = fij(1) + swderiv*d_grp(1)
866     fij(2) = fij(2) + swderiv*d_grp(2)
867     fij(3) = fij(3) + swderiv*d_grp(3)
868 chrisfen 2229
869 gezelter 1610 do ia=groupStartRow(i), groupStartRow(i+1)-1
870     atom1=groupListRow(ia)
871     mf = mfactRow(atom1)
872     #ifdef IS_MPI
873     f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
874     f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
875     f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
876     #else
877     f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
878     f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
879     f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
880     #endif
881     enddo
882 chrisfen 2229
883 gezelter 1610 do jb=groupStartCol(j), groupStartCol(j+1)-1
884     atom2=groupListCol(jb)
885     mf = mfactCol(atom2)
886     #ifdef IS_MPI
887     f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
888     f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
889     f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
890     #else
891     f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
892     f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
893     f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
894     #endif
895     enddo
896     endif
897 chrisfen 2229
898 gezelter 1610 if (do_stress) call add_stress_tensor(d_grp, fij)
899     endif
900     end if
901     enddo
902     enddo outer
903 chrisfen 2229
904 gezelter 1610 if (update_nlist) then
905     #ifdef IS_MPI
906     point(nGroupsInRow + 1) = nlist + 1
907     #else
908     point(nGroups) = nlist + 1
909     #endif
910     if (loop .eq. PREPAIR_LOOP) then
911     ! we just did the neighbor list update on the first
912     ! pass, so we don't need to do it
913     ! again on the second pass
914     update_nlist = .false.
915     endif
916     endif
917 chrisfen 2229
918 gezelter 1610 if (loop .eq. PREPAIR_LOOP) then
919     call do_preforce(nlocal, pot)
920     endif
921 chrisfen 2229
922 gezelter 1610 enddo
923 chrisfen 2229
924 gezelter 1610 !! Do timing
925     #ifdef PROFILE
926     call cpu_time(forceTimeFinal)
927     forceTime = forceTime + forceTimeFinal - forceTimeInitial
928     #endif
929 chrisfen 2229
930 gezelter 1610 #ifdef IS_MPI
931     !!distribute forces
932 chrisfen 2229
933 gezelter 1610 f_temp = 0.0_dp
934     call scatter(f_Row,f_temp,plan_atom_row_3d)
935     do i = 1,nlocal
936     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
937     end do
938 chrisfen 2229
939 gezelter 1610 f_temp = 0.0_dp
940     call scatter(f_Col,f_temp,plan_atom_col_3d)
941     do i = 1,nlocal
942     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
943     end do
944 chrisfen 2229
945 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
946 gezelter 1610 t_temp = 0.0_dp
947     call scatter(t_Row,t_temp,plan_atom_row_3d)
948     do i = 1,nlocal
949     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
950     end do
951     t_temp = 0.0_dp
952     call scatter(t_Col,t_temp,plan_atom_col_3d)
953 chrisfen 2229
954 gezelter 1610 do i = 1,nlocal
955     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
956     end do
957     endif
958 chrisfen 2229
959 gezelter 1610 if (do_pot) then
960     ! scatter/gather pot_row into the members of my column
961     call scatter(pot_Row, pot_Temp, plan_atom_row)
962 chrisfen 2229
963 gezelter 1610 ! scatter/gather pot_local into all other procs
964     ! add resultant to get total pot
965     do i = 1, nlocal
966     pot_local = pot_local + pot_Temp(i)
967     enddo
968 chrisfen 2229
969 gezelter 1610 pot_Temp = 0.0_DP
970 chrisfen 2229
971 gezelter 1610 call scatter(pot_Col, pot_Temp, plan_atom_col)
972     do i = 1, nlocal
973     pot_local = pot_local + pot_Temp(i)
974     enddo
975 chrisfen 2229
976 gezelter 1610 endif
977     #endif
978 chrisfen 2229
979 gezelter 1610 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
980 chrisfen 2229
981 chrisfen 2295 if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then
982 chrisfen 2229
983 gezelter 1610 #ifdef IS_MPI
984     call scatter(rf_Row,rf,plan_atom_row_3d)
985     call scatter(rf_Col,rf_Temp,plan_atom_col_3d)
986     do i = 1,nlocal
987     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
988     end do
989     #endif
990 chrisfen 2229
991 gezelter 1610 do i = 1, nLocal
992 chrisfen 2229
993 gezelter 1610 rfpot = 0.0_DP
994     #ifdef IS_MPI
995     me_i = atid_row(i)
996     #else
997     me_i = atid(i)
998     #endif
999 gezelter 2270 iHash = InteractionHash(me_i,me_j)
1000 gezelter 2261
1001 gezelter 2270 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1002 chrisfen 2229
1003 gezelter 1634 mu_i = getDipoleMoment(me_i)
1004 chrisfen 2229
1005 gezelter 1610 !! The reaction field needs to include a self contribution
1006     !! to the field:
1007 gezelter 1930 call accumulate_self_rf(i, mu_i, eFrame)
1008 gezelter 1610 !! Get the reaction field contribution to the
1009     !! potential and torques:
1010 gezelter 1930 call reaction_field_final(i, mu_i, eFrame, rfpot, t, do_pot)
1011 gezelter 1610 #ifdef IS_MPI
1012     pot_local = pot_local + rfpot
1013     #else
1014     pot = pot + rfpot
1015 chrisfen 2229
1016 gezelter 1610 #endif
1017 chrisfen 2229 endif
1018 gezelter 1610 enddo
1019     endif
1020     endif
1021 chrisfen 2229
1022    
1023 gezelter 1610 #ifdef IS_MPI
1024 chrisfen 2229
1025 gezelter 1610 if (do_pot) then
1026     pot = pot + pot_local
1027     !! we assume the c code will do the allreduce to get the total potential
1028     !! we could do it right here if we needed to...
1029     endif
1030 chrisfen 2229
1031 gezelter 1610 if (do_stress) then
1032     call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
1033     mpi_comm_world,mpi_err)
1034     call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
1035     mpi_comm_world,mpi_err)
1036     endif
1037 chrisfen 2229
1038 gezelter 1610 #else
1039 chrisfen 2229
1040 gezelter 1610 if (do_stress) then
1041     tau = tau_Temp
1042     virial = virial_Temp
1043     endif
1044 chrisfen 2229
1045 gezelter 1610 #endif
1046 chrisfen 2229
1047 gezelter 1610 end subroutine do_force_loop
1048 chrisfen 2229
1049 gezelter 1610 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1050 gezelter 1930 eFrame, A, f, t, pot, vpair, fpair)
1051 gezelter 1610
1052     real( kind = dp ) :: pot, vpair, sw
1053     real( kind = dp ), dimension(3) :: fpair
1054     real( kind = dp ), dimension(nLocal) :: mfact
1055 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
1056 gezelter 1610 real( kind = dp ), dimension(9,nLocal) :: A
1057     real( kind = dp ), dimension(3,nLocal) :: f
1058     real( kind = dp ), dimension(3,nLocal) :: t
1059    
1060     logical, intent(inout) :: do_pot
1061     integer, intent(in) :: i, j
1062     real ( kind = dp ), intent(inout) :: rijsq
1063     real ( kind = dp ) :: r
1064     real ( kind = dp ), intent(inout) :: d(3)
1065     integer :: me_i, me_j
1066    
1067 gezelter 2270 integer :: iHash
1068 gezelter 2259
1069 gezelter 1610 r = sqrt(rijsq)
1070     vpair = 0.0d0
1071     fpair(1:3) = 0.0d0
1072    
1073     #ifdef IS_MPI
1074     me_i = atid_row(i)
1075     me_j = atid_col(j)
1076     #else
1077     me_i = atid(i)
1078     me_j = atid(j)
1079     #endif
1080 gezelter 1706
1081 gezelter 2270 iHash = InteractionHash(me_i, me_j)
1082 chrisfen 2229
1083 gezelter 2270 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1084 gezelter 2259 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
1085 gezelter 1610 endif
1086 chrisfen 2229
1087 gezelter 2270 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1088 gezelter 2259 call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, &
1089 chrisfen 2295 pot, eFrame, f, t, do_pot, corrMethod, rcuti)
1090 chrisfen 2229
1091 chrisfen 2295 if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then
1092 gezelter 2261
1093     ! CHECK ME (RF needs to know about all electrostatic types)
1094     call accumulate_rf(i, j, r, eFrame, sw)
1095     call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair)
1096 gezelter 1610 endif
1097 gezelter 2261
1098 gezelter 1610 endif
1099    
1100 gezelter 2270 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1101 gezelter 2259 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1102     pot, A, f, t, do_pot)
1103     endif
1104 gezelter 2085
1105 gezelter 2270 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1106 gezelter 2259 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1107     pot, A, f, t, do_pot)
1108 gezelter 1610 endif
1109    
1110 gezelter 2270 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1111 gezelter 2259 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1112     pot, A, f, t, do_pot)
1113 chrisfen 2229 endif
1114    
1115 gezelter 2270 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1116 chuckv 2262 ! call do_gblj_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1117     ! pot, A, f, t, do_pot)
1118 gezelter 2259 endif
1119 kdaily 2226
1120 gezelter 2270 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1121 gezelter 2259 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
1122     do_pot)
1123 gezelter 1610 endif
1124 chrisfen 2229
1125 gezelter 2270 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1126 gezelter 2259 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1127     pot, A, f, t, do_pot)
1128 gezelter 1610 endif
1129 gezelter 1634
1130 gezelter 2270 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1131 gezelter 2259 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1132     pot, A, f, t, do_pot)
1133 gezelter 1634 endif
1134 gezelter 2259
1135 gezelter 1610 end subroutine do_pair
1136    
1137     subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
1138 gezelter 1930 do_pot, do_stress, eFrame, A, f, t, pot)
1139 gezelter 1610
1140 chrisfen 2229 real( kind = dp ) :: pot, sw
1141     real( kind = dp ), dimension(9,nLocal) :: eFrame
1142     real (kind=dp), dimension(9,nLocal) :: A
1143     real (kind=dp), dimension(3,nLocal) :: f
1144     real (kind=dp), dimension(3,nLocal) :: t
1145 gezelter 1610
1146 chrisfen 2229 logical, intent(inout) :: do_pot, do_stress
1147     integer, intent(in) :: i, j
1148     real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1149     real ( kind = dp ) :: r, rc
1150     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1151    
1152 gezelter 2270 integer :: me_i, me_j, iHash
1153 chrisfen 2229
1154 gezelter 2290 r = sqrt(rijsq)
1155    
1156 gezelter 1610 #ifdef IS_MPI
1157 chrisfen 2229 me_i = atid_row(i)
1158     me_j = atid_col(j)
1159 gezelter 1610 #else
1160 chrisfen 2229 me_i = atid(i)
1161     me_j = atid(j)
1162 gezelter 1610 #endif
1163 chrisfen 2229
1164 gezelter 2270 iHash = InteractionHash(me_i, me_j)
1165 chrisfen 2229
1166 gezelter 2270 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1167 chrisfen 2229 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1168     endif
1169 gezelter 2259
1170 chrisfen 2229 end subroutine do_prepair
1171    
1172    
1173     subroutine do_preforce(nlocal,pot)
1174     integer :: nlocal
1175     real( kind = dp ) :: pot
1176    
1177     if (FF_uses_EAM .and. SIM_uses_EAM) then
1178     call calc_EAM_preforce_Frho(nlocal,pot)
1179     endif
1180    
1181    
1182     end subroutine do_preforce
1183    
1184    
1185     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1186    
1187     real (kind = dp), dimension(3) :: q_i
1188     real (kind = dp), dimension(3) :: q_j
1189     real ( kind = dp ), intent(out) :: r_sq
1190     real( kind = dp ) :: d(3), scaled(3)
1191     integer i
1192    
1193     d(1:3) = q_j(1:3) - q_i(1:3)
1194    
1195     ! Wrap back into periodic box if necessary
1196     if ( SIM_uses_PBC ) then
1197    
1198     if( .not.boxIsOrthorhombic ) then
1199     ! calc the scaled coordinates.
1200    
1201     scaled = matmul(HmatInv, d)
1202    
1203     ! wrap the scaled coordinates
1204    
1205     scaled = scaled - anint(scaled)
1206    
1207    
1208     ! calc the wrapped real coordinates from the wrapped scaled
1209     ! coordinates
1210    
1211     d = matmul(Hmat,scaled)
1212    
1213     else
1214     ! calc the scaled coordinates.
1215    
1216     do i = 1, 3
1217     scaled(i) = d(i) * HmatInv(i,i)
1218    
1219     ! wrap the scaled coordinates
1220    
1221     scaled(i) = scaled(i) - anint(scaled(i))
1222    
1223     ! calc the wrapped real coordinates from the wrapped scaled
1224     ! coordinates
1225    
1226     d(i) = scaled(i)*Hmat(i,i)
1227     enddo
1228     endif
1229    
1230     endif
1231    
1232     r_sq = dot_product(d,d)
1233    
1234     end subroutine get_interatomic_vector
1235    
1236     subroutine zero_work_arrays()
1237    
1238 gezelter 1610 #ifdef IS_MPI
1239    
1240 chrisfen 2229 q_Row = 0.0_dp
1241     q_Col = 0.0_dp
1242    
1243     q_group_Row = 0.0_dp
1244     q_group_Col = 0.0_dp
1245    
1246     eFrame_Row = 0.0_dp
1247     eFrame_Col = 0.0_dp
1248    
1249     A_Row = 0.0_dp
1250     A_Col = 0.0_dp
1251    
1252     f_Row = 0.0_dp
1253     f_Col = 0.0_dp
1254     f_Temp = 0.0_dp
1255    
1256     t_Row = 0.0_dp
1257     t_Col = 0.0_dp
1258     t_Temp = 0.0_dp
1259    
1260     pot_Row = 0.0_dp
1261     pot_Col = 0.0_dp
1262     pot_Temp = 0.0_dp
1263    
1264     rf_Row = 0.0_dp
1265     rf_Col = 0.0_dp
1266     rf_Temp = 0.0_dp
1267    
1268 gezelter 1610 #endif
1269 chrisfen 2229
1270     if (FF_uses_EAM .and. SIM_uses_EAM) then
1271     call clean_EAM()
1272     endif
1273    
1274     rf = 0.0_dp
1275     tau_Temp = 0.0_dp
1276     virial_Temp = 0.0_dp
1277     end subroutine zero_work_arrays
1278    
1279     function skipThisPair(atom1, atom2) result(skip_it)
1280     integer, intent(in) :: atom1
1281     integer, intent(in), optional :: atom2
1282     logical :: skip_it
1283     integer :: unique_id_1, unique_id_2
1284     integer :: me_i,me_j
1285     integer :: i
1286    
1287     skip_it = .false.
1288    
1289     !! there are a number of reasons to skip a pair or a particle
1290     !! mostly we do this to exclude atoms who are involved in short
1291     !! range interactions (bonds, bends, torsions), but we also need
1292     !! to exclude some overcounted interactions that result from
1293     !! the parallel decomposition
1294    
1295 gezelter 1610 #ifdef IS_MPI
1296 chrisfen 2229 !! in MPI, we have to look up the unique IDs for each atom
1297     unique_id_1 = AtomRowToGlobal(atom1)
1298 gezelter 1610 #else
1299 chrisfen 2229 !! in the normal loop, the atom numbers are unique
1300     unique_id_1 = atom1
1301 gezelter 1610 #endif
1302 chrisfen 2229
1303     !! We were called with only one atom, so just check the global exclude
1304     !! list for this atom
1305     if (.not. present(atom2)) then
1306     do i = 1, nExcludes_global
1307     if (excludesGlobal(i) == unique_id_1) then
1308     skip_it = .true.
1309     return
1310     end if
1311     end do
1312     return
1313     end if
1314    
1315 gezelter 1610 #ifdef IS_MPI
1316 chrisfen 2229 unique_id_2 = AtomColToGlobal(atom2)
1317 gezelter 1610 #else
1318 chrisfen 2229 unique_id_2 = atom2
1319 gezelter 1610 #endif
1320 chrisfen 2229
1321 gezelter 1610 #ifdef IS_MPI
1322 chrisfen 2229 !! this situation should only arise in MPI simulations
1323     if (unique_id_1 == unique_id_2) then
1324     skip_it = .true.
1325     return
1326     end if
1327    
1328     !! this prevents us from doing the pair on multiple processors
1329     if (unique_id_1 < unique_id_2) then
1330     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1331     skip_it = .true.
1332     return
1333     endif
1334     else
1335     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1336     skip_it = .true.
1337     return
1338     endif
1339     endif
1340 gezelter 1610 #endif
1341 chrisfen 2229
1342     !! the rest of these situations can happen in all simulations:
1343     do i = 1, nExcludes_global
1344     if ((excludesGlobal(i) == unique_id_1) .or. &
1345     (excludesGlobal(i) == unique_id_2)) then
1346     skip_it = .true.
1347     return
1348     endif
1349     enddo
1350    
1351     do i = 1, nSkipsForAtom(atom1)
1352     if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1353     skip_it = .true.
1354     return
1355     endif
1356     end do
1357    
1358     return
1359     end function skipThisPair
1360    
1361     function FF_UsesDirectionalAtoms() result(doesit)
1362     logical :: doesit
1363 gezelter 2270 doesit = FF_uses_DirectionalAtoms
1364 chrisfen 2229 end function FF_UsesDirectionalAtoms
1365    
1366     function FF_RequiresPrepairCalc() result(doesit)
1367     logical :: doesit
1368     doesit = FF_uses_EAM
1369     end function FF_RequiresPrepairCalc
1370    
1371     function FF_RequiresPostpairCalc() result(doesit)
1372     logical :: doesit
1373     doesit = FF_uses_RF
1374 chrisfen 2295 if (corrMethod == 3) doesit = .true.
1375 chrisfen 2229 end function FF_RequiresPostpairCalc
1376    
1377 gezelter 1610 #ifdef PROFILE
1378 chrisfen 2229 function getforcetime() result(totalforcetime)
1379     real(kind=dp) :: totalforcetime
1380     totalforcetime = forcetime
1381     end function getforcetime
1382 gezelter 1610 #endif
1383    
1384 chrisfen 2229 !! This cleans componets of force arrays belonging only to fortran
1385    
1386     subroutine add_stress_tensor(dpair, fpair)
1387    
1388     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1389    
1390     ! because the d vector is the rj - ri vector, and
1391     ! because fx, fy, fz are the force on atom i, we need a
1392     ! negative sign here:
1393    
1394     tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1395     tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1396     tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1397     tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1398     tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1399     tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1400     tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1401     tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1402     tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1403    
1404     virial_Temp = virial_Temp + &
1405     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1406    
1407     end subroutine add_stress_tensor
1408    
1409 gezelter 1610 end module doForces