ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1386
Committed: Fri Oct 23 18:41:09 2009 UTC (15 years, 6 months ago) by gezelter
Original Path: trunk/src/UseTheForce/doForces.F90
File size: 63517 byte(s)
Log Message:
removing MPI responsibilities from the lowest level force routines.  This is
in preparation for migrating these routines (LJ, electrostatic, eam, suttonchen,
gay-berne, sticky) to C++

File Contents

# User Rev Content
1 gezelter 246 !!
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 gezelter 1313 !! A. Meineke, Charles F. Vardeman I, Teng Lin, Christopher
14 gezelter 246 !! 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 117 !! doForces.F90
43     !! module doForces
44     !! Calculates Long Range forces.
45    
46     !! @author Charles F. Vardeman II
47     !! @author Matthew Meineke
48 gezelter 1386 !! @version $Id: doForces.F90,v 1.103 2009-10-23 18:41:08 gezelter Exp $, $Date: 2009-10-23 18:41:08 $, $Name: not supported by cvs2svn $, $Revision: 1.103 $
49 gezelter 117
50 gezelter 246
51 gezelter 117 module doForces
52     use force_globals
53 gezelter 1286 use fForceOptions
54 gezelter 117 use simulation
55     use definitions
56     use atype_module
57     use switcheroo
58     use neighborLists
59     use lj
60 gezelter 246 use sticky
61 gezelter 401 use electrostatic_module
62 gezelter 676 use gayberne
63 chrisfen 143 use shapes
64 gezelter 117 use vector_class
65     use eam
66 chuckv 1162 use MetalNonMetal
67 chuckv 733 use suttonchen
68 gezelter 117 use status
69     #ifdef IS_MPI
70     use mpiSimulation
71     #endif
72    
73     implicit none
74     PRIVATE
75    
76     #define __FORTRAN90
77 gezelter 574 #include "UseTheForce/fCutoffPolicy.h"
78 gezelter 560 #include "UseTheForce/DarkSide/fInteractionMap.h"
79 chrisfen 611 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
80 gezelter 117
81     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 571 logical, save :: haveInteractionHash = .false.
88     logical, save :: haveGtypeCutoffMap = .false.
89 chrisfen 618 logical, save :: haveDefaultCutoffs = .false.
90 gezelter 762 logical, save :: haveSkinThickness = .false.
91     logical, save :: haveElectrostaticSummationMethod = .false.
92     logical, save :: haveCutoffPolicy = .false.
93     logical, save :: VisitCutoffsAfterComputing = .false.
94 chrisfen 998 logical, save :: do_box_dipole = .false.
95 chrisfen 532
96 gezelter 141 logical, save :: FF_uses_DirectionalAtoms
97 gezelter 401 logical, save :: FF_uses_Dipoles
98 gezelter 141 logical, save :: FF_uses_GayBerne
99     logical, save :: FF_uses_EAM
100 chuckv 733 logical, save :: FF_uses_SC
101 chuckv 1162 logical, save :: FF_uses_MNM
102 chuckv 733
103 gezelter 141
104     logical, save :: SIM_uses_DirectionalAtoms
105     logical, save :: SIM_uses_EAM
106 chuckv 733 logical, save :: SIM_uses_SC
107 chuckv 1162 logical, save :: SIM_uses_MNM
108 gezelter 117 logical, save :: SIM_requires_postpair_calc
109     logical, save :: SIM_requires_prepair_calc
110     logical, save :: SIM_uses_PBC
111 gezelter 1126 logical, save :: SIM_uses_AtomicVirial
112 gezelter 117
113 chrisfen 607 integer, save :: electrostaticSummationMethod
114 gezelter 762 integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
115 chrisfen 580
116 gezelter 762 real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut
117     real(kind=dp), save :: skinThickness
118 chrisfen 1129 logical, save :: defaultDoShiftPot
119     logical, save :: defaultDoShiftFrc
120 gezelter 762
121 gezelter 117 public :: init_FF
122 gezelter 762 public :: setCutoffs
123     public :: cWasLame
124     public :: setElectrostaticMethod
125 chrisfen 998 public :: setBoxDipole
126     public :: getBoxDipole
127 gezelter 762 public :: setCutoffPolicy
128     public :: setSkinThickness
129 gezelter 117 public :: do_force_loop
130    
131     #ifdef PROFILE
132     public :: getforcetime
133     real, save :: forceTime = 0
134     real :: forceTimeInitial, forceTimeFinal
135     integer :: nLoops
136     #endif
137 chuckv 561
138 gezelter 571 !! Variables for cutoff mapping and interaction mapping
139     ! Bit hash to determine pair-pair interactions.
140     integer, dimension(:,:), allocatable :: InteractionHash
141     real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff
142 chuckv 651 real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow
143     real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol
144    
145     integer, dimension(:), allocatable, target :: groupToGtypeRow
146     integer, dimension(:), pointer :: groupToGtypeCol => null()
147    
148     real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow
149     real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol
150 gezelter 571 type ::gtypeCutoffs
151     real(kind=dp) :: rcut
152     real(kind=dp) :: rcutsq
153     real(kind=dp) :: rlistsq
154     end type gtypeCutoffs
155     type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap
156 gezelter 574
157 chrisfen 998 real(kind=dp), dimension(3) :: boxDipole
158 gezelter 939
159 gezelter 117 contains
160    
161 gezelter 762 subroutine createInteractionHash()
162 chuckv 561 integer :: nAtypes
163     integer :: i
164     integer :: j
165 gezelter 571 integer :: iHash
166 tim 568 !! Test Types
167 chuckv 561 logical :: i_is_LJ
168     logical :: i_is_Elect
169     logical :: i_is_Sticky
170     logical :: i_is_StickyP
171     logical :: i_is_GB
172     logical :: i_is_EAM
173     logical :: i_is_Shape
174 chuckv 733 logical :: i_is_SC
175 chuckv 561 logical :: j_is_LJ
176     logical :: j_is_Elect
177     logical :: j_is_Sticky
178     logical :: j_is_StickyP
179     logical :: j_is_GB
180     logical :: j_is_EAM
181     logical :: j_is_Shape
182 chuckv 733 logical :: j_is_SC
183 gezelter 576 real(kind=dp) :: myRcut
184    
185 chuckv 561 if (.not. associated(atypes)) then
186 gezelter 762 call handleError("doForces", "atypes was not present before call of createInteractionHash!")
187 chuckv 561 return
188     endif
189    
190     nAtypes = getSize(atypes)
191    
192     if (nAtypes == 0) then
193 gezelter 762 call handleError("doForces", "nAtypes was zero during call of createInteractionHash!")
194 chuckv 561 return
195     end if
196 chrisfen 532
197 chuckv 570 if (.not. allocated(InteractionHash)) then
198     allocate(InteractionHash(nAtypes,nAtypes))
199 chuckv 655 else
200     deallocate(InteractionHash)
201     allocate(InteractionHash(nAtypes,nAtypes))
202 chuckv 561 endif
203 gezelter 571
204     if (.not. allocated(atypeMaxCutoff)) then
205     allocate(atypeMaxCutoff(nAtypes))
206 chuckv 655 else
207     deallocate(atypeMaxCutoff)
208     allocate(atypeMaxCutoff(nAtypes))
209 gezelter 571 endif
210 chuckv 561
211     do i = 1, nAtypes
212     call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
213     call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
214     call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
215     call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
216     call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
217     call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
218     call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
219 chuckv 733 call getElementProperty(atypes, i, "is_SC", i_is_SC)
220 gezelter 117
221 chuckv 561 do j = i, nAtypes
222 chrisfen 532
223 chuckv 561 iHash = 0
224     myRcut = 0.0_dp
225 gezelter 117
226 chuckv 561 call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ)
227     call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect)
228     call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky)
229     call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP)
230     call getElementProperty(atypes, j, "is_GayBerne", j_is_GB)
231     call getElementProperty(atypes, j, "is_EAM", j_is_EAM)
232     call getElementProperty(atypes, j, "is_Shape", j_is_Shape)
233 chuckv 733 call getElementProperty(atypes, j, "is_SC", j_is_SC)
234 gezelter 117
235 chuckv 561 if (i_is_LJ .and. j_is_LJ) then
236 gezelter 562 iHash = ior(iHash, LJ_PAIR)
237     endif
238    
239     if (i_is_Elect .and. j_is_Elect) then
240     iHash = ior(iHash, ELECTROSTATIC_PAIR)
241     endif
242    
243     if (i_is_Sticky .and. j_is_Sticky) then
244     iHash = ior(iHash, STICKY_PAIR)
245     endif
246 chuckv 561
247 gezelter 562 if (i_is_StickyP .and. j_is_StickyP) then
248     iHash = ior(iHash, STICKYPOWER_PAIR)
249     endif
250 chuckv 561
251 gezelter 562 if (i_is_EAM .and. j_is_EAM) then
252     iHash = ior(iHash, EAM_PAIR)
253 chuckv 561 endif
254    
255 chuckv 733 if (i_is_SC .and. j_is_SC) then
256     iHash = ior(iHash, SC_PAIR)
257     endif
258    
259 chuckv 561 if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR)
260     if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ)
261     if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ)
262 chuckv 1162
263     if ((i_is_EAM.or.i_is_SC).and.(.not.(j_is_EAM.or.j_is_SC))) iHash = ior(iHash, MNM_PAIR)
264     if ((j_is_EAM.or.j_is_SC).and.(.not.(i_is_EAM.or.i_is_SC))) iHash = ior(iHash, MNM_PAIR)
265 chuckv 561
266     if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR)
267     if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ)
268     if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ)
269    
270    
271 chuckv 570 InteractionHash(i,j) = iHash
272     InteractionHash(j,i) = iHash
273 chuckv 561
274     end do
275    
276     end do
277 tim 568
278 gezelter 571 haveInteractionHash = .true.
279     end subroutine createInteractionHash
280 chuckv 561
281 gezelter 762 subroutine createGtypeCutoffMap()
282 gezelter 569
283 gezelter 574 logical :: i_is_LJ
284     logical :: i_is_Elect
285     logical :: i_is_Sticky
286     logical :: i_is_StickyP
287     logical :: i_is_GB
288     logical :: i_is_EAM
289     logical :: i_is_Shape
290 chuckv 831 logical :: i_is_SC
291 gezelter 587 logical :: GtypeFound
292 chuckv 561
293 gezelter 576 integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend
294 chuckv 652 integer :: n_in_i, me_i, ia, g, atom1, ja, n_in_j,me_j
295 chuckv 589 integer :: nGroupsInRow
296 chuckv 651 integer :: nGroupsInCol
297     integer :: nGroupTypesRow,nGroupTypesCol
298 gezelter 762 real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol
299 gezelter 576 real(kind=dp) :: biggestAtypeCutoff
300 gezelter 571
301     if (.not. haveInteractionHash) then
302 gezelter 762 call createInteractionHash()
303 chuckv 567 endif
304 chuckv 589 #ifdef IS_MPI
305     nGroupsInRow = getNgroupsInRow(plan_group_row)
306 chuckv 651 nGroupsInCol = getNgroupsInCol(plan_group_col)
307 chuckv 589 #endif
308 chuckv 563 nAtypes = getSize(atypes)
309 chuckv 599 ! Set all of the initial cutoffs to zero.
310     atypeMaxCutoff = 0.0_dp
311 gezelter 1313 biggestAtypeCutoff = 0.0_dp
312 gezelter 571 do i = 1, nAtypes
313 gezelter 582 if (SimHasAtype(i)) then
314 gezelter 575 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
315     call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
316     call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
317     call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
318     call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
319     call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
320     call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
321 chuckv 831 call getElementProperty(atypes, i, "is_SC", i_is_SC)
322 chuckv 599
323 chrisfen 618 if (haveDefaultCutoffs) then
324     atypeMaxCutoff(i) = defaultRcut
325     else
326     if (i_is_LJ) then
327     thisRcut = getSigma(i) * 2.5_dp
328     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
329     endif
330     if (i_is_Elect) then
331     thisRcut = defaultRcut
332     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
333     endif
334     if (i_is_Sticky) then
335     thisRcut = getStickyCut(i)
336     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
337     endif
338     if (i_is_StickyP) then
339     thisRcut = getStickyPowerCut(i)
340     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
341     endif
342     if (i_is_GB) then
343     thisRcut = getGayBerneCut(i)
344     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
345     endif
346     if (i_is_EAM) then
347     thisRcut = getEAMCut(i)
348     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
349     endif
350     if (i_is_Shape) then
351     thisRcut = getShapeCut(i)
352     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
353     endif
354 chuckv 831 if (i_is_SC) then
355     thisRcut = getSCCut(i)
356     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
357     endif
358 gezelter 575 endif
359 gezelter 762
360 gezelter 575 if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
361     biggestAtypeCutoff = atypeMaxCutoff(i)
362     endif
363 chrisfen 618
364 gezelter 574 endif
365 gezelter 575 enddo
366 gezelter 581
367 gezelter 575 istart = 1
368 chuckv 651 jstart = 1
369 gezelter 575 #ifdef IS_MPI
370     iend = nGroupsInRow
371 chuckv 651 jend = nGroupsInCol
372 gezelter 575 #else
373     iend = nGroups
374 chuckv 651 jend = nGroups
375 gezelter 575 #endif
376 gezelter 582
377 gezelter 581 !! allocate the groupToGtype and gtypeMaxCutoff here.
378 chuckv 651 if(.not.allocated(groupToGtypeRow)) then
379     ! allocate(groupToGtype(iend))
380     allocate(groupToGtypeRow(iend))
381     else
382     deallocate(groupToGtypeRow)
383     allocate(groupToGtypeRow(iend))
384 chuckv 583 endif
385 chuckv 651 if(.not.allocated(groupMaxCutoffRow)) then
386     allocate(groupMaxCutoffRow(iend))
387     else
388     deallocate(groupMaxCutoffRow)
389     allocate(groupMaxCutoffRow(iend))
390     end if
391    
392     if(.not.allocated(gtypeMaxCutoffRow)) then
393     allocate(gtypeMaxCutoffRow(iend))
394     else
395     deallocate(gtypeMaxCutoffRow)
396     allocate(gtypeMaxCutoffRow(iend))
397     endif
398    
399    
400     #ifdef IS_MPI
401     ! We only allocate new storage if we are in MPI because Ncol /= Nrow
402 chuckv 652 if(.not.associated(groupToGtypeCol)) then
403 chuckv 651 allocate(groupToGtypeCol(jend))
404     else
405     deallocate(groupToGtypeCol)
406     allocate(groupToGtypeCol(jend))
407     end if
408    
409 tim 833 if(.not.associated(groupMaxCutoffCol)) then
410     allocate(groupMaxCutoffCol(jend))
411 chuckv 651 else
412 tim 833 deallocate(groupMaxCutoffCol)
413     allocate(groupMaxCutoffCol(jend))
414 chuckv 651 end if
415 chuckv 652 if(.not.associated(gtypeMaxCutoffCol)) then
416 chuckv 651 allocate(gtypeMaxCutoffCol(jend))
417     else
418     deallocate(gtypeMaxCutoffCol)
419     allocate(gtypeMaxCutoffCol(jend))
420     end if
421    
422     groupMaxCutoffCol = 0.0_dp
423     gtypeMaxCutoffCol = 0.0_dp
424    
425     #endif
426     groupMaxCutoffRow = 0.0_dp
427     gtypeMaxCutoffRow = 0.0_dp
428    
429    
430 gezelter 582 !! first we do a single loop over the cutoff groups to find the
431     !! largest cutoff for any atypes present in this group. We also
432     !! create gtypes at this point.
433    
434 gezelter 960 tol = 1.0e-6_dp
435 chuckv 651 nGroupTypesRow = 0
436 tim 833 nGroupTypesCol = 0
437 gezelter 581 do i = istart, iend
438 gezelter 575 n_in_i = groupStartRow(i+1) - groupStartRow(i)
439 chuckv 651 groupMaxCutoffRow(i) = 0.0_dp
440 gezelter 581 do ia = groupStartRow(i), groupStartRow(i+1)-1
441     atom1 = groupListRow(ia)
442 gezelter 575 #ifdef IS_MPI
443 gezelter 581 me_i = atid_row(atom1)
444 gezelter 575 #else
445 gezelter 581 me_i = atid(atom1)
446     #endif
447 chuckv 651 if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
448     groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
449 gezelter 587 endif
450 gezelter 581 enddo
451 chuckv 651 if (nGroupTypesRow.eq.0) then
452     nGroupTypesRow = nGroupTypesRow + 1
453     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
454     groupToGtypeRow(i) = nGroupTypesRow
455 gezelter 581 else
456 gezelter 587 GtypeFound = .false.
457 chuckv 651 do g = 1, nGroupTypesRow
458     if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
459     groupToGtypeRow(i) = g
460 gezelter 587 GtypeFound = .true.
461 gezelter 581 endif
462     enddo
463 gezelter 587 if (.not.GtypeFound) then
464 chuckv 651 nGroupTypesRow = nGroupTypesRow + 1
465     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
466     groupToGtypeRow(i) = nGroupTypesRow
467 gezelter 587 endif
468 gezelter 581 endif
469 gezelter 587 enddo
470    
471 chuckv 651 #ifdef IS_MPI
472     do j = jstart, jend
473     n_in_j = groupStartCol(j+1) - groupStartCol(j)
474     groupMaxCutoffCol(j) = 0.0_dp
475     do ja = groupStartCol(j), groupStartCol(j+1)-1
476     atom1 = groupListCol(ja)
477    
478     me_j = atid_col(atom1)
479    
480     if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then
481     groupMaxCutoffCol(j)=atypeMaxCutoff(me_j)
482     endif
483     enddo
484    
485     if (nGroupTypesCol.eq.0) then
486     nGroupTypesCol = nGroupTypesCol + 1
487     gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
488     groupToGtypeCol(j) = nGroupTypesCol
489     else
490     GtypeFound = .false.
491     do g = 1, nGroupTypesCol
492     if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then
493     groupToGtypeCol(j) = g
494     GtypeFound = .true.
495     endif
496     enddo
497     if (.not.GtypeFound) then
498     nGroupTypesCol = nGroupTypesCol + 1
499     gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
500     groupToGtypeCol(j) = nGroupTypesCol
501     endif
502     endif
503     enddo
504    
505     #else
506     ! Set pointers to information we just found
507     nGroupTypesCol = nGroupTypesRow
508     groupToGtypeCol => groupToGtypeRow
509     gtypeMaxCutoffCol => gtypeMaxCutoffRow
510     groupMaxCutoffCol => groupMaxCutoffRow
511     #endif
512    
513 gezelter 581 !! allocate the gtypeCutoffMap here.
514 chuckv 651 allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
515 gezelter 581 !! then we do a double loop over all the group TYPES to find the cutoff
516     !! map between groups of two types
517 chuckv 651 tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
518    
519 gezelter 762 do i = 1, nGroupTypesRow
520 chuckv 651 do j = 1, nGroupTypesCol
521 gezelter 576
522 gezelter 581 select case(cutoffPolicy)
523 gezelter 582 case(TRADITIONAL_CUTOFF_POLICY)
524 chuckv 651 thisRcut = tradRcut
525 gezelter 582 case(MIX_CUTOFF_POLICY)
526 chuckv 651 thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
527 gezelter 582 case(MAX_CUTOFF_POLICY)
528 chuckv 651 thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
529 gezelter 582 case default
530     call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
531     return
532     end select
533     gtypeCutoffMap(i,j)%rcut = thisRcut
534 gezelter 762
535     if (thisRcut.gt.largestRcut) largestRcut = thisRcut
536    
537 gezelter 582 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
538 gezelter 585
539 gezelter 762 if (.not.haveSkinThickness) then
540     skinThickness = 1.0_dp
541     endif
542    
543     gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2
544    
545 chrisfen 618 ! sanity check
546    
547     if (haveDefaultCutoffs) then
548     if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then
549     call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff")
550     endif
551     endif
552 gezelter 581 enddo
553     enddo
554 gezelter 762
555 chuckv 651 if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow)
556     if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow)
557     if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff)
558     #ifdef IS_MPI
559     if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol)
560     if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol)
561     #endif
562     groupMaxCutoffCol => null()
563     gtypeMaxCutoffCol => null()
564    
565 gezelter 581 haveGtypeCutoffMap = .true.
566 chrisfen 596 end subroutine createGtypeCutoffMap
567 chrisfen 578
568 chrisfen 1129 subroutine setCutoffs(defRcut, defRsw, defSP, defSF)
569 chrisfen 596
570 gezelter 762 real(kind=dp),intent(in) :: defRcut, defRsw
571 gezelter 1386 integer, intent(in) :: defSP, defSF
572 gezelter 762 character(len = statusMsgSize) :: errMsg
573     integer :: localError
574    
575 chrisfen 596 defaultRcut = defRcut
576     defaultRsw = defRsw
577 gezelter 1386
578     if (defSP .ne. 0) then
579     defaultDoShiftPot = .true.
580     else
581     defaultDoShiftPot = .false.
582     endif
583     if (defSF .ne. 0) then
584     defaultDoShiftFrc = .true.
585     else
586     defaultDoShiftFrc = .false.
587     endif
588 chrisfen 1129
589 gezelter 762 if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then
590 chrisfen 1129 if (defaultDoShiftFrc) then
591     write(errMsg, *) &
592     'cutoffRadius and switchingRadius are set to the', newline &
593     // tab, 'same value. OOPSE will use shifted force', newline &
594     // tab, 'potentials instead of switching functions.'
595    
596     call handleInfo("setCutoffs", errMsg)
597     else
598     write(errMsg, *) &
599     'cutoffRadius and switchingRadius are set to the', newline &
600     // tab, 'same value. OOPSE will use shifted', newline &
601     // tab, 'potentials instead of switching functions.'
602    
603     call handleInfo("setCutoffs", errMsg)
604    
605     defaultDoShiftPot = .true.
606     endif
607    
608 gezelter 762 endif
609 gezelter 939
610 gezelter 762 localError = 0
611 chrisfen 1129 call setLJDefaultCutoff( defaultRcut, defaultDoShiftPot, &
612     defaultDoShiftFrc )
613 gezelter 813 call setElectrostaticCutoffRadius( defaultRcut, defaultRsw )
614 gezelter 938 call setCutoffEAM( defaultRcut )
615     call setCutoffSC( defaultRcut )
616 chuckv 1162 call setMnMDefaultCutoff( defaultRcut, defaultDoShiftPot, &
617     defaultDoShiftFrc )
618 gezelter 939 call set_switch(defaultRsw, defaultRcut)
619 gezelter 889 call setHmatDangerousRcutValue(defaultRcut)
620 gezelter 939
621 chrisfen 618 haveDefaultCutoffs = .true.
622 gezelter 813 haveGtypeCutoffMap = .false.
623 gezelter 939
624 gezelter 762 end subroutine setCutoffs
625 chrisfen 596
626 gezelter 762 subroutine cWasLame()
627    
628     VisitCutoffsAfterComputing = .true.
629     return
630    
631     end subroutine cWasLame
632    
633 chrisfen 596 subroutine setCutoffPolicy(cutPolicy)
634 gezelter 762
635 chrisfen 596 integer, intent(in) :: cutPolicy
636 gezelter 762
637 chrisfen 596 cutoffPolicy = cutPolicy
638 gezelter 762 haveCutoffPolicy = .true.
639 gezelter 813 haveGtypeCutoffMap = .false.
640 gezelter 762
641 gezelter 576 end subroutine setCutoffPolicy
642 gezelter 1126
643 chrisfen 998 subroutine setBoxDipole()
644    
645     do_box_dipole = .true.
646    
647     end subroutine setBoxDipole
648    
649     subroutine getBoxDipole( box_dipole )
650    
651     real(kind=dp), intent(inout), dimension(3) :: box_dipole
652    
653     box_dipole = boxDipole
654    
655     end subroutine getBoxDipole
656    
657 gezelter 762 subroutine setElectrostaticMethod( thisESM )
658    
659     integer, intent(in) :: thisESM
660    
661     electrostaticSummationMethod = thisESM
662     haveElectrostaticSummationMethod = .true.
663 gezelter 574
664 gezelter 762 end subroutine setElectrostaticMethod
665    
666     subroutine setSkinThickness( thisSkin )
667 gezelter 574
668 gezelter 762 real(kind=dp), intent(in) :: thisSkin
669    
670     skinThickness = thisSkin
671 gezelter 813 haveSkinThickness = .true.
672     haveGtypeCutoffMap = .false.
673 gezelter 762
674     end subroutine setSkinThickness
675    
676     subroutine setSimVariables()
677     SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
678     SIM_uses_EAM = SimUsesEAM()
679     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
680     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
681     SIM_uses_PBC = SimUsesPBC()
682 chuckv 841 SIM_uses_SC = SimUsesSC()
683 gezelter 1126 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
684 chrisfen 998
685 gezelter 762 haveSIMvariables = .true.
686    
687     return
688     end subroutine setSimVariables
689 gezelter 117
690     subroutine doReadyCheck(error)
691     integer, intent(out) :: error
692     integer :: myStatus
693    
694     error = 0
695 chrisfen 532
696 gezelter 571 if (.not. haveInteractionHash) then
697 gezelter 762 call createInteractionHash()
698 gezelter 117 endif
699    
700 gezelter 571 if (.not. haveGtypeCutoffMap) then
701 gezelter 762 call createGtypeCutoffMap()
702 gezelter 571 endif
703    
704 gezelter 762 if (VisitCutoffsAfterComputing) then
705 gezelter 939 call set_switch(largestRcut, largestRcut)
706 gezelter 889 call setHmatDangerousRcutValue(largestRcut)
707 gezelter 938 call setCutoffEAM(largestRcut)
708     call setCutoffSC(largestRcut)
709     VisitCutoffsAfterComputing = .false.
710 gezelter 762 endif
711    
712 gezelter 117 if (.not. haveSIMvariables) then
713     call setSimVariables()
714     endif
715    
716     if (.not. haveNeighborList) then
717     write(default_error, *) 'neighbor list has not been initialized in doForces!'
718     error = -1
719     return
720     end if
721 gezelter 939
722 gezelter 117 if (.not. haveSaneForceField) then
723     write(default_error, *) 'Force Field is not sane in doForces!'
724     error = -1
725     return
726     end if
727 gezelter 939
728 gezelter 117 #ifdef IS_MPI
729     if (.not. isMPISimSet()) then
730     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
731     error = -1
732     return
733     endif
734     #endif
735     return
736     end subroutine doReadyCheck
737    
738 chrisfen 532
739 gezelter 762 subroutine init_FF(thisStat)
740 gezelter 117
741     integer, intent(out) :: thisStat
742     integer :: my_status, nMatches
743     integer, pointer :: MatchList(:) => null()
744    
745     !! assume things are copacetic, unless they aren't
746     thisStat = 0
747    
748     !! init_FF is called *after* all of the atom types have been
749     !! defined in atype_module using the new_atype subroutine.
750     !!
751     !! this will scan through the known atypes and figure out what
752     !! interactions are used by the force field.
753 chrisfen 532
754 gezelter 141 FF_uses_DirectionalAtoms = .false.
755     FF_uses_Dipoles = .false.
756     FF_uses_GayBerne = .false.
757 gezelter 117 FF_uses_EAM = .false.
758 chuckv 834 FF_uses_SC = .false.
759 chrisfen 532
760 gezelter 141 call getMatchingElementList(atypes, "is_Directional", .true., &
761     nMatches, MatchList)
762     if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
763    
764     call getMatchingElementList(atypes, "is_Dipole", .true., &
765     nMatches, MatchList)
766 gezelter 571 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
767 chrisfen 523
768 gezelter 141 call getMatchingElementList(atypes, "is_GayBerne", .true., &
769     nMatches, MatchList)
770 gezelter 571 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
771 chrisfen 532
772 gezelter 117 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
773     if (nMatches .gt. 0) FF_uses_EAM = .true.
774 chrisfen 532
775 chuckv 834 call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList)
776     if (nMatches .gt. 0) FF_uses_SC = .true.
777 gezelter 141
778 chuckv 834
779 gezelter 117 haveSaneForceField = .true.
780 chrisfen 532
781 gezelter 117 if (FF_uses_EAM) then
782 chrisfen 532 call init_EAM_FF(my_status)
783 gezelter 117 if (my_status /= 0) then
784     write(default_error, *) "init_EAM_FF returned a bad status"
785     thisStat = -1
786     haveSaneForceField = .false.
787     return
788     end if
789     endif
790    
791     if (.not. haveNeighborList) then
792     !! Create neighbor lists
793     call expandNeighborList(nLocal, my_status)
794     if (my_Status /= 0) then
795     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
796     thisStat = -1
797     return
798     endif
799     haveNeighborList = .true.
800 chrisfen 532 endif
801    
802 gezelter 117 end subroutine init_FF
803    
804 chrisfen 532
805 gezelter 117 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
806     !------------------------------------------------------------->
807 gezelter 1285 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
808 gezelter 117 do_pot_c, do_stress_c, error)
809     !! Position array provided by C, dimensioned by getNlocal
810     real ( kind = dp ), dimension(3, nLocal) :: q
811     !! molecular center-of-mass position array
812     real ( kind = dp ), dimension(3, nGroups) :: q_group
813     !! Rotation Matrix for each long range particle in simulation.
814     real( kind = dp), dimension(9, nLocal) :: A
815     !! Unit vectors for dipoles (lab frame)
816 gezelter 246 real( kind = dp ), dimension(9,nLocal) :: eFrame
817 gezelter 117 !! Force array provided by C, dimensioned by getNlocal
818     real ( kind = dp ), dimension(3,nLocal) :: f
819     !! Torsion array provided by C, dimensioned by getNlocal
820     real( kind = dp ), dimension(3,nLocal) :: t
821    
822     !! Stress Tensor
823     real( kind = dp), dimension(9) :: tau
824 gezelter 662 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
825 chuckv 1245 real( kind = dp ), dimension(nLocal) :: particle_pot
826 gezelter 117 logical ( kind = 2) :: do_pot_c, do_stress_c
827     logical :: do_pot
828     logical :: do_stress
829     logical :: in_switching_region
830     #ifdef IS_MPI
831 gezelter 662 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
832 gezelter 117 integer :: nAtomsInRow
833     integer :: nAtomsInCol
834     integer :: nprocs
835     integer :: nGroupsInRow
836     integer :: nGroupsInCol
837     #endif
838     integer :: natoms
839     logical :: update_nlist
840     integer :: i, j, jstart, jend, jnab
841     integer :: istart, iend
842     integer :: ia, jb, atom1, atom2
843     integer :: nlist
844 gezelter 1126 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
845 gezelter 117 real( kind = DP ) :: sw, dswdr, swderiv, mf
846 chrisfen 699 real( kind = DP ) :: rVal
847 gezelter 1126 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
848     real(kind=dp) :: rfpot, mu_i
849 gezelter 762 real(kind=dp):: rCut
850 gezelter 1345 integer :: me_i, me_j, n_in_i, n_in_j, iG, j1
851 gezelter 117 logical :: is_dp_i
852     integer :: neighborListSize
853     integer :: listerror, error
854     integer :: localError
855     integer :: propPack_i, propPack_j
856     integer :: loopStart, loopEnd, loop
857 gezelter 1352 integer :: iHash, jHash
858 gezelter 1286 integer :: i1, topoDist
859 chrisfen 532
860 chrisfen 998 !! the variables for the box dipole moment
861     #ifdef IS_MPI
862     integer :: pChgCount_local
863     integer :: nChgCount_local
864     real(kind=dp) :: pChg_local
865     real(kind=dp) :: nChg_local
866     real(kind=dp), dimension(3) :: pChgPos_local
867     real(kind=dp), dimension(3) :: nChgPos_local
868     real(kind=dp), dimension(3) :: dipVec_local
869     #endif
870     integer :: pChgCount
871     integer :: nChgCount
872     real(kind=dp) :: pChg
873     real(kind=dp) :: nChg
874     real(kind=dp) :: chg_value
875     real(kind=dp), dimension(3) :: pChgPos
876     real(kind=dp), dimension(3) :: nChgPos
877     real(kind=dp), dimension(3) :: dipVec
878     real(kind=dp), dimension(3) :: chgVec
879 gezelter 1340 real(kind=dp) :: skch
880 chrisfen 998
881     !! initialize box dipole variables
882     if (do_box_dipole) then
883     #ifdef IS_MPI
884     pChg_local = 0.0_dp
885     nChg_local = 0.0_dp
886     pChgCount_local = 0
887     nChgCount_local = 0
888     do i=1, 3
889     pChgPos_local = 0.0_dp
890     nChgPos_local = 0.0_dp
891     dipVec_local = 0.0_dp
892     enddo
893     #endif
894     pChg = 0.0_dp
895     nChg = 0.0_dp
896     pChgCount = 0
897     nChgCount = 0
898     chg_value = 0.0_dp
899    
900     do i=1, 3
901     pChgPos(i) = 0.0_dp
902     nChgPos(i) = 0.0_dp
903     dipVec(i) = 0.0_dp
904     chgVec(i) = 0.0_dp
905     boxDipole(i) = 0.0_dp
906     enddo
907     endif
908    
909 gezelter 117 !! initialize local variables
910 chrisfen 532
911 gezelter 117 #ifdef IS_MPI
912     pot_local = 0.0_dp
913     nAtomsInRow = getNatomsInRow(plan_atom_row)
914     nAtomsInCol = getNatomsInCol(plan_atom_col)
915     nGroupsInRow = getNgroupsInRow(plan_group_row)
916     nGroupsInCol = getNgroupsInCol(plan_group_col)
917     #else
918     natoms = nlocal
919     #endif
920 chrisfen 532
921 gezelter 117 call doReadyCheck(localError)
922     if ( localError .ne. 0 ) then
923     call handleError("do_force_loop", "Not Initialized")
924     error = -1
925     return
926     end if
927     call zero_work_arrays()
928 chrisfen 532
929 gezelter 117 do_pot = do_pot_c
930     do_stress = do_stress_c
931 chrisfen 532
932 gezelter 117 ! Gather all information needed by all force loops:
933 chrisfen 532
934 gezelter 117 #ifdef IS_MPI
935 chrisfen 532
936 gezelter 117 call gather(q, q_Row, plan_atom_row_3d)
937     call gather(q, q_Col, plan_atom_col_3d)
938    
939     call gather(q_group, q_group_Row, plan_group_row_3d)
940     call gather(q_group, q_group_Col, plan_group_col_3d)
941 chrisfen 532
942 gezelter 141 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
943 gezelter 246 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
944     call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
945 chrisfen 532
946 gezelter 117 call gather(A, A_Row, plan_atom_row_rotation)
947     call gather(A, A_Col, plan_atom_col_rotation)
948     endif
949 chrisfen 532
950 gezelter 117 #endif
951 chrisfen 532
952 gezelter 117 !! Begin force loop timing:
953     #ifdef PROFILE
954     call cpu_time(forceTimeInitial)
955     nloops = nloops + 1
956     #endif
957 chrisfen 532
958 gezelter 117 loopEnd = PAIR_LOOP
959     if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
960     loopStart = PREPAIR_LOOP
961     else
962     loopStart = PAIR_LOOP
963     endif
964    
965     do loop = loopStart, loopEnd
966    
967     ! See if we need to update neighbor lists
968     ! (but only on the first time through):
969     if (loop .eq. loopStart) then
970     #ifdef IS_MPI
971 gezelter 762 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
972 chrisfen 532 update_nlist)
973 gezelter 117 #else
974 gezelter 762 call checkNeighborList(nGroups, q_group, skinThickness, &
975 chrisfen 532 update_nlist)
976 gezelter 117 #endif
977     endif
978 chrisfen 532
979 gezelter 117 if (update_nlist) then
980     !! save current configuration and construct neighbor list
981     #ifdef IS_MPI
982     call saveNeighborList(nGroupsInRow, q_group_row)
983     #else
984     call saveNeighborList(nGroups, q_group)
985     #endif
986     neighborListSize = size(list)
987     nlist = 0
988     endif
989 chrisfen 532
990 gezelter 117 istart = 1
991     #ifdef IS_MPI
992     iend = nGroupsInRow
993     #else
994     iend = nGroups - 1
995     #endif
996     outer: do i = istart, iend
997    
998     if (update_nlist) point(i) = nlist + 1
999 chrisfen 532
1000 gezelter 117 n_in_i = groupStartRow(i+1) - groupStartRow(i)
1001 chrisfen 532
1002 gezelter 117 if (update_nlist) then
1003     #ifdef IS_MPI
1004     jstart = 1
1005     jend = nGroupsInCol
1006     #else
1007     jstart = i+1
1008     jend = nGroups
1009     #endif
1010     else
1011     jstart = point(i)
1012     jend = point(i+1) - 1
1013     ! make sure group i has neighbors
1014     if (jstart .gt. jend) cycle outer
1015     endif
1016 chrisfen 532
1017 gezelter 117 do jnab = jstart, jend
1018     if (update_nlist) then
1019     j = jnab
1020     else
1021     j = list(jnab)
1022     endif
1023    
1024     #ifdef IS_MPI
1025 chuckv 567 me_j = atid_col(j)
1026 gezelter 117 call get_interatomic_vector(q_group_Row(:,i), &
1027     q_group_Col(:,j), d_grp, rgrpsq)
1028     #else
1029 chuckv 567 me_j = atid(j)
1030 gezelter 117 call get_interatomic_vector(q_group(:,i), &
1031     q_group(:,j), d_grp, rgrpsq)
1032 chrisfen 618 #endif
1033 gezelter 117
1034 chuckv 651 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1035 gezelter 117 if (update_nlist) then
1036     nlist = nlist + 1
1037 chrisfen 532
1038 gezelter 117 if (nlist > neighborListSize) then
1039     #ifdef IS_MPI
1040     call expandNeighborList(nGroupsInRow, listerror)
1041     #else
1042     call expandNeighborList(nGroups, listerror)
1043     #endif
1044     if (listerror /= 0) then
1045     error = -1
1046     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1047     return
1048     end if
1049     neighborListSize = size(list)
1050     endif
1051 chrisfen 532
1052 gezelter 117 list(nlist) = j
1053     endif
1054 gezelter 939
1055 chrisfen 708 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1056 chrisfen 532
1057 gezelter 762 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1058 chrisfen 708 if (loop .eq. PAIR_LOOP) then
1059 gezelter 960 vij = 0.0_dp
1060 gezelter 938 fij(1) = 0.0_dp
1061     fij(2) = 0.0_dp
1062     fij(3) = 0.0_dp
1063 chrisfen 708 endif
1064    
1065 gezelter 939 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1066 chrisfen 708
1067     n_in_j = groupStartCol(j+1) - groupStartCol(j)
1068    
1069     do ia = groupStartRow(i), groupStartRow(i+1)-1
1070 chrisfen 703
1071 chrisfen 708 atom1 = groupListRow(ia)
1072    
1073     inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1074    
1075     atom2 = groupListCol(jb)
1076    
1077     if (skipThisPair(atom1, atom2)) cycle inner
1078    
1079     if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1080 gezelter 938 d_atm(1) = d_grp(1)
1081     d_atm(2) = d_grp(2)
1082     d_atm(3) = d_grp(3)
1083 chrisfen 708 ratmsq = rgrpsq
1084     else
1085 gezelter 117 #ifdef IS_MPI
1086 chrisfen 708 call get_interatomic_vector(q_Row(:,atom1), &
1087     q_Col(:,atom2), d_atm, ratmsq)
1088 gezelter 117 #else
1089 chrisfen 708 call get_interatomic_vector(q(:,atom1), &
1090     q(:,atom2), d_atm, ratmsq)
1091 gezelter 117 #endif
1092 gezelter 1286 endif
1093    
1094     topoDist = getTopoDistance(atom1, atom2)
1095    
1096 chrisfen 708 if (loop .eq. PREPAIR_LOOP) then
1097 gezelter 117 #ifdef IS_MPI
1098 chrisfen 708 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1099 gezelter 762 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1100 chrisfen 708 eFrame, A, f, t, pot_local)
1101 gezelter 117 #else
1102 chrisfen 708 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1103 gezelter 762 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1104 chrisfen 708 eFrame, A, f, t, pot)
1105 gezelter 117 #endif
1106 chrisfen 708 else
1107 gezelter 117 #ifdef IS_MPI
1108 chrisfen 708 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1109 gezelter 1309 do_pot, eFrame, A, f, t, pot_local, particle_pot, vpair, &
1110 gezelter 1286 fpair, d_grp, rgrp, rCut, topoDist)
1111 chuckv 1245 ! particle_pot will be accumulated from row & column
1112     ! arrays later
1113 gezelter 117 #else
1114 chrisfen 708 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1115 gezelter 1309 do_pot, eFrame, A, f, t, pot, particle_pot, vpair, &
1116 gezelter 1286 fpair, d_grp, rgrp, rCut, topoDist)
1117 gezelter 117 #endif
1118 chrisfen 708 vij = vij + vpair
1119 gezelter 938 fij(1) = fij(1) + fpair(1)
1120     fij(2) = fij(2) + fpair(2)
1121     fij(3) = fij(3) + fpair(3)
1122 gezelter 1127 if (do_stress) then
1123 gezelter 1126 call add_stress_tensor(d_atm, fpair, tau)
1124     endif
1125 chrisfen 708 endif
1126     enddo inner
1127     enddo
1128 gezelter 117
1129 chrisfen 708 if (loop .eq. PAIR_LOOP) then
1130     if (in_switching_region) then
1131     swderiv = vij*dswdr/rgrp
1132 chrisfen 1131 fg = swderiv*d_grp
1133    
1134     fij(1) = fij(1) + fg(1)
1135     fij(2) = fij(2) + fg(2)
1136     fij(3) = fij(3) + fg(3)
1137 chrisfen 708
1138 chrisfen 1132 if (do_stress .and. (n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1139 chrisfen 1131 call add_stress_tensor(d_atm, fg, tau)
1140     endif
1141    
1142 chrisfen 708 do ia=groupStartRow(i), groupStartRow(i+1)-1
1143     atom1=groupListRow(ia)
1144     mf = mfactRow(atom1)
1145 gezelter 1126 ! fg is the force on atom ia due to cutoff group's
1146     ! presence in switching region
1147     fg = swderiv*d_grp*mf
1148 gezelter 117 #ifdef IS_MPI
1149 gezelter 1126 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1150     f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1151     f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1152 gezelter 117 #else
1153 gezelter 1126 f(1,atom1) = f(1,atom1) + fg(1)
1154     f(2,atom1) = f(2,atom1) + fg(2)
1155     f(3,atom1) = f(3,atom1) + fg(3)
1156 gezelter 117 #endif
1157 gezelter 1127 if (n_in_i .gt. 1) then
1158     if (do_stress.and.SIM_uses_AtomicVirial) then
1159     ! find the distance between the atom and the center of
1160     ! the cutoff group:
1161 gezelter 1126 #ifdef IS_MPI
1162 gezelter 1127 call get_interatomic_vector(q_Row(:,atom1), &
1163     q_group_Row(:,i), dag, rag)
1164 gezelter 1126 #else
1165 gezelter 1127 call get_interatomic_vector(q(:,atom1), &
1166     q_group(:,i), dag, rag)
1167 gezelter 1126 #endif
1168 gezelter 1127 call add_stress_tensor(dag,fg,tau)
1169     endif
1170 gezelter 1126 endif
1171 chrisfen 708 enddo
1172    
1173     do jb=groupStartCol(j), groupStartCol(j+1)-1
1174     atom2=groupListCol(jb)
1175     mf = mfactCol(atom2)
1176 gezelter 1126 ! fg is the force on atom jb due to cutoff group's
1177     ! presence in switching region
1178     fg = -swderiv*d_grp*mf
1179 gezelter 117 #ifdef IS_MPI
1180 gezelter 1126 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1181     f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1182     f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1183 gezelter 117 #else
1184 gezelter 1126 f(1,atom2) = f(1,atom2) + fg(1)
1185     f(2,atom2) = f(2,atom2) + fg(2)
1186     f(3,atom2) = f(3,atom2) + fg(3)
1187 gezelter 117 #endif
1188 gezelter 1127 if (n_in_j .gt. 1) then
1189     if (do_stress.and.SIM_uses_AtomicVirial) then
1190     ! find the distance between the atom and the center of
1191     ! the cutoff group:
1192 gezelter 1126 #ifdef IS_MPI
1193 gezelter 1127 call get_interatomic_vector(q_Col(:,atom2), &
1194     q_group_Col(:,j), dag, rag)
1195 gezelter 1126 #else
1196 gezelter 1127 call get_interatomic_vector(q(:,atom2), &
1197     q_group(:,j), dag, rag)
1198 gezelter 1126 #endif
1199 gezelter 1127 call add_stress_tensor(dag,fg,tau)
1200     endif
1201     endif
1202 chrisfen 708 enddo
1203     endif
1204 gezelter 1174 !if (do_stress.and.(.not.SIM_uses_AtomicVirial)) then
1205     ! call add_stress_tensor(d_grp, fij, tau)
1206     !endif
1207 gezelter 117 endif
1208     endif
1209 chrisfen 708 endif
1210 gezelter 117 enddo
1211 chrisfen 708
1212 gezelter 117 enddo outer
1213 chrisfen 532
1214 gezelter 117 if (update_nlist) then
1215     #ifdef IS_MPI
1216     point(nGroupsInRow + 1) = nlist + 1
1217     #else
1218     point(nGroups) = nlist + 1
1219     #endif
1220     if (loop .eq. PREPAIR_LOOP) then
1221     ! we just did the neighbor list update on the first
1222     ! pass, so we don't need to do it
1223     ! again on the second pass
1224     update_nlist = .false.
1225     endif
1226     endif
1227 chrisfen 532
1228 gezelter 117 if (loop .eq. PREPAIR_LOOP) then
1229 chuckv 1133 #ifdef IS_MPI
1230 gezelter 1285 call do_preforce(nlocal, pot_local, particle_pot)
1231 chuckv 1133 #else
1232 gezelter 1285 call do_preforce(nlocal, pot, particle_pot)
1233 chuckv 1133 #endif
1234 gezelter 117 endif
1235 chrisfen 532
1236 gezelter 117 enddo
1237 chrisfen 532
1238 gezelter 117 !! Do timing
1239     #ifdef PROFILE
1240     call cpu_time(forceTimeFinal)
1241     forceTime = forceTime + forceTimeFinal - forceTimeInitial
1242     #endif
1243 chrisfen 532
1244 gezelter 117 #ifdef IS_MPI
1245     !!distribute forces
1246 chrisfen 532
1247 gezelter 117 f_temp = 0.0_dp
1248     call scatter(f_Row,f_temp,plan_atom_row_3d)
1249     do i = 1,nlocal
1250     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1251     end do
1252 chrisfen 532
1253 gezelter 117 f_temp = 0.0_dp
1254     call scatter(f_Col,f_temp,plan_atom_col_3d)
1255     do i = 1,nlocal
1256     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1257     end do
1258 chrisfen 532
1259 gezelter 141 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1260 gezelter 117 t_temp = 0.0_dp
1261     call scatter(t_Row,t_temp,plan_atom_row_3d)
1262     do i = 1,nlocal
1263     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1264     end do
1265     t_temp = 0.0_dp
1266     call scatter(t_Col,t_temp,plan_atom_col_3d)
1267 chrisfen 532
1268 gezelter 117 do i = 1,nlocal
1269     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1270     end do
1271     endif
1272 chrisfen 532
1273 gezelter 117 if (do_pot) then
1274     ! scatter/gather pot_row into the members of my column
1275 gezelter 662 do i = 1,LR_POT_TYPES
1276 chuckv 657 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1277     end do
1278 gezelter 117 ! scatter/gather pot_local into all other procs
1279     ! add resultant to get total pot
1280     do i = 1, nlocal
1281 gezelter 662 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1282     + pot_Temp(1:LR_POT_TYPES,i)
1283 gezelter 117 enddo
1284 chrisfen 532
1285 chuckv 1245 do i = 1,LR_POT_TYPES
1286     particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1287     enddo
1288    
1289 gezelter 117 pot_Temp = 0.0_DP
1290 chuckv 1245
1291 gezelter 662 do i = 1,LR_POT_TYPES
1292 chuckv 657 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1293     end do
1294 chuckv 1245
1295 gezelter 117 do i = 1, nlocal
1296 gezelter 662 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1297     + pot_Temp(1:LR_POT_TYPES,i)
1298 gezelter 117 enddo
1299 chrisfen 532
1300 chuckv 1245 do i = 1,LR_POT_TYPES
1301     particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1302     enddo
1303 gezelter 1309
1304     ppot_Temp = 0.0_DP
1305    
1306     call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
1307     do i = 1, nlocal
1308     particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1309     enddo
1310 chuckv 1245
1311 gezelter 1309 ppot_Temp = 0.0_DP
1312 chuckv 1245
1313 gezelter 1309 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
1314     do i = 1, nlocal
1315     particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1316     enddo
1317    
1318    
1319 gezelter 117 endif
1320     #endif
1321 chrisfen 532
1322 chrisfen 691 if (SIM_requires_postpair_calc) then
1323 chrisfen 695 do i = 1, nlocal
1324    
1325     ! we loop only over the local atoms, so we don't need row and column
1326     ! lookups for the types
1327 gezelter 1346
1328 chrisfen 691 me_i = atid(i)
1329    
1330 chrisfen 695 ! is the atom electrostatic? See if it would have an
1331     ! electrostatic interaction with itself
1332     iHash = InteractionHash(me_i,me_i)
1333 chrisfen 699
1334 chrisfen 691 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1335 gezelter 1340
1336     ! loop over the excludes to accumulate charge in the
1337     ! cutoff sphere that we've left out of the normal pair loop
1338     skch = 0.0_dp
1339 gezelter 1345
1340 gezelter 1346 do i1 = 1, nSkipsForLocalAtom(i)
1341     j = skipsForLocalAtom(i, i1)
1342 gezelter 1340 me_j = atid(j)
1343 gezelter 1352 jHash = InteractionHash(me_i,me_j)
1344     if ( iand(jHash, ELECTROSTATIC_PAIR).ne.0 ) then
1345     skch = skch + getCharge(me_j)
1346     endif
1347 gezelter 1340 enddo
1348 gezelter 1346
1349 gezelter 117 #ifdef IS_MPI
1350 gezelter 1340 call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), &
1351 chrisfen 695 t, do_pot)
1352 gezelter 117 #else
1353 gezelter 1340 call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), &
1354 chrisfen 695 t, do_pot)
1355 gezelter 117 #endif
1356 chrisfen 691 endif
1357 chrisfen 699
1358 chrisfen 703
1359 chrisfen 708 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1360 chrisfen 699
1361 chrisfen 703 ! loop over the excludes to accumulate RF stuff we've
1362     ! left out of the normal pair loop
1363    
1364 gezelter 1346 do i1 = 1, nSkipsForLocalAtom(i)
1365     j = skipsForLocalAtom(i, i1)
1366 chrisfen 703
1367     ! prevent overcounting of the skips
1368     if (i.lt.j) then
1369 gezelter 939 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1370 gezelter 960 rVal = sqrt(ratmsq)
1371 gezelter 939 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1372 chrisfen 699 #ifdef IS_MPI
1373 gezelter 1286 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1374 chrisfen 703 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1375 chrisfen 699 #else
1376 gezelter 1286 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1377 chrisfen 703 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1378 chrisfen 699 #endif
1379 chrisfen 703 endif
1380     enddo
1381 chrisfen 708 endif
1382 chrisfen 998
1383     if (do_box_dipole) then
1384     #ifdef IS_MPI
1385     call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1386     nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1387     pChgCount_local, nChgCount_local)
1388     #else
1389     call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1390     pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1391     #endif
1392     endif
1393 chrisfen 703 enddo
1394 gezelter 117 endif
1395 chrisfen 998
1396 gezelter 117 #ifdef IS_MPI
1397     if (do_pot) then
1398 gezelter 962 #ifdef SINGLE_PRECISION
1399     call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1400     mpi_comm_world,mpi_err)
1401     #else
1402 chrisfen 998 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1403     mpi_sum, mpi_comm_world,mpi_err)
1404 gezelter 962 #endif
1405 gezelter 117 endif
1406 gezelter 1126
1407 chrisfen 998 if (do_box_dipole) then
1408    
1409     #ifdef SINGLE_PRECISION
1410     call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1411     mpi_comm_world, mpi_err)
1412     call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1413     mpi_comm_world, mpi_err)
1414     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1415     mpi_comm_world, mpi_err)
1416     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1417     mpi_comm_world, mpi_err)
1418     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1419     mpi_comm_world, mpi_err)
1420     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1421     mpi_comm_world, mpi_err)
1422     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1423     mpi_comm_world, mpi_err)
1424 gezelter 117 #else
1425 chrisfen 998 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1426     mpi_comm_world, mpi_err)
1427     call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1428     mpi_comm_world, mpi_err)
1429     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1430     mpi_sum, mpi_comm_world, mpi_err)
1431     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1432     mpi_sum, mpi_comm_world, mpi_err)
1433     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1434     mpi_sum, mpi_comm_world, mpi_err)
1435     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1436     mpi_sum, mpi_comm_world, mpi_err)
1437     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1438     mpi_sum, mpi_comm_world, mpi_err)
1439     #endif
1440    
1441     endif
1442 chrisfen 695
1443 gezelter 117 #endif
1444 chrisfen 998
1445     if (do_box_dipole) then
1446     ! first load the accumulated dipole moment (if dipoles were present)
1447     boxDipole(1) = dipVec(1)
1448     boxDipole(2) = dipVec(2)
1449     boxDipole(3) = dipVec(3)
1450    
1451     ! now include the dipole moment due to charges
1452     ! use the lesser of the positive and negative charge totals
1453     if (nChg .le. pChg) then
1454     chg_value = nChg
1455     else
1456     chg_value = pChg
1457     endif
1458    
1459     ! find the average positions
1460     if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1461     pChgPos = pChgPos / pChgCount
1462     nChgPos = nChgPos / nChgCount
1463     endif
1464    
1465     ! dipole is from the negative to the positive (physics notation)
1466     chgVec(1) = pChgPos(1) - nChgPos(1)
1467     chgVec(2) = pChgPos(2) - nChgPos(2)
1468     chgVec(3) = pChgPos(3) - nChgPos(3)
1469    
1470     boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1471     boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1472     boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1473    
1474     endif
1475    
1476 gezelter 117 end subroutine do_force_loop
1477 chrisfen 532
1478 gezelter 117 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1479 gezelter 1309 eFrame, A, f, t, pot, particle_pot, vpair, &
1480     fpair, d_grp, r_grp, rCut, topoDist)
1481 gezelter 117
1482 chuckv 656 real( kind = dp ) :: vpair, sw
1483 gezelter 662 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1484 chuckv 1245 real( kind = dp ), dimension(nLocal) :: particle_pot
1485 gezelter 117 real( kind = dp ), dimension(3) :: fpair
1486     real( kind = dp ), dimension(nLocal) :: mfact
1487 gezelter 246 real( kind = dp ), dimension(9,nLocal) :: eFrame
1488 gezelter 117 real( kind = dp ), dimension(9,nLocal) :: A
1489     real( kind = dp ), dimension(3,nLocal) :: f
1490     real( kind = dp ), dimension(3,nLocal) :: t
1491    
1492     logical, intent(inout) :: do_pot
1493     integer, intent(in) :: i, j
1494     real ( kind = dp ), intent(inout) :: rijsq
1495 chrisfen 695 real ( kind = dp ), intent(inout) :: r_grp
1496 gezelter 117 real ( kind = dp ), intent(inout) :: d(3)
1497 chrisfen 695 real ( kind = dp ), intent(inout) :: d_grp(3)
1498 gezelter 762 real ( kind = dp ), intent(inout) :: rCut
1499 gezelter 1286 integer, intent(inout) :: topoDist
1500     real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1501 gezelter 939 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1502 gezelter 1386
1503     real( kind = dp), dimension(3) :: f1, t1, t2
1504     real( kind = dp), dimension(9) :: A1, A2, eF1, eF2
1505     real( kind = dp) :: p_vdw, p_elect, p_hb, p_met
1506     integer :: atid_i, atid_j, id1, id2, idx
1507 gezelter 939 integer :: k
1508 gezelter 117
1509 gezelter 571 integer :: iHash
1510 gezelter 560
1511 gezelter 1313 !!$ write(*,*) i, j, rijsq, d(1), d(2), d(3), sw, do_pot
1512     !!$ write(*,*) particle_pot(1), vpair, fpair(1), fpair(2), fpair(3)
1513     !!$ write(*,*) rCut
1514    
1515 chrisfen 942 r = sqrt(rijsq)
1516    
1517 gezelter 960 vpair = 0.0_dp
1518     fpair(1:3) = 0.0_dp
1519 gezelter 117
1520 gezelter 1386 p_vdw = 0.0
1521     p_elect = 0.0
1522     p_hb = 0.0
1523     p_met = 0.0
1524    
1525     f1(1:3) = 0.0
1526     t1(1:3) = 0.0
1527     t2(1:3) = 0.0
1528    
1529 gezelter 117 #ifdef IS_MPI
1530 gezelter 1386 atid_i = atid_row(i)
1531     atid_j = atid_col(j)
1532    
1533     do idx = 1, 9
1534     A1(idx) = A_Row(idx, i)
1535     A2(idx) = A_Col(idx, j)
1536     eF1(idx) = eFrame_Row(idx, i)
1537     eF2(idx) = eFrame_Col(idx, j)
1538     enddo
1539    
1540 gezelter 117 #else
1541 gezelter 1386 atid_i = atid(i)
1542     atid_j = atid(j)
1543     do idx = 1, 9
1544     A1(idx) = A(idx, i)
1545     A2(idx) = A(idx, j)
1546     eF1(idx) = eFrame(idx, i)
1547     eF2(idx) = eFrame(idx, j)
1548     enddo
1549    
1550 gezelter 117 #endif
1551 gezelter 202
1552 gezelter 1386 iHash = InteractionHash(atid_i, atid_j)
1553 cli2 1289
1554 gezelter 1286 vdwMult = vdwScale(topoDist)
1555     electroMult = electrostaticScale(topoDist)
1556 cli2 1289
1557 chrisfen 703 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1558 gezelter 1386 call do_lj_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1559     p_vdw, f1, do_pot)
1560 gezelter 117 endif
1561 chrisfen 532
1562 chrisfen 703 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1563 gezelter 1386 call doElectrostaticPair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, &
1564     vpair, fpair, p_elect, eF1, eF2, f1, t1, t2, do_pot)
1565 chrisfen 703 endif
1566    
1567     if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1568 gezelter 1386 call do_sticky_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1569     p_hb, A1, A2, f1, t1, t2, do_pot)
1570 chrisfen 703 endif
1571    
1572     if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1573 gezelter 1386 call do_sticky_power_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1574     p_hb, A1, A2, f1, t1, t2, do_pot)
1575 chrisfen 703 endif
1576    
1577     if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1578 gezelter 1386 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1579     p_vdw, A1, A2, f1, t1, t2, do_pot)
1580 chrisfen 703 endif
1581    
1582     if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1583 gezelter 1386 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1584     p_vdw, A1, A2, f1, t1, t2, do_pot)
1585 chrisfen 703 endif
1586    
1587     if ( iand(iHash, EAM_PAIR).ne.0 ) then
1588 gezelter 1386 call do_eam_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, particle_pot, &
1589     fpair, p_met, f1, do_pot)
1590 chrisfen 703 endif
1591    
1592     if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1593 gezelter 1386 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1594     p_vdw, A1, A2, f1, t1, t2, do_pot)
1595 chrisfen 703 endif
1596    
1597     if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1598 gezelter 1386 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1599     p_vdw, A1, A2, f1, t1, t2, do_pot)
1600 chrisfen 703 endif
1601 chuckv 733
1602     if ( iand(iHash, SC_PAIR).ne.0 ) then
1603 gezelter 1386 call do_SC_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vpair, particle_pot, &
1604     fpair, p_met, f1, do_pot)
1605 chuckv 733 endif
1606 chrisfen 703
1607 gezelter 1174 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1608 gezelter 1386 call do_mnm_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1609     p_vdw, A1, A2, f1, t1, t2, do_pot)
1610 gezelter 1174 endif
1611 gezelter 1386
1612    
1613     #ifdef IS_MPI
1614     id1 = AtomRowToGlobal(i)
1615     id2 = AtomColToGlobal(j)
1616    
1617     pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*p_vdw
1618     pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*p_vdw
1619     pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*p_elect
1620     pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*p_elect
1621     pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*p_hb
1622     pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*p_hb
1623     pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*p_met
1624     pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*p_met
1625    
1626     do idx = 1, 3
1627     f_Row(idx,i) = f_Row(idx,i) + f1(idx)
1628     f_Col(idx,j) = f_Col(idx,j) - f1(idx)
1629    
1630     t_Row(idx,i) = t_Row(idx,i) + t1(idx)
1631     t_Col(idx,j) = t_Col(idx,j) + t2(idx)
1632     enddo
1633     #else
1634     id1 = i
1635     id2 = j
1636    
1637     pot(VDW_POT) = pot(VDW_POT) + p_vdw
1638     pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + p_elect
1639     pot(HB_POT) = pot(HB_POT) + p_hb
1640     pot(METALLIC_POT) = pot(METALLIC_POT) + p_met
1641    
1642     do idx = 1, 3
1643     f(idx,i) = f(idx,i) + f1(idx)
1644     f(idx,j) = f(idx,j) - f1(idx)
1645    
1646     t(idx,i) = t(idx,i) + t1(idx)
1647     t(idx,j) = t(idx,j) + t2(idx)
1648     enddo
1649     #endif
1650    
1651     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1652    
1653     fpair(1) = fpair(1) + f1(1)
1654     fpair(2) = fpair(2) + f1(2)
1655     fpair(3) = fpair(3) + f1(3)
1656    
1657     endif
1658    
1659    
1660 gezelter 1313 !!$
1661     !!$ particle_pot(i) = particle_pot(i) + vpair*sw
1662     !!$ particle_pot(j) = particle_pot(j) + vpair*sw
1663 gezelter 1174
1664 gezelter 117 end subroutine do_pair
1665    
1666 gezelter 762 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1667 gezelter 246 do_pot, do_stress, eFrame, A, f, t, pot)
1668 gezelter 117
1669 chuckv 656 real( kind = dp ) :: sw
1670 gezelter 662 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1671 chrisfen 532 real( kind = dp ), dimension(9,nLocal) :: eFrame
1672     real (kind=dp), dimension(9,nLocal) :: A
1673     real (kind=dp), dimension(3,nLocal) :: f
1674     real (kind=dp), dimension(3,nLocal) :: t
1675 gezelter 117
1676 chrisfen 532 logical, intent(inout) :: do_pot, do_stress
1677     integer, intent(in) :: i, j
1678 gezelter 762 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1679 chrisfen 532 real ( kind = dp ) :: r, rc
1680     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1681    
1682 gezelter 1386 integer :: atid_i, atid_j, iHash
1683 chrisfen 532
1684 chrisfen 942 r = sqrt(rijsq)
1685    
1686 gezelter 117 #ifdef IS_MPI
1687 gezelter 1386 atid_i = atid_row(i)
1688     atid_j = atid_col(j)
1689 gezelter 117 #else
1690 gezelter 1386 atid_i = atid(i)
1691     atid_j = atid(j)
1692 gezelter 117 #endif
1693 chrisfen 532
1694 gezelter 1386 iHash = InteractionHash(atid_i, atid_j)
1695 chrisfen 532
1696 gezelter 571 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1697 gezelter 1386 call calc_EAM_prepair_rho(i, j, atid_i, atid_j, d, r, rijsq)
1698 chrisfen 532 endif
1699 chuckv 733
1700     if ( iand(iHash, SC_PAIR).ne.0 ) then
1701 gezelter 1386 call calc_SC_prepair_rho(i, j, atid_i, atid_j, d, r, rijsq, rcut )
1702 chuckv 733 endif
1703 gezelter 560
1704 chrisfen 532 end subroutine do_prepair
1705    
1706    
1707 gezelter 1285 subroutine do_preforce(nlocal, pot, particle_pot)
1708 chrisfen 532 integer :: nlocal
1709 gezelter 662 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1710 gezelter 1285 real( kind = dp ),dimension(nlocal) :: particle_pot
1711 chrisfen 532
1712     if (FF_uses_EAM .and. SIM_uses_EAM) then
1713 gezelter 1285 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1714 chrisfen 532 endif
1715 chuckv 733 if (FF_uses_SC .and. SIM_uses_SC) then
1716 gezelter 1285 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1717 chuckv 733 endif
1718 chrisfen 532 end subroutine do_preforce
1719    
1720    
1721     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1722    
1723     real (kind = dp), dimension(3) :: q_i
1724     real (kind = dp), dimension(3) :: q_j
1725     real ( kind = dp ), intent(out) :: r_sq
1726     real( kind = dp ) :: d(3), scaled(3)
1727     integer i
1728    
1729 gezelter 938 d(1) = q_j(1) - q_i(1)
1730     d(2) = q_j(2) - q_i(2)
1731     d(3) = q_j(3) - q_i(3)
1732 chrisfen 532
1733     ! Wrap back into periodic box if necessary
1734     if ( SIM_uses_PBC ) then
1735    
1736     if( .not.boxIsOrthorhombic ) then
1737     ! calc the scaled coordinates.
1738 gezelter 939 ! scaled = matmul(HmatInv, d)
1739 chrisfen 532
1740 gezelter 939 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1741     scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1742     scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1743    
1744 chrisfen 532 ! wrap the scaled coordinates
1745    
1746 gezelter 960 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1747     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1748     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1749 chrisfen 532
1750     ! calc the wrapped real coordinates from the wrapped scaled
1751     ! coordinates
1752 gezelter 939 ! d = matmul(Hmat,scaled)
1753     d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1754     d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1755     d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1756 chrisfen 532
1757     else
1758     ! calc the scaled coordinates.
1759    
1760 gezelter 938 scaled(1) = d(1) * HmatInv(1,1)
1761     scaled(2) = d(2) * HmatInv(2,2)
1762     scaled(3) = d(3) * HmatInv(3,3)
1763    
1764     ! wrap the scaled coordinates
1765    
1766 gezelter 960 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1767     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1768     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1769 chrisfen 532
1770 gezelter 938 ! calc the wrapped real coordinates from the wrapped scaled
1771     ! coordinates
1772 chrisfen 532
1773 gezelter 938 d(1) = scaled(1)*Hmat(1,1)
1774     d(2) = scaled(2)*Hmat(2,2)
1775     d(3) = scaled(3)*Hmat(3,3)
1776 chrisfen 532
1777     endif
1778    
1779     endif
1780    
1781 gezelter 938 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1782 chrisfen 532
1783     end subroutine get_interatomic_vector
1784    
1785     subroutine zero_work_arrays()
1786    
1787 gezelter 117 #ifdef IS_MPI
1788    
1789 chrisfen 532 q_Row = 0.0_dp
1790     q_Col = 0.0_dp
1791    
1792     q_group_Row = 0.0_dp
1793     q_group_Col = 0.0_dp
1794    
1795     eFrame_Row = 0.0_dp
1796     eFrame_Col = 0.0_dp
1797    
1798     A_Row = 0.0_dp
1799     A_Col = 0.0_dp
1800    
1801     f_Row = 0.0_dp
1802     f_Col = 0.0_dp
1803     f_Temp = 0.0_dp
1804    
1805     t_Row = 0.0_dp
1806     t_Col = 0.0_dp
1807     t_Temp = 0.0_dp
1808    
1809     pot_Row = 0.0_dp
1810     pot_Col = 0.0_dp
1811     pot_Temp = 0.0_dp
1812 gezelter 1309 ppot_Temp = 0.0_dp
1813 chrisfen 532
1814 gezelter 117 #endif
1815 chrisfen 532
1816     if (FF_uses_EAM .and. SIM_uses_EAM) then
1817     call clean_EAM()
1818     endif
1819    
1820     end subroutine zero_work_arrays
1821    
1822     function skipThisPair(atom1, atom2) result(skip_it)
1823     integer, intent(in) :: atom1
1824     integer, intent(in), optional :: atom2
1825     logical :: skip_it
1826     integer :: unique_id_1, unique_id_2
1827     integer :: me_i,me_j
1828     integer :: i
1829    
1830     skip_it = .false.
1831    
1832     !! there are a number of reasons to skip a pair or a particle
1833     !! mostly we do this to exclude atoms who are involved in short
1834     !! range interactions (bonds, bends, torsions), but we also need
1835     !! to exclude some overcounted interactions that result from
1836     !! the parallel decomposition
1837    
1838 gezelter 117 #ifdef IS_MPI
1839 chrisfen 532 !! in MPI, we have to look up the unique IDs for each atom
1840     unique_id_1 = AtomRowToGlobal(atom1)
1841     unique_id_2 = AtomColToGlobal(atom2)
1842     !! this situation should only arise in MPI simulations
1843     if (unique_id_1 == unique_id_2) then
1844     skip_it = .true.
1845     return
1846     end if
1847    
1848     !! this prevents us from doing the pair on multiple processors
1849     if (unique_id_1 < unique_id_2) then
1850     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1851     skip_it = .true.
1852     return
1853     endif
1854     else
1855     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1856     skip_it = .true.
1857     return
1858     endif
1859     endif
1860 gezelter 1286 #else
1861     !! in the normal loop, the atom numbers are unique
1862     unique_id_1 = atom1
1863     unique_id_2 = atom2
1864 gezelter 117 #endif
1865 gezelter 1346
1866     #ifdef IS_MPI
1867     do i = 1, nSkipsForRowAtom(atom1)
1868     if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then
1869 chrisfen 532 skip_it = .true.
1870     return
1871     endif
1872     end do
1873 gezelter 1346 #else
1874     do i = 1, nSkipsForLocalAtom(atom1)
1875     if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then
1876     skip_it = .true.
1877     return
1878     endif
1879     end do
1880     #endif
1881 chrisfen 532
1882     return
1883     end function skipThisPair
1884    
1885 gezelter 1286 function getTopoDistance(atom1, atom2) result(topoDist)
1886     integer, intent(in) :: atom1
1887     integer, intent(in) :: atom2
1888     integer :: topoDist
1889     integer :: unique_id_2
1890     integer :: i
1891    
1892     #ifdef IS_MPI
1893     unique_id_2 = AtomColToGlobal(atom2)
1894     #else
1895     unique_id_2 = atom2
1896     #endif
1897    
1898     ! zero is default for unconnected (i.e. normal) pair interactions
1899    
1900     topoDist = 0
1901    
1902     do i = 1, nTopoPairsForAtom(atom1)
1903     if (toposForAtom(atom1, i) .eq. unique_id_2) then
1904     topoDist = topoDistance(atom1, i)
1905     return
1906     endif
1907     end do
1908    
1909     return
1910     end function getTopoDistance
1911    
1912 chrisfen 532 function FF_UsesDirectionalAtoms() result(doesit)
1913     logical :: doesit
1914 gezelter 571 doesit = FF_uses_DirectionalAtoms
1915 chrisfen 532 end function FF_UsesDirectionalAtoms
1916    
1917     function FF_RequiresPrepairCalc() result(doesit)
1918     logical :: doesit
1919 chuckv 1162 doesit = FF_uses_EAM .or. FF_uses_SC
1920 chrisfen 532 end function FF_RequiresPrepairCalc
1921    
1922 gezelter 117 #ifdef PROFILE
1923 chrisfen 532 function getforcetime() result(totalforcetime)
1924     real(kind=dp) :: totalforcetime
1925     totalforcetime = forcetime
1926     end function getforcetime
1927 gezelter 117 #endif
1928    
1929 chrisfen 532 !! This cleans componets of force arrays belonging only to fortran
1930    
1931 gezelter 1126 subroutine add_stress_tensor(dpair, fpair, tau)
1932 chrisfen 532
1933     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1934 gezelter 1126 real( kind = dp ), dimension(9), intent(inout) :: tau
1935 chrisfen 532
1936     ! because the d vector is the rj - ri vector, and
1937     ! because fx, fy, fz are the force on atom i, we need a
1938     ! negative sign here:
1939    
1940 gezelter 1126 tau(1) = tau(1) - dpair(1) * fpair(1)
1941     tau(2) = tau(2) - dpair(1) * fpair(2)
1942     tau(3) = tau(3) - dpair(1) * fpair(3)
1943     tau(4) = tau(4) - dpair(2) * fpair(1)
1944     tau(5) = tau(5) - dpair(2) * fpair(2)
1945     tau(6) = tau(6) - dpair(2) * fpair(3)
1946     tau(7) = tau(7) - dpair(3) * fpair(1)
1947     tau(8) = tau(8) - dpair(3) * fpair(2)
1948     tau(9) = tau(9) - dpair(3) * fpair(3)
1949 chrisfen 532
1950     end subroutine add_stress_tensor
1951    
1952 gezelter 117 end module doForces