ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2342
Committed: Tue Oct 4 19:32:58 2005 UTC (18 years, 9 months ago) by chrisfen
File size: 42763 byte(s)
Log Message:
just some random changes when testing

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