ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2298
Committed: Thu Sep 15 02:48:43 2005 UTC (18 years, 9 months ago) by chuckv
File size: 42759 byte(s)
Log Message:
Fixed bug where gtypeMaxCutoff was not initialized after creation. When maxval(gtypeMaxCutoff) was called, the largest random garbage value was returned from the array.

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