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