ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2394
Committed: Sun Oct 23 21:08:08 2005 UTC (18 years, 8 months ago) by chrisfen
File size: 45465 byte(s)
Log Message:
streamlined reaction field for dipoles (now a good bit faster) and added reaction field for charges - still need to do charge-dipole RF

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