ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2592
Committed: Thu Feb 16 21:40:20 2006 UTC (18 years, 5 months ago) by gezelter
File size: 49977 byte(s)
Log Message:
fixed a problem with cutoff radii being set larger than 1/2 the
length of the shortest box dimension.  added a few fortran utility
routines

File Contents

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