ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2268
Committed: Fri Jul 29 19:38:27 2005 UTC (18 years, 11 months ago) by gezelter
File size: 42654 byte(s)
Log Message:
fixes in progress

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