ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 3561
Committed: Fri Oct 30 16:38:48 2009 UTC (15 years, 9 months ago) by chuckv
File size: 67648 byte(s)
Log Message:
Removed remaining MPI from metallic potentials. Bug in MPI calculation of energies in Sutton-Chen.

File Contents

# User Rev Content
1 gezelter 1930 !!
2 chuckv 3561 !! Copyright (c) 2005, 2009 The University of Notre Dame. All Rights Reserved.
3 gezelter 1930 !!
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 3470 !! A. Meineke, Charles F. Vardeman I, Teng Lin, Christopher
14 gezelter 1930 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42 gezelter 1610 !! doForces.F90
43     !! module doForces
44     !! Calculates Long Range forces.
45    
46     !! @author Charles F. Vardeman II
47     !! @author Matthew Meineke
48 chuckv 3561 !! @version $Id: doForces.F90,v 1.104 2009-10-30 16:38:48 chuckv Exp $, $Date: 2009-10-30 16:38:48 $, $Name: not supported by cvs2svn $, $Revision: 1.104 $
49 gezelter 1610
50 gezelter 1930
51 gezelter 1610 module doForces
52     use force_globals
53 gezelter 3441 use fForceOptions
54 gezelter 1610 use simulation
55     use definitions
56     use atype_module
57     use switcheroo
58     use neighborLists
59     use lj
60 gezelter 1930 use sticky
61 gezelter 2085 use electrostatic_module
62 gezelter 2375 use gayberne
63 chrisfen 1636 use shapes
64 gezelter 1610 use vector_class
65     use eam
66 chuckv 3164 use MetalNonMetal
67 chuckv 2432 use suttonchen
68 gezelter 1610 use status
69     #ifdef IS_MPI
70     use mpiSimulation
71     #endif
72    
73     implicit none
74     PRIVATE
75    
76     #define __FORTRAN90
77 gezelter 2273 #include "UseTheForce/fCutoffPolicy.h"
78 gezelter 2259 #include "UseTheForce/DarkSide/fInteractionMap.h"
79 chrisfen 2310 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
80 gezelter 1610
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 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 2917 logical, save :: do_box_dipole = .false.
95 chrisfen 2229
96 gezelter 1634 logical, save :: FF_uses_DirectionalAtoms
97 gezelter 2085 logical, save :: FF_uses_Dipoles
98 gezelter 1634 logical, save :: FF_uses_GayBerne
99     logical, save :: FF_uses_EAM
100 chuckv 2432 logical, save :: FF_uses_SC
101 chuckv 3164 logical, save :: FF_uses_MNM
102 chuckv 2432
103 gezelter 1634
104     logical, save :: SIM_uses_DirectionalAtoms
105     logical, save :: SIM_uses_EAM
106 chuckv 2432 logical, save :: SIM_uses_SC
107 chuckv 3164 logical, save :: SIM_uses_MNM
108 gezelter 1610 logical, save :: SIM_requires_postpair_calc
109     logical, save :: SIM_requires_prepair_calc
110     logical, save :: SIM_uses_PBC
111 gezelter 3126 logical, save :: SIM_uses_AtomicVirial
112 gezelter 1610
113 chrisfen 2306 integer, save :: electrostaticSummationMethod
114 gezelter 2461 integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
115 chrisfen 2279
116 gezelter 2461 real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut
117     real(kind=dp), save :: skinThickness
118 chrisfen 3129 logical, save :: defaultDoShiftPot
119     logical, save :: defaultDoShiftFrc
120 gezelter 2461
121 gezelter 1610 public :: init_FF
122 gezelter 2461 public :: setCutoffs
123     public :: cWasLame
124     public :: setElectrostaticMethod
125 chrisfen 2917 public :: setBoxDipole
126     public :: getBoxDipole
127 gezelter 2461 public :: setCutoffPolicy
128     public :: setSkinThickness
129 gezelter 1610 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 2260
138 gezelter 2270 !! 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 2350 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 2270 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 2273
157 chrisfen 2917 real(kind=dp), dimension(3) :: boxDipole
158 gezelter 2722
159 gezelter 1610 contains
160    
161 gezelter 2461 subroutine createInteractionHash()
162 chuckv 2260 integer :: nAtypes
163     integer :: i
164     integer :: j
165 gezelter 2270 integer :: iHash
166 tim 2267 !! Test Types
167 chuckv 2260 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 2432 logical :: i_is_SC
175 chuckv 2260 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 2432 logical :: j_is_SC
183 gezelter 2275 real(kind=dp) :: myRcut
184    
185 chuckv 2260 if (.not. associated(atypes)) then
186 gezelter 2461 call handleError("doForces", "atypes was not present before call of createInteractionHash!")
187 chuckv 2260 return
188     endif
189    
190     nAtypes = getSize(atypes)
191    
192     if (nAtypes == 0) then
193 gezelter 2461 call handleError("doForces", "nAtypes was zero during call of createInteractionHash!")
194 chuckv 2260 return
195     end if
196 chrisfen 2229
197 chuckv 2269 if (.not. allocated(InteractionHash)) then
198     allocate(InteractionHash(nAtypes,nAtypes))
199 chuckv 2354 else
200     deallocate(InteractionHash)
201     allocate(InteractionHash(nAtypes,nAtypes))
202 chuckv 2260 endif
203 gezelter 2270
204     if (.not. allocated(atypeMaxCutoff)) then
205     allocate(atypeMaxCutoff(nAtypes))
206 chuckv 2354 else
207     deallocate(atypeMaxCutoff)
208     allocate(atypeMaxCutoff(nAtypes))
209 gezelter 2270 endif
210 chuckv 2260
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 2432 call getElementProperty(atypes, i, "is_SC", i_is_SC)
220 gezelter 1610
221 chuckv 2260 do j = i, nAtypes
222 chrisfen 2229
223 chuckv 2260 iHash = 0
224     myRcut = 0.0_dp
225 gezelter 1610
226 chuckv 2260 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 2432 call getElementProperty(atypes, j, "is_SC", j_is_SC)
234 gezelter 1610
235 chuckv 2260 if (i_is_LJ .and. j_is_LJ) then
236 gezelter 2261 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 2260
247 gezelter 2261 if (i_is_StickyP .and. j_is_StickyP) then
248     iHash = ior(iHash, STICKYPOWER_PAIR)
249     endif
250 chuckv 2260
251 gezelter 2261 if (i_is_EAM .and. j_is_EAM) then
252     iHash = ior(iHash, EAM_PAIR)
253 chuckv 2260 endif
254    
255 chuckv 2432 if (i_is_SC .and. j_is_SC) then
256     iHash = ior(iHash, SC_PAIR)
257     endif
258    
259 chuckv 2260 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 3164
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 2260
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 2269 InteractionHash(i,j) = iHash
272     InteractionHash(j,i) = iHash
273 chuckv 2260
274     end do
275    
276     end do
277 tim 2267
278 gezelter 2270 haveInteractionHash = .true.
279     end subroutine createInteractionHash
280 chuckv 2260
281 gezelter 2461 subroutine createGtypeCutoffMap()
282 gezelter 2268
283 gezelter 2273 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 2530 logical :: i_is_SC
291 gezelter 2286 logical :: GtypeFound
292 chuckv 2260
293 gezelter 2275 integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend
294 chuckv 2351 integer :: n_in_i, me_i, ia, g, atom1, ja, n_in_j,me_j
295 chuckv 2288 integer :: nGroupsInRow
296 chuckv 2350 integer :: nGroupsInCol
297     integer :: nGroupTypesRow,nGroupTypesCol
298 gezelter 2461 real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol
299 gezelter 2275 real(kind=dp) :: biggestAtypeCutoff
300 gezelter 2270
301     if (.not. haveInteractionHash) then
302 gezelter 2461 call createInteractionHash()
303 chuckv 2266 endif
304 chuckv 2288 #ifdef IS_MPI
305     nGroupsInRow = getNgroupsInRow(plan_group_row)
306 chuckv 2350 nGroupsInCol = getNgroupsInCol(plan_group_col)
307 chuckv 2288 #endif
308 chuckv 2262 nAtypes = getSize(atypes)
309 chuckv 2298 ! Set all of the initial cutoffs to zero.
310     atypeMaxCutoff = 0.0_dp
311 gezelter 3470 biggestAtypeCutoff = 0.0_dp
312 gezelter 2270 do i = 1, nAtypes
313 gezelter 2281 if (SimHasAtype(i)) then
314 gezelter 2274 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 2530 call getElementProperty(atypes, i, "is_SC", i_is_SC)
322 chuckv 2298
323 chrisfen 2317 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 2530 if (i_is_SC) then
355     thisRcut = getSCCut(i)
356     if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
357     endif
358 gezelter 2274 endif
359 gezelter 2461
360 gezelter 2274 if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
361     biggestAtypeCutoff = atypeMaxCutoff(i)
362     endif
363 chrisfen 2317
364 gezelter 2273 endif
365 gezelter 2274 enddo
366 gezelter 2280
367 gezelter 2274 istart = 1
368 chuckv 2350 jstart = 1
369 gezelter 2274 #ifdef IS_MPI
370     iend = nGroupsInRow
371 chuckv 2350 jend = nGroupsInCol
372 gezelter 2274 #else
373     iend = nGroups
374 chuckv 2350 jend = nGroups
375 gezelter 2274 #endif
376 gezelter 2281
377 gezelter 2280 !! allocate the groupToGtype and gtypeMaxCutoff here.
378 chuckv 2350 if(.not.allocated(groupToGtypeRow)) then
379     ! allocate(groupToGtype(iend))
380     allocate(groupToGtypeRow(iend))
381     else
382     deallocate(groupToGtypeRow)
383     allocate(groupToGtypeRow(iend))
384 chuckv 2282 endif
385 chuckv 2350 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 2351 if(.not.associated(groupToGtypeCol)) then
403 chuckv 2350 allocate(groupToGtypeCol(jend))
404     else
405     deallocate(groupToGtypeCol)
406     allocate(groupToGtypeCol(jend))
407     end if
408    
409 tim 2532 if(.not.associated(groupMaxCutoffCol)) then
410     allocate(groupMaxCutoffCol(jend))
411 chuckv 2350 else
412 tim 2532 deallocate(groupMaxCutoffCol)
413     allocate(groupMaxCutoffCol(jend))
414 chuckv 2350 end if
415 chuckv 2351 if(.not.associated(gtypeMaxCutoffCol)) then
416 chuckv 2350 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 2281 !! 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 2756 tol = 1.0e-6_dp
435 chuckv 2350 nGroupTypesRow = 0
436 tim 2532 nGroupTypesCol = 0
437 gezelter 2280 do i = istart, iend
438 gezelter 2274 n_in_i = groupStartRow(i+1) - groupStartRow(i)
439 chuckv 2350 groupMaxCutoffRow(i) = 0.0_dp
440 gezelter 2280 do ia = groupStartRow(i), groupStartRow(i+1)-1
441     atom1 = groupListRow(ia)
442 gezelter 2274 #ifdef IS_MPI
443 gezelter 2280 me_i = atid_row(atom1)
444 gezelter 2274 #else
445 gezelter 2280 me_i = atid(atom1)
446     #endif
447 chuckv 2350 if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
448     groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
449 gezelter 2286 endif
450 gezelter 2280 enddo
451 chuckv 2350 if (nGroupTypesRow.eq.0) then
452     nGroupTypesRow = nGroupTypesRow + 1
453     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
454     groupToGtypeRow(i) = nGroupTypesRow
455 gezelter 2280 else
456 gezelter 2286 GtypeFound = .false.
457 chuckv 2350 do g = 1, nGroupTypesRow
458     if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
459     groupToGtypeRow(i) = g
460 gezelter 2286 GtypeFound = .true.
461 gezelter 2280 endif
462     enddo
463 gezelter 2286 if (.not.GtypeFound) then
464 chuckv 2350 nGroupTypesRow = nGroupTypesRow + 1
465     gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
466     groupToGtypeRow(i) = nGroupTypesRow
467 gezelter 2286 endif
468 gezelter 2280 endif
469 gezelter 2286 enddo
470    
471 chuckv 2350 #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 2280 !! allocate the gtypeCutoffMap here.
514 chuckv 2350 allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
515 gezelter 2280 !! then we do a double loop over all the group TYPES to find the cutoff
516     !! map between groups of two types
517 chuckv 2350 tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
518    
519 gezelter 2461 do i = 1, nGroupTypesRow
520 chuckv 2350 do j = 1, nGroupTypesCol
521 gezelter 2275
522 gezelter 2280 select case(cutoffPolicy)
523 gezelter 2281 case(TRADITIONAL_CUTOFF_POLICY)
524 chuckv 2350 thisRcut = tradRcut
525 gezelter 2281 case(MIX_CUTOFF_POLICY)
526 chuckv 2350 thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
527 gezelter 2281 case(MAX_CUTOFF_POLICY)
528 chuckv 2350 thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
529 gezelter 2281 case default
530     call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
531     return
532     end select
533     gtypeCutoffMap(i,j)%rcut = thisRcut
534 gezelter 2461
535     if (thisRcut.gt.largestRcut) largestRcut = thisRcut
536    
537 gezelter 2281 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
538 gezelter 2284
539 gezelter 2461 if (.not.haveSkinThickness) then
540     skinThickness = 1.0_dp
541     endif
542    
543     gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2
544    
545 chrisfen 2317 ! 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 2280 enddo
553     enddo
554 gezelter 2461
555 chuckv 2350 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 2280 haveGtypeCutoffMap = .true.
566 chrisfen 2295 end subroutine createGtypeCutoffMap
567 chrisfen 2277
568 chrisfen 3129 subroutine setCutoffs(defRcut, defRsw, defSP, defSF)
569 chrisfen 2295
570 gezelter 2461 real(kind=dp),intent(in) :: defRcut, defRsw
571 gezelter 3559 integer, intent(in) :: defSP, defSF
572 gezelter 2461 character(len = statusMsgSize) :: errMsg
573     integer :: localError
574    
575 chrisfen 2295 defaultRcut = defRcut
576     defaultRsw = defRsw
577 gezelter 3559
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 3129
589 gezelter 2461 if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then
590 chrisfen 3129 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 2461 endif
609 gezelter 2722
610 gezelter 2461 localError = 0
611 chrisfen 3129 call setLJDefaultCutoff( defaultRcut, defaultDoShiftPot, &
612     defaultDoShiftFrc )
613 gezelter 2512 call setElectrostaticCutoffRadius( defaultRcut, defaultRsw )
614 gezelter 2717 call setCutoffEAM( defaultRcut )
615     call setCutoffSC( defaultRcut )
616 chuckv 3164 call setMnMDefaultCutoff( defaultRcut, defaultDoShiftPot, &
617     defaultDoShiftFrc )
618 gezelter 2722 call set_switch(defaultRsw, defaultRcut)
619 gezelter 2592 call setHmatDangerousRcutValue(defaultRcut)
620 gezelter 2722
621 chrisfen 2317 haveDefaultCutoffs = .true.
622 gezelter 2512 haveGtypeCutoffMap = .false.
623 gezelter 2722
624 gezelter 2461 end subroutine setCutoffs
625 chrisfen 2295
626 gezelter 2461 subroutine cWasLame()
627    
628     VisitCutoffsAfterComputing = .true.
629     return
630    
631     end subroutine cWasLame
632    
633 chrisfen 2295 subroutine setCutoffPolicy(cutPolicy)
634 gezelter 2461
635 chrisfen 2295 integer, intent(in) :: cutPolicy
636 gezelter 2461
637 chrisfen 2295 cutoffPolicy = cutPolicy
638 gezelter 2461 haveCutoffPolicy = .true.
639 gezelter 2512 haveGtypeCutoffMap = .false.
640 gezelter 2461
641 gezelter 2275 end subroutine setCutoffPolicy
642 gezelter 3126
643 chrisfen 2917 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 2461 subroutine setElectrostaticMethod( thisESM )
658    
659     integer, intent(in) :: thisESM
660    
661     electrostaticSummationMethod = thisESM
662     haveElectrostaticSummationMethod = .true.
663 gezelter 2273
664 gezelter 2461 end subroutine setElectrostaticMethod
665    
666     subroutine setSkinThickness( thisSkin )
667 gezelter 2273
668 gezelter 2461 real(kind=dp), intent(in) :: thisSkin
669    
670     skinThickness = thisSkin
671 gezelter 2512 haveSkinThickness = .true.
672     haveGtypeCutoffMap = .false.
673 gezelter 2461
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 2540 SIM_uses_SC = SimUsesSC()
683 gezelter 3126 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
684 chrisfen 2917
685 gezelter 2461 haveSIMvariables = .true.
686    
687     return
688     end subroutine setSimVariables
689 gezelter 1610
690     subroutine doReadyCheck(error)
691     integer, intent(out) :: error
692     integer :: myStatus
693    
694     error = 0
695 chrisfen 2229
696 gezelter 2270 if (.not. haveInteractionHash) then
697 gezelter 2461 call createInteractionHash()
698 gezelter 1610 endif
699    
700 gezelter 2270 if (.not. haveGtypeCutoffMap) then
701 gezelter 2461 call createGtypeCutoffMap()
702 gezelter 2270 endif
703    
704 gezelter 2461 if (VisitCutoffsAfterComputing) then
705 gezelter 2722 call set_switch(largestRcut, largestRcut)
706 gezelter 2592 call setHmatDangerousRcutValue(largestRcut)
707 gezelter 2717 call setCutoffEAM(largestRcut)
708     call setCutoffSC(largestRcut)
709     VisitCutoffsAfterComputing = .false.
710 gezelter 2461 endif
711    
712 gezelter 1610 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 2722
722 gezelter 1610 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 2722
728 gezelter 1610 #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 2229
739 gezelter 2461 subroutine init_FF(thisStat)
740 gezelter 1610
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 2229
754 gezelter 1634 FF_uses_DirectionalAtoms = .false.
755     FF_uses_Dipoles = .false.
756     FF_uses_GayBerne = .false.
757 gezelter 1610 FF_uses_EAM = .false.
758 chuckv 2533 FF_uses_SC = .false.
759 chrisfen 2229
760 gezelter 1634 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 2270 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
767 chrisfen 2220
768 gezelter 1634 call getMatchingElementList(atypes, "is_GayBerne", .true., &
769     nMatches, MatchList)
770 gezelter 2270 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
771 chrisfen 2229
772 gezelter 1610 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
773     if (nMatches .gt. 0) FF_uses_EAM = .true.
774 chrisfen 2229
775 chuckv 2533 call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList)
776     if (nMatches .gt. 0) FF_uses_SC = .true.
777 gezelter 1634
778 chuckv 2533
779 gezelter 1610 haveSaneForceField = .true.
780 chrisfen 2229
781 gezelter 1610
782     if (.not. haveNeighborList) then
783     !! Create neighbor lists
784     call expandNeighborList(nLocal, my_status)
785     if (my_Status /= 0) then
786     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
787     thisStat = -1
788     return
789     endif
790     haveNeighborList = .true.
791 chrisfen 2229 endif
792    
793 gezelter 1610 end subroutine init_FF
794    
795 chrisfen 2229
796 gezelter 1610 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
797     !------------------------------------------------------------->
798 gezelter 3440 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
799 gezelter 1610 do_pot_c, do_stress_c, error)
800     !! Position array provided by C, dimensioned by getNlocal
801     real ( kind = dp ), dimension(3, nLocal) :: q
802     !! molecular center-of-mass position array
803     real ( kind = dp ), dimension(3, nGroups) :: q_group
804     !! Rotation Matrix for each long range particle in simulation.
805     real( kind = dp), dimension(9, nLocal) :: A
806     !! Unit vectors for dipoles (lab frame)
807 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
808 gezelter 1610 !! Force array provided by C, dimensioned by getNlocal
809     real ( kind = dp ), dimension(3,nLocal) :: f
810     !! Torsion array provided by C, dimensioned by getNlocal
811     real( kind = dp ), dimension(3,nLocal) :: t
812    
813     !! Stress Tensor
814     real( kind = dp), dimension(9) :: tau
815 gezelter 2361 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
816 chuckv 3397 real( kind = dp ), dimension(nLocal) :: particle_pot
817 gezelter 1610 logical ( kind = 2) :: do_pot_c, do_stress_c
818     logical :: do_pot
819     logical :: do_stress
820     logical :: in_switching_region
821     #ifdef IS_MPI
822 gezelter 2361 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
823 gezelter 1610 integer :: nAtomsInRow
824     integer :: nAtomsInCol
825     integer :: nprocs
826     integer :: nGroupsInRow
827     integer :: nGroupsInCol
828     #endif
829     integer :: natoms
830     logical :: update_nlist
831     integer :: i, j, jstart, jend, jnab
832     integer :: istart, iend
833     integer :: ia, jb, atom1, atom2
834     integer :: nlist
835 gezelter 3126 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
836 gezelter 1610 real( kind = DP ) :: sw, dswdr, swderiv, mf
837 chrisfen 2398 real( kind = DP ) :: rVal
838 gezelter 3126 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
839     real(kind=dp) :: rfpot, mu_i
840 gezelter 2461 real(kind=dp):: rCut
841 gezelter 3505 integer :: me_i, me_j, n_in_i, n_in_j, iG, j1
842 gezelter 1610 logical :: is_dp_i
843     integer :: neighborListSize
844     integer :: listerror, error
845     integer :: localError
846     integer :: propPack_i, propPack_j
847     integer :: loopStart, loopEnd, loop
848 gezelter 3512 integer :: iHash, jHash
849 gezelter 3441 integer :: i1, topoDist
850 chrisfen 2229
851 chrisfen 2917 !! the variables for the box dipole moment
852     #ifdef IS_MPI
853     integer :: pChgCount_local
854     integer :: nChgCount_local
855     real(kind=dp) :: pChg_local
856     real(kind=dp) :: nChg_local
857     real(kind=dp), dimension(3) :: pChgPos_local
858     real(kind=dp), dimension(3) :: nChgPos_local
859     real(kind=dp), dimension(3) :: dipVec_local
860     #endif
861     integer :: pChgCount
862     integer :: nChgCount
863     real(kind=dp) :: pChg
864     real(kind=dp) :: nChg
865     real(kind=dp) :: chg_value
866     real(kind=dp), dimension(3) :: pChgPos
867     real(kind=dp), dimension(3) :: nChgPos
868     real(kind=dp), dimension(3) :: dipVec
869     real(kind=dp), dimension(3) :: chgVec
870 gezelter 3500 real(kind=dp) :: skch
871 chrisfen 2917
872     !! initialize box dipole variables
873     if (do_box_dipole) then
874     #ifdef IS_MPI
875     pChg_local = 0.0_dp
876     nChg_local = 0.0_dp
877     pChgCount_local = 0
878     nChgCount_local = 0
879     do i=1, 3
880     pChgPos_local = 0.0_dp
881     nChgPos_local = 0.0_dp
882     dipVec_local = 0.0_dp
883     enddo
884     #endif
885     pChg = 0.0_dp
886     nChg = 0.0_dp
887     pChgCount = 0
888     nChgCount = 0
889     chg_value = 0.0_dp
890    
891     do i=1, 3
892     pChgPos(i) = 0.0_dp
893     nChgPos(i) = 0.0_dp
894     dipVec(i) = 0.0_dp
895     chgVec(i) = 0.0_dp
896     boxDipole(i) = 0.0_dp
897     enddo
898     endif
899    
900 gezelter 1610 !! initialize local variables
901 chrisfen 2229
902 gezelter 1610 #ifdef IS_MPI
903     pot_local = 0.0_dp
904     nAtomsInRow = getNatomsInRow(plan_atom_row)
905     nAtomsInCol = getNatomsInCol(plan_atom_col)
906     nGroupsInRow = getNgroupsInRow(plan_group_row)
907     nGroupsInCol = getNgroupsInCol(plan_group_col)
908     #else
909     natoms = nlocal
910     #endif
911 chrisfen 2229
912 gezelter 1610 call doReadyCheck(localError)
913     if ( localError .ne. 0 ) then
914     call handleError("do_force_loop", "Not Initialized")
915     error = -1
916     return
917     end if
918     call zero_work_arrays()
919 chrisfen 2229
920 gezelter 1610 do_pot = do_pot_c
921     do_stress = do_stress_c
922 chrisfen 2229
923 gezelter 1610 ! Gather all information needed by all force loops:
924 chrisfen 2229
925 gezelter 1610 #ifdef IS_MPI
926 chrisfen 2229
927 gezelter 1610 call gather(q, q_Row, plan_atom_row_3d)
928     call gather(q, q_Col, plan_atom_col_3d)
929    
930     call gather(q_group, q_group_Row, plan_group_row_3d)
931     call gather(q_group, q_group_Col, plan_group_col_3d)
932 chrisfen 2229
933 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
934 gezelter 1930 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
935     call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
936 chrisfen 2229
937 gezelter 1610 call gather(A, A_Row, plan_atom_row_rotation)
938     call gather(A, A_Col, plan_atom_col_rotation)
939     endif
940 chrisfen 2229
941 gezelter 1610 #endif
942 chrisfen 2229
943 gezelter 1610 !! Begin force loop timing:
944     #ifdef PROFILE
945     call cpu_time(forceTimeInitial)
946     nloops = nloops + 1
947     #endif
948 chrisfen 2229
949 gezelter 1610 loopEnd = PAIR_LOOP
950     if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
951     loopStart = PREPAIR_LOOP
952     else
953     loopStart = PAIR_LOOP
954     endif
955    
956     do loop = loopStart, loopEnd
957    
958     ! See if we need to update neighbor lists
959     ! (but only on the first time through):
960     if (loop .eq. loopStart) then
961     #ifdef IS_MPI
962 gezelter 2461 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
963 chrisfen 2229 update_nlist)
964 gezelter 1610 #else
965 gezelter 2461 call checkNeighborList(nGroups, q_group, skinThickness, &
966 chrisfen 2229 update_nlist)
967 gezelter 1610 #endif
968     endif
969 chrisfen 2229
970 gezelter 1610 if (update_nlist) then
971     !! save current configuration and construct neighbor list
972     #ifdef IS_MPI
973     call saveNeighborList(nGroupsInRow, q_group_row)
974     #else
975     call saveNeighborList(nGroups, q_group)
976     #endif
977     neighborListSize = size(list)
978     nlist = 0
979     endif
980 chrisfen 2229
981 gezelter 1610 istart = 1
982     #ifdef IS_MPI
983     iend = nGroupsInRow
984     #else
985     iend = nGroups - 1
986     #endif
987     outer: do i = istart, iend
988    
989     if (update_nlist) point(i) = nlist + 1
990 chrisfen 2229
991 gezelter 1610 n_in_i = groupStartRow(i+1) - groupStartRow(i)
992 chrisfen 2229
993 gezelter 1610 if (update_nlist) then
994     #ifdef IS_MPI
995     jstart = 1
996     jend = nGroupsInCol
997     #else
998     jstart = i+1
999     jend = nGroups
1000     #endif
1001     else
1002     jstart = point(i)
1003     jend = point(i+1) - 1
1004     ! make sure group i has neighbors
1005     if (jstart .gt. jend) cycle outer
1006     endif
1007 chrisfen 2229
1008 gezelter 1610 do jnab = jstart, jend
1009     if (update_nlist) then
1010     j = jnab
1011     else
1012     j = list(jnab)
1013     endif
1014    
1015     #ifdef IS_MPI
1016 chuckv 2266 me_j = atid_col(j)
1017 gezelter 1610 call get_interatomic_vector(q_group_Row(:,i), &
1018     q_group_Col(:,j), d_grp, rgrpsq)
1019     #else
1020 chuckv 2266 me_j = atid(j)
1021 gezelter 1610 call get_interatomic_vector(q_group(:,i), &
1022     q_group(:,j), d_grp, rgrpsq)
1023 chrisfen 2317 #endif
1024 gezelter 1610
1025 chuckv 2350 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
1026 gezelter 1610 if (update_nlist) then
1027     nlist = nlist + 1
1028 chrisfen 2229
1029 gezelter 1610 if (nlist > neighborListSize) then
1030     #ifdef IS_MPI
1031     call expandNeighborList(nGroupsInRow, listerror)
1032     #else
1033     call expandNeighborList(nGroups, listerror)
1034     #endif
1035     if (listerror /= 0) then
1036     error = -1
1037     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
1038     return
1039     end if
1040     neighborListSize = size(list)
1041     endif
1042 chrisfen 2229
1043 gezelter 1610 list(nlist) = j
1044     endif
1045 gezelter 2722
1046 chrisfen 2407 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
1047 chrisfen 2229
1048 gezelter 2461 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
1049 chrisfen 2407 if (loop .eq. PAIR_LOOP) then
1050 gezelter 2756 vij = 0.0_dp
1051 gezelter 2717 fij(1) = 0.0_dp
1052     fij(2) = 0.0_dp
1053     fij(3) = 0.0_dp
1054 chrisfen 2407 endif
1055    
1056 gezelter 2722 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
1057 chrisfen 2407
1058     n_in_j = groupStartCol(j+1) - groupStartCol(j)
1059    
1060     do ia = groupStartRow(i), groupStartRow(i+1)-1
1061 chrisfen 2402
1062 chrisfen 2407 atom1 = groupListRow(ia)
1063    
1064     inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
1065    
1066     atom2 = groupListCol(jb)
1067    
1068     if (skipThisPair(atom1, atom2)) cycle inner
1069    
1070     if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1071 gezelter 2717 d_atm(1) = d_grp(1)
1072     d_atm(2) = d_grp(2)
1073     d_atm(3) = d_grp(3)
1074 chrisfen 2407 ratmsq = rgrpsq
1075     else
1076 gezelter 1610 #ifdef IS_MPI
1077 chrisfen 2407 call get_interatomic_vector(q_Row(:,atom1), &
1078     q_Col(:,atom2), d_atm, ratmsq)
1079 gezelter 1610 #else
1080 chrisfen 2407 call get_interatomic_vector(q(:,atom1), &
1081     q(:,atom2), d_atm, ratmsq)
1082 gezelter 1610 #endif
1083 gezelter 3441 endif
1084    
1085     topoDist = getTopoDistance(atom1, atom2)
1086    
1087 chrisfen 2407 if (loop .eq. PREPAIR_LOOP) then
1088 gezelter 1610 #ifdef IS_MPI
1089 chrisfen 2407 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1090 gezelter 2461 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1091 chrisfen 2407 eFrame, A, f, t, pot_local)
1092 gezelter 1610 #else
1093 chrisfen 2407 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1094 gezelter 2461 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1095 chrisfen 2407 eFrame, A, f, t, pot)
1096 gezelter 1610 #endif
1097 chrisfen 2407 else
1098 gezelter 1610 #ifdef IS_MPI
1099 chrisfen 2407 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1100 gezelter 3466 do_pot, eFrame, A, f, t, pot_local, particle_pot, vpair, &
1101 gezelter 3441 fpair, d_grp, rgrp, rCut, topoDist)
1102 chuckv 3397 ! particle_pot will be accumulated from row & column
1103     ! arrays later
1104 gezelter 1610 #else
1105 chrisfen 2407 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1106 gezelter 3466 do_pot, eFrame, A, f, t, pot, particle_pot, vpair, &
1107 gezelter 3441 fpair, d_grp, rgrp, rCut, topoDist)
1108 gezelter 1610 #endif
1109 chrisfen 2407 vij = vij + vpair
1110 gezelter 2717 fij(1) = fij(1) + fpair(1)
1111     fij(2) = fij(2) + fpair(2)
1112     fij(3) = fij(3) + fpair(3)
1113 gezelter 3127 if (do_stress) then
1114 gezelter 3126 call add_stress_tensor(d_atm, fpair, tau)
1115     endif
1116 chrisfen 2407 endif
1117     enddo inner
1118     enddo
1119 gezelter 1610
1120 chrisfen 2407 if (loop .eq. PAIR_LOOP) then
1121     if (in_switching_region) then
1122     swderiv = vij*dswdr/rgrp
1123 chrisfen 3131 fg = swderiv*d_grp
1124    
1125     fij(1) = fij(1) + fg(1)
1126     fij(2) = fij(2) + fg(2)
1127     fij(3) = fij(3) + fg(3)
1128 chrisfen 2407
1129 chrisfen 3132 if (do_stress .and. (n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
1130 chrisfen 3131 call add_stress_tensor(d_atm, fg, tau)
1131     endif
1132    
1133 chrisfen 2407 do ia=groupStartRow(i), groupStartRow(i+1)-1
1134     atom1=groupListRow(ia)
1135     mf = mfactRow(atom1)
1136 gezelter 3126 ! fg is the force on atom ia due to cutoff group's
1137     ! presence in switching region
1138     fg = swderiv*d_grp*mf
1139 gezelter 1610 #ifdef IS_MPI
1140 gezelter 3126 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
1141     f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
1142     f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
1143 gezelter 1610 #else
1144 gezelter 3126 f(1,atom1) = f(1,atom1) + fg(1)
1145     f(2,atom1) = f(2,atom1) + fg(2)
1146     f(3,atom1) = f(3,atom1) + fg(3)
1147 gezelter 1610 #endif
1148 gezelter 3127 if (n_in_i .gt. 1) then
1149     if (do_stress.and.SIM_uses_AtomicVirial) then
1150     ! find the distance between the atom and the center of
1151     ! the cutoff group:
1152 gezelter 3126 #ifdef IS_MPI
1153 gezelter 3127 call get_interatomic_vector(q_Row(:,atom1), &
1154     q_group_Row(:,i), dag, rag)
1155 gezelter 3126 #else
1156 gezelter 3127 call get_interatomic_vector(q(:,atom1), &
1157     q_group(:,i), dag, rag)
1158 gezelter 3126 #endif
1159 gezelter 3127 call add_stress_tensor(dag,fg,tau)
1160     endif
1161 gezelter 3126 endif
1162 chrisfen 2407 enddo
1163    
1164     do jb=groupStartCol(j), groupStartCol(j+1)-1
1165     atom2=groupListCol(jb)
1166     mf = mfactCol(atom2)
1167 gezelter 3126 ! fg is the force on atom jb due to cutoff group's
1168     ! presence in switching region
1169     fg = -swderiv*d_grp*mf
1170 gezelter 1610 #ifdef IS_MPI
1171 gezelter 3126 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
1172     f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
1173     f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
1174 gezelter 1610 #else
1175 gezelter 3126 f(1,atom2) = f(1,atom2) + fg(1)
1176     f(2,atom2) = f(2,atom2) + fg(2)
1177     f(3,atom2) = f(3,atom2) + fg(3)
1178 gezelter 1610 #endif
1179 gezelter 3127 if (n_in_j .gt. 1) then
1180     if (do_stress.and.SIM_uses_AtomicVirial) then
1181     ! find the distance between the atom and the center of
1182     ! the cutoff group:
1183 gezelter 3126 #ifdef IS_MPI
1184 gezelter 3127 call get_interatomic_vector(q_Col(:,atom2), &
1185     q_group_Col(:,j), dag, rag)
1186 gezelter 3126 #else
1187 gezelter 3127 call get_interatomic_vector(q(:,atom2), &
1188     q_group(:,j), dag, rag)
1189 gezelter 3126 #endif
1190 gezelter 3127 call add_stress_tensor(dag,fg,tau)
1191     endif
1192     endif
1193 chrisfen 2407 enddo
1194     endif
1195 gezelter 3177 !if (do_stress.and.(.not.SIM_uses_AtomicVirial)) then
1196     ! call add_stress_tensor(d_grp, fij, tau)
1197     !endif
1198 gezelter 1610 endif
1199     endif
1200 chrisfen 2407 endif
1201 gezelter 1610 enddo
1202 chrisfen 2407
1203 gezelter 1610 enddo outer
1204 chrisfen 2229
1205 gezelter 1610 if (update_nlist) then
1206     #ifdef IS_MPI
1207     point(nGroupsInRow + 1) = nlist + 1
1208     #else
1209     point(nGroups) = nlist + 1
1210     #endif
1211     if (loop .eq. PREPAIR_LOOP) then
1212     ! we just did the neighbor list update on the first
1213     ! pass, so we don't need to do it
1214     ! again on the second pass
1215     update_nlist = .false.
1216     endif
1217     endif
1218 chrisfen 2229
1219 gezelter 1610 if (loop .eq. PREPAIR_LOOP) then
1220 chuckv 3133 #ifdef IS_MPI
1221 gezelter 3440 call do_preforce(nlocal, pot_local, particle_pot)
1222 chuckv 3133 #else
1223 gezelter 3440 call do_preforce(nlocal, pot, particle_pot)
1224 chuckv 3133 #endif
1225 gezelter 1610 endif
1226 chrisfen 2229
1227 gezelter 1610 enddo
1228 chrisfen 2229
1229 gezelter 1610 !! Do timing
1230     #ifdef PROFILE
1231     call cpu_time(forceTimeFinal)
1232     forceTime = forceTime + forceTimeFinal - forceTimeInitial
1233     #endif
1234 chrisfen 2229
1235 gezelter 1610 #ifdef IS_MPI
1236     !!distribute forces
1237 chrisfen 2229
1238 gezelter 1610 f_temp = 0.0_dp
1239     call scatter(f_Row,f_temp,plan_atom_row_3d)
1240     do i = 1,nlocal
1241     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1242     end do
1243 chrisfen 2229
1244 gezelter 1610 f_temp = 0.0_dp
1245     call scatter(f_Col,f_temp,plan_atom_col_3d)
1246     do i = 1,nlocal
1247     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1248     end do
1249 chrisfen 2229
1250 gezelter 1634 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1251 gezelter 1610 t_temp = 0.0_dp
1252     call scatter(t_Row,t_temp,plan_atom_row_3d)
1253     do i = 1,nlocal
1254     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1255     end do
1256     t_temp = 0.0_dp
1257     call scatter(t_Col,t_temp,plan_atom_col_3d)
1258 chrisfen 2229
1259 gezelter 1610 do i = 1,nlocal
1260     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1261     end do
1262     endif
1263 chrisfen 2229
1264 gezelter 1610 if (do_pot) then
1265     ! scatter/gather pot_row into the members of my column
1266 gezelter 2361 do i = 1,LR_POT_TYPES
1267 chuckv 2356 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1268     end do
1269 gezelter 1610 ! scatter/gather pot_local into all other procs
1270     ! add resultant to get total pot
1271     do i = 1, nlocal
1272 gezelter 2361 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1273     + pot_Temp(1:LR_POT_TYPES,i)
1274 gezelter 1610 enddo
1275 chrisfen 2229
1276 chuckv 3397 do i = 1,LR_POT_TYPES
1277     particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1278     enddo
1279    
1280 gezelter 1610 pot_Temp = 0.0_DP
1281 chuckv 3397
1282 gezelter 2361 do i = 1,LR_POT_TYPES
1283 chuckv 2356 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1284     end do
1285 chuckv 3397
1286 gezelter 1610 do i = 1, nlocal
1287 gezelter 2361 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1288     + pot_Temp(1:LR_POT_TYPES,i)
1289 gezelter 1610 enddo
1290 chrisfen 2229
1291 chuckv 3397 do i = 1,LR_POT_TYPES
1292     particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
1293     enddo
1294 gezelter 3466
1295     ppot_Temp = 0.0_DP
1296    
1297     call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
1298     do i = 1, nlocal
1299     particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1300     enddo
1301 chuckv 3397
1302 gezelter 3466 ppot_Temp = 0.0_DP
1303 chuckv 3397
1304 gezelter 3466 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
1305     do i = 1, nlocal
1306     particle_pot(i) = particle_pot(i) + ppot_Temp(i)
1307     enddo
1308    
1309    
1310 gezelter 1610 endif
1311     #endif
1312 chrisfen 2229
1313 chrisfen 2390 if (SIM_requires_postpair_calc) then
1314 chrisfen 2394 do i = 1, nlocal
1315    
1316     ! we loop only over the local atoms, so we don't need row and column
1317     ! lookups for the types
1318 gezelter 3506
1319 chrisfen 2390 me_i = atid(i)
1320    
1321 chrisfen 2394 ! is the atom electrostatic? See if it would have an
1322     ! electrostatic interaction with itself
1323     iHash = InteractionHash(me_i,me_i)
1324 chrisfen 2398
1325 chrisfen 2390 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1326 gezelter 3500
1327     ! loop over the excludes to accumulate charge in the
1328     ! cutoff sphere that we've left out of the normal pair loop
1329     skch = 0.0_dp
1330 gezelter 3505
1331 gezelter 3506 do i1 = 1, nSkipsForLocalAtom(i)
1332     j = skipsForLocalAtom(i, i1)
1333 gezelter 3500 me_j = atid(j)
1334 gezelter 3512 jHash = InteractionHash(me_i,me_j)
1335     if ( iand(jHash, ELECTROSTATIC_PAIR).ne.0 ) then
1336     skch = skch + getCharge(me_j)
1337     endif
1338 gezelter 3500 enddo
1339 gezelter 3506
1340 gezelter 1610 #ifdef IS_MPI
1341 gezelter 3500 call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), &
1342 chrisfen 2394 t, do_pot)
1343 gezelter 1610 #else
1344 gezelter 3500 call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), &
1345 chrisfen 2394 t, do_pot)
1346 gezelter 1610 #endif
1347 chrisfen 2390 endif
1348 chrisfen 2398
1349 chrisfen 2402
1350 chrisfen 2407 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1351 chrisfen 2398
1352 chrisfen 2402 ! loop over the excludes to accumulate RF stuff we've
1353     ! left out of the normal pair loop
1354    
1355 gezelter 3506 do i1 = 1, nSkipsForLocalAtom(i)
1356     j = skipsForLocalAtom(i, i1)
1357 chrisfen 2402
1358     ! prevent overcounting of the skips
1359     if (i.lt.j) then
1360 gezelter 2722 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
1361 gezelter 2756 rVal = sqrt(ratmsq)
1362 gezelter 2722 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
1363 chrisfen 2398 #ifdef IS_MPI
1364 gezelter 3441 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1365 chrisfen 2402 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1366 chrisfen 2398 #else
1367 gezelter 3441 call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, &
1368 chrisfen 2402 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1369 chrisfen 2398 #endif
1370 chrisfen 2402 endif
1371     enddo
1372 chrisfen 2407 endif
1373 chrisfen 2917
1374     if (do_box_dipole) then
1375     #ifdef IS_MPI
1376     call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, &
1377     nChg_local, pChgPos_local, nChgPos_local, dipVec_local, &
1378     pChgCount_local, nChgCount_local)
1379     #else
1380     call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, &
1381     pChgPos, nChgPos, dipVec, pChgCount, nChgCount)
1382     #endif
1383     endif
1384 chrisfen 2402 enddo
1385 gezelter 1610 endif
1386 chrisfen 2917
1387 gezelter 1610 #ifdef IS_MPI
1388     if (do_pot) then
1389 gezelter 2758 #ifdef SINGLE_PRECISION
1390     call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
1391     mpi_comm_world,mpi_err)
1392     #else
1393 chrisfen 2917 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
1394     mpi_sum, mpi_comm_world,mpi_err)
1395 gezelter 2758 #endif
1396 gezelter 1610 endif
1397 gezelter 3126
1398 chrisfen 2917 if (do_box_dipole) then
1399    
1400     #ifdef SINGLE_PRECISION
1401     call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, &
1402     mpi_comm_world, mpi_err)
1403     call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, &
1404     mpi_comm_world, mpi_err)
1405     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,&
1406     mpi_comm_world, mpi_err)
1407     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,&
1408     mpi_comm_world, mpi_err)
1409     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, &
1410     mpi_comm_world, mpi_err)
1411     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, &
1412     mpi_comm_world, mpi_err)
1413     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, &
1414     mpi_comm_world, mpi_err)
1415 gezelter 1610 #else
1416 chrisfen 2917 call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, &
1417     mpi_comm_world, mpi_err)
1418     call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, &
1419     mpi_comm_world, mpi_err)
1420     call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,&
1421     mpi_sum, mpi_comm_world, mpi_err)
1422     call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,&
1423     mpi_sum, mpi_comm_world, mpi_err)
1424     call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, &
1425     mpi_sum, mpi_comm_world, mpi_err)
1426     call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, &
1427     mpi_sum, mpi_comm_world, mpi_err)
1428     call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, &
1429     mpi_sum, mpi_comm_world, mpi_err)
1430     #endif
1431    
1432     endif
1433 chrisfen 2394
1434 gezelter 1610 #endif
1435 chrisfen 2917
1436     if (do_box_dipole) then
1437     ! first load the accumulated dipole moment (if dipoles were present)
1438     boxDipole(1) = dipVec(1)
1439     boxDipole(2) = dipVec(2)
1440     boxDipole(3) = dipVec(3)
1441    
1442     ! now include the dipole moment due to charges
1443     ! use the lesser of the positive and negative charge totals
1444     if (nChg .le. pChg) then
1445     chg_value = nChg
1446     else
1447     chg_value = pChg
1448     endif
1449    
1450     ! find the average positions
1451     if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then
1452     pChgPos = pChgPos / pChgCount
1453     nChgPos = nChgPos / nChgCount
1454     endif
1455    
1456     ! dipole is from the negative to the positive (physics notation)
1457     chgVec(1) = pChgPos(1) - nChgPos(1)
1458     chgVec(2) = pChgPos(2) - nChgPos(2)
1459     chgVec(3) = pChgPos(3) - nChgPos(3)
1460    
1461     boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value
1462     boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value
1463     boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value
1464    
1465     endif
1466    
1467 gezelter 1610 end subroutine do_force_loop
1468 chrisfen 2229
1469 gezelter 1610 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1470 gezelter 3466 eFrame, A, f, t, pot, particle_pot, vpair, &
1471     fpair, d_grp, r_grp, rCut, topoDist)
1472 gezelter 1610
1473 chuckv 2355 real( kind = dp ) :: vpair, sw
1474 gezelter 2361 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1475 chuckv 3397 real( kind = dp ), dimension(nLocal) :: particle_pot
1476 gezelter 1610 real( kind = dp ), dimension(3) :: fpair
1477     real( kind = dp ), dimension(nLocal) :: mfact
1478 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
1479 gezelter 1610 real( kind = dp ), dimension(9,nLocal) :: A
1480     real( kind = dp ), dimension(3,nLocal) :: f
1481     real( kind = dp ), dimension(3,nLocal) :: t
1482    
1483     logical, intent(inout) :: do_pot
1484     integer, intent(in) :: i, j
1485     real ( kind = dp ), intent(inout) :: rijsq
1486 chrisfen 2394 real ( kind = dp ), intent(inout) :: r_grp
1487 gezelter 1610 real ( kind = dp ), intent(inout) :: d(3)
1488 chrisfen 2394 real ( kind = dp ), intent(inout) :: d_grp(3)
1489 gezelter 2461 real ( kind = dp ), intent(inout) :: rCut
1490 gezelter 3441 integer, intent(inout) :: topoDist
1491     real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
1492 gezelter 2722 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
1493 gezelter 3559
1494     real( kind = dp), dimension(3) :: f1, t1, t2
1495     real( kind = dp), dimension(9) :: A1, A2, eF1, eF2
1496 chuckv 3561 real( kind = dp) :: dfrhodrho_i, dfrhodrho_j
1497     real( kind = dp) :: rho_i, rho_j
1498     real( kind = dp) :: fshift_i, fshift_j
1499 gezelter 3559 real( kind = dp) :: p_vdw, p_elect, p_hb, p_met
1500     integer :: atid_i, atid_j, id1, id2, idx
1501 gezelter 2722 integer :: k
1502 gezelter 1610
1503 gezelter 2270 integer :: iHash
1504 gezelter 2259
1505 gezelter 3470 !!$ write(*,*) i, j, rijsq, d(1), d(2), d(3), sw, do_pot
1506     !!$ write(*,*) particle_pot(1), vpair, fpair(1), fpair(2), fpair(3)
1507     !!$ write(*,*) rCut
1508    
1509 chrisfen 2727 r = sqrt(rijsq)
1510    
1511 gezelter 2756 vpair = 0.0_dp
1512     fpair(1:3) = 0.0_dp
1513 gezelter 1610
1514 gezelter 3559 p_vdw = 0.0
1515     p_elect = 0.0
1516     p_hb = 0.0
1517     p_met = 0.0
1518    
1519     f1(1:3) = 0.0
1520     t1(1:3) = 0.0
1521     t2(1:3) = 0.0
1522    
1523 gezelter 1610 #ifdef IS_MPI
1524 gezelter 3559 atid_i = atid_row(i)
1525     atid_j = atid_col(j)
1526    
1527     do idx = 1, 9
1528     A1(idx) = A_Row(idx, i)
1529     A2(idx) = A_Col(idx, j)
1530     eF1(idx) = eFrame_Row(idx, i)
1531     eF2(idx) = eFrame_Col(idx, j)
1532     enddo
1533    
1534 gezelter 1610 #else
1535 gezelter 3559 atid_i = atid(i)
1536     atid_j = atid(j)
1537     do idx = 1, 9
1538     A1(idx) = A(idx, i)
1539     A2(idx) = A(idx, j)
1540     eF1(idx) = eFrame(idx, i)
1541     eF2(idx) = eFrame(idx, j)
1542     enddo
1543    
1544 gezelter 1610 #endif
1545 gezelter 1706
1546 chuckv 3561
1547 gezelter 3559 iHash = InteractionHash(atid_i, atid_j)
1548 cli2 3444
1549 chuckv 3561 !! For the metallic potentials, we need to pass dF[rho]/drho since the pair calculation routines no longer are aware of parallel.
1550     if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1551     #ifdef IS_MPI
1552     dfrhodrho_i = dfrhodrho_row(i)
1553     dfrhodrho_j = dfrhodrho_col(j)
1554     rho_i = rho_row(i)
1555     rho_j = rho_col(j)
1556     #else
1557     dfrhodrho_i = dfrhodrho(i)
1558     dfrhodrho_j = dfrhodrho(j)
1559     rho_i = rho(i)
1560     rho_j = rho(j)
1561     #endif
1562     end if
1563    
1564 gezelter 3441 vdwMult = vdwScale(topoDist)
1565     electroMult = electrostaticScale(topoDist)
1566 cli2 3444
1567 chrisfen 2402 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1568 gezelter 3559 call do_lj_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1569     p_vdw, f1, do_pot)
1570 gezelter 1610 endif
1571 chrisfen 2229
1572 chrisfen 2402 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1573 gezelter 3559 call doElectrostaticPair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, &
1574     vpair, fpair, p_elect, eF1, eF2, f1, t1, t2, do_pot)
1575 chrisfen 2402 endif
1576    
1577     if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1578 gezelter 3559 call do_sticky_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1579     p_hb, A1, A2, f1, t1, t2, do_pot)
1580 chrisfen 2402 endif
1581    
1582     if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1583 gezelter 3559 call do_sticky_power_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1584     p_hb, A1, A2, f1, t1, t2, do_pot)
1585 chrisfen 2402 endif
1586    
1587     if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1588 gezelter 3559 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1589     p_vdw, A1, A2, f1, t1, t2, do_pot)
1590 chrisfen 2402 endif
1591    
1592     if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1593 gezelter 3559 call do_gb_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, fpair, &
1594     p_vdw, A1, A2, f1, t1, t2, do_pot)
1595 chrisfen 2402 endif
1596    
1597     if ( iand(iHash, EAM_PAIR).ne.0 ) then
1598 chuckv 3561 call do_eam_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, &
1599     fpair, p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j,fshift_i,fshift_j, do_pot)
1600 chrisfen 2402 endif
1601    
1602     if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1603 gezelter 3559 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1604     p_vdw, A1, A2, f1, t1, t2, do_pot)
1605 chrisfen 2402 endif
1606    
1607     if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1608 gezelter 3559 call do_shape_pair(i, j, atid_i, atid_j, d, r, rijsq, sw, vpair, fpair, &
1609     p_vdw, A1, A2, f1, t1, t2, do_pot)
1610 chrisfen 2402 endif
1611 chuckv 2432
1612     if ( iand(iHash, SC_PAIR).ne.0 ) then
1613 chuckv 3561 call do_SC_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vpair, &
1614     fpair, p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j, do_pot)
1615 chuckv 2432 endif
1616 chrisfen 2402
1617 gezelter 3177 if ( iand(iHash, MNM_PAIR).ne.0 ) then
1618 gezelter 3559 call do_mnm_pair(i, j, atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, fpair, &
1619     p_vdw, A1, A2, f1, t1, t2, do_pot)
1620 gezelter 3177 endif
1621 gezelter 3559
1622    
1623     #ifdef IS_MPI
1624     id1 = AtomRowToGlobal(i)
1625     id2 = AtomColToGlobal(j)
1626    
1627     pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*p_vdw
1628     pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*p_vdw
1629     pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*p_elect
1630     pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*p_elect
1631     pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*p_hb
1632     pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*p_hb
1633     pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*p_met
1634     pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*p_met
1635    
1636     do idx = 1, 3
1637     f_Row(idx,i) = f_Row(idx,i) + f1(idx)
1638     f_Col(idx,j) = f_Col(idx,j) - f1(idx)
1639    
1640     t_Row(idx,i) = t_Row(idx,i) + t1(idx)
1641     t_Col(idx,j) = t_Col(idx,j) + t2(idx)
1642     enddo
1643 chuckv 3561 ! particle_pot is the difference between the full potential
1644     ! and the full potential without the presence of a particular
1645     ! particle (atom1).
1646     !
1647     ! This reduces the density at other particle locations, so
1648     ! we need to recompute the density at atom2 assuming atom1
1649     ! didn't contribute. This then requires recomputing the
1650     ! density functional for atom2 as well.
1651     !
1652     ! Most of the particle_pot heavy lifting comes from the
1653     ! pair interaction, and will be handled by vpair. Parallel version.
1654    
1655     if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1656     ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j
1657     ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i
1658     end if
1659    
1660 gezelter 3559 #else
1661     id1 = i
1662     id2 = j
1663    
1664     pot(VDW_POT) = pot(VDW_POT) + p_vdw
1665     pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + p_elect
1666     pot(HB_POT) = pot(HB_POT) + p_hb
1667     pot(METALLIC_POT) = pot(METALLIC_POT) + p_met
1668    
1669     do idx = 1, 3
1670     f(idx,i) = f(idx,i) + f1(idx)
1671     f(idx,j) = f(idx,j) - f1(idx)
1672    
1673     t(idx,i) = t(idx,i) + t1(idx)
1674     t(idx,j) = t(idx,j) + t2(idx)
1675     enddo
1676 chuckv 3561 ! particle_pot is the difference between the full potential
1677     ! and the full potential without the presence of a particular
1678     ! particle (atom1).
1679     !
1680     ! This reduces the density at other particle locations, so
1681     ! we need to recompute the density at atom2 assuming atom1
1682     ! didn't contribute. This then requires recomputing the
1683     ! density functional for atom2 as well.
1684     !
1685     ! Most of the particle_pot heavy lifting comes from the
1686     ! pair interaction, and will be handled by vpair. NonParallel version.
1687    
1688     if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
1689     particle_pot(i) = particle_pot(i) - frho(j) + fshift_j
1690     particle_pot(j) = particle_pot(j) - frho(i) + fshift_i
1691     end if
1692    
1693    
1694 gezelter 3559 #endif
1695    
1696     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
1697    
1698     fpair(1) = fpair(1) + f1(1)
1699     fpair(2) = fpair(2) + f1(2)
1700     fpair(3) = fpair(3) + f1(3)
1701    
1702     endif
1703    
1704    
1705 gezelter 3470 !!$
1706     !!$ particle_pot(i) = particle_pot(i) + vpair*sw
1707     !!$ particle_pot(j) = particle_pot(j) + vpair*sw
1708 gezelter 3177
1709 gezelter 1610 end subroutine do_pair
1710    
1711 gezelter 2461 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1712 gezelter 1930 do_pot, do_stress, eFrame, A, f, t, pot)
1713 gezelter 1610
1714 chuckv 2355 real( kind = dp ) :: sw
1715 gezelter 2361 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1716 chrisfen 2229 real( kind = dp ), dimension(9,nLocal) :: eFrame
1717     real (kind=dp), dimension(9,nLocal) :: A
1718     real (kind=dp), dimension(3,nLocal) :: f
1719     real (kind=dp), dimension(3,nLocal) :: t
1720 gezelter 1610
1721 chrisfen 2229 logical, intent(inout) :: do_pot, do_stress
1722     integer, intent(in) :: i, j
1723 gezelter 2461 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1724 chrisfen 2229 real ( kind = dp ) :: r, rc
1725     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1726 chuckv 3561 real ( kind = dp ), rho_i_at_j, rho_j_at_i
1727 gezelter 3559 integer :: atid_i, atid_j, iHash
1728 chrisfen 2229
1729 chrisfen 2727 r = sqrt(rijsq)
1730    
1731 gezelter 1610 #ifdef IS_MPI
1732 gezelter 3559 atid_i = atid_row(i)
1733     atid_j = atid_col(j)
1734 gezelter 1610 #else
1735 gezelter 3559 atid_i = atid(i)
1736     atid_j = atid(j)
1737 gezelter 1610 #endif
1738 chuckv 3561 rho_i_at_j = 0.0_dp
1739     rho_j_at_i = 0.0_dp
1740    
1741 gezelter 3559 iHash = InteractionHash(atid_i, atid_j)
1742 chrisfen 2229
1743 gezelter 2270 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1744 chuckv 3561 call calc_EAM_prepair_rho(i, j, atid_i, atid_j, d, r, rho_i_at_j, rho_j_at_i, rijsq)
1745 chrisfen 2229 endif
1746 chuckv 2432
1747     if ( iand(iHash, SC_PAIR).ne.0 ) then
1748 chuckv 3561 call calc_SC_prepair_rho(i, j, atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i, rcut )
1749 chuckv 2432 endif
1750 chuckv 3561
1751    
1752    
1753     if ( iand(iHash, EAM_PAIR).ne.0 .or. iand(iHash, SC_PAIR).ne.0 ) then
1754     #ifdef IS_MPI
1755     rho_col(j) = rho_col(j) + rho_i_at_j
1756     rho_row(i) = rho_row(i) + rho_j_at_i
1757     #else
1758     rho(j) = rho(j) + rho_i_at_j
1759     rho(i) = rho(i) + rho_j_at_i
1760     #endif
1761     endif
1762 gezelter 2259
1763 chuckv 3561
1764    
1765 chrisfen 2229 end subroutine do_prepair
1766    
1767    
1768 gezelter 3440 subroutine do_preforce(nlocal, pot, particle_pot)
1769 chrisfen 2229 integer :: nlocal
1770 gezelter 2361 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1771 gezelter 3440 real( kind = dp ),dimension(nlocal) :: particle_pot
1772 chuckv 3561 integer :: sc_err = 0
1773 chrisfen 2229
1774 chuckv 3561 #ifdef IS_MPI
1775     if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_EAM .and. SIM_uses_EAM)) then
1776     call scatter(rho_row,rho,plan_atom_row,sc_err)
1777     if (sc_err /= 0 ) then
1778     call handleError("do_preforce()", "Error scattering rho_row into rho")
1779     endif
1780     call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
1781     if (sc_err /= 0 ) then
1782     call handleError("do_preforce()", "Error scattering rho_col into rho")
1783     endif
1784     rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
1785     end if
1786     #endif
1787    
1788    
1789    
1790 chrisfen 2229 if (FF_uses_EAM .and. SIM_uses_EAM) then
1791 gezelter 3440 call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1792 chrisfen 2229 endif
1793 chuckv 2432 if (FF_uses_SC .and. SIM_uses_SC) then
1794 gezelter 3440 call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot)
1795 chuckv 2432 endif
1796 chuckv 3561
1797     #ifdef IS_MPI
1798     if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_EAM .and. SIM_uses_EAM)) then
1799     !! communicate f(rho) and derivatives back into row and column arrays
1800     call gather(frho,frho_row,plan_atom_row, sc_err)
1801     if (sc_err /= 0) then
1802     call handleError("do_preforce()","MPI gather frho_row failure")
1803     endif
1804     call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
1805     if (sc_err /= 0) then
1806     call handleError("do_preforce()","MPI gather dfrhodrho_row failure")
1807     endif
1808     call gather(frho,frho_col,plan_atom_col, sc_err)
1809     if (sc_err /= 0) then
1810     call handleError("do_preforce()","MPI gather frho_col failure")
1811     endif
1812     call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
1813     if (sc_err /= 0) then
1814     call handleError("do_preforce()","MPI gather dfrhodrho_col failure")
1815     endif
1816     end if
1817     #endif
1818    
1819    
1820 chrisfen 2229 end subroutine do_preforce
1821    
1822    
1823     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1824    
1825     real (kind = dp), dimension(3) :: q_i
1826     real (kind = dp), dimension(3) :: q_j
1827     real ( kind = dp ), intent(out) :: r_sq
1828     real( kind = dp ) :: d(3), scaled(3)
1829     integer i
1830    
1831 gezelter 2717 d(1) = q_j(1) - q_i(1)
1832     d(2) = q_j(2) - q_i(2)
1833     d(3) = q_j(3) - q_i(3)
1834 chrisfen 2229
1835     ! Wrap back into periodic box if necessary
1836     if ( SIM_uses_PBC ) then
1837    
1838     if( .not.boxIsOrthorhombic ) then
1839     ! calc the scaled coordinates.
1840 gezelter 2722 ! scaled = matmul(HmatInv, d)
1841 chrisfen 2229
1842 gezelter 2722 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1843     scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1844     scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1845    
1846 chrisfen 2229 ! wrap the scaled coordinates
1847    
1848 gezelter 2756 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1849     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1850     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1851 chrisfen 2229
1852     ! calc the wrapped real coordinates from the wrapped scaled
1853     ! coordinates
1854 gezelter 2722 ! d = matmul(Hmat,scaled)
1855     d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1856     d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1857     d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1858 chrisfen 2229
1859     else
1860     ! calc the scaled coordinates.
1861    
1862 gezelter 2717 scaled(1) = d(1) * HmatInv(1,1)
1863     scaled(2) = d(2) * HmatInv(2,2)
1864     scaled(3) = d(3) * HmatInv(3,3)
1865    
1866     ! wrap the scaled coordinates
1867    
1868 gezelter 2756 scaled(1) = scaled(1) - anint(scaled(1), kind=dp)
1869     scaled(2) = scaled(2) - anint(scaled(2), kind=dp)
1870     scaled(3) = scaled(3) - anint(scaled(3), kind=dp)
1871 chrisfen 2229
1872 gezelter 2717 ! calc the wrapped real coordinates from the wrapped scaled
1873     ! coordinates
1874 chrisfen 2229
1875 gezelter 2717 d(1) = scaled(1)*Hmat(1,1)
1876     d(2) = scaled(2)*Hmat(2,2)
1877     d(3) = scaled(3)*Hmat(3,3)
1878 chrisfen 2229
1879     endif
1880    
1881     endif
1882    
1883 gezelter 2717 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1884 chrisfen 2229
1885     end subroutine get_interatomic_vector
1886    
1887     subroutine zero_work_arrays()
1888    
1889 gezelter 1610 #ifdef IS_MPI
1890    
1891 chrisfen 2229 q_Row = 0.0_dp
1892     q_Col = 0.0_dp
1893    
1894     q_group_Row = 0.0_dp
1895     q_group_Col = 0.0_dp
1896    
1897     eFrame_Row = 0.0_dp
1898     eFrame_Col = 0.0_dp
1899    
1900     A_Row = 0.0_dp
1901     A_Col = 0.0_dp
1902    
1903     f_Row = 0.0_dp
1904     f_Col = 0.0_dp
1905     f_Temp = 0.0_dp
1906    
1907     t_Row = 0.0_dp
1908     t_Col = 0.0_dp
1909     t_Temp = 0.0_dp
1910    
1911     pot_Row = 0.0_dp
1912     pot_Col = 0.0_dp
1913     pot_Temp = 0.0_dp
1914 gezelter 3466 ppot_Temp = 0.0_dp
1915 chrisfen 2229
1916 chuckv 3561 frho_row = 0.0_dp
1917     frho_col = 0.0_dp
1918     rho_row = 0.0_dp
1919     rho_col = 0.0_dp
1920     rho_tmp = 0.0_dp
1921     dfrhodrho_row = 0.0_dp
1922     dfrhodrho_col = 0.0_dp
1923    
1924 gezelter 1610 #endif
1925 chuckv 3561 rho = 0.0_dp
1926     frho = 0.0_dp
1927     dfrhodrho = 0.0_dp
1928 chrisfen 2229
1929     end subroutine zero_work_arrays
1930    
1931     function skipThisPair(atom1, atom2) result(skip_it)
1932     integer, intent(in) :: atom1
1933     integer, intent(in), optional :: atom2
1934     logical :: skip_it
1935     integer :: unique_id_1, unique_id_2
1936     integer :: me_i,me_j
1937     integer :: i
1938    
1939     skip_it = .false.
1940    
1941     !! there are a number of reasons to skip a pair or a particle
1942     !! mostly we do this to exclude atoms who are involved in short
1943     !! range interactions (bonds, bends, torsions), but we also need
1944     !! to exclude some overcounted interactions that result from
1945     !! the parallel decomposition
1946    
1947 gezelter 1610 #ifdef IS_MPI
1948 chrisfen 2229 !! in MPI, we have to look up the unique IDs for each atom
1949     unique_id_1 = AtomRowToGlobal(atom1)
1950     unique_id_2 = AtomColToGlobal(atom2)
1951     !! this situation should only arise in MPI simulations
1952     if (unique_id_1 == unique_id_2) then
1953     skip_it = .true.
1954     return
1955     end if
1956    
1957     !! this prevents us from doing the pair on multiple processors
1958     if (unique_id_1 < unique_id_2) then
1959     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1960     skip_it = .true.
1961     return
1962     endif
1963     else
1964     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1965     skip_it = .true.
1966     return
1967     endif
1968     endif
1969 gezelter 3441 #else
1970     !! in the normal loop, the atom numbers are unique
1971     unique_id_1 = atom1
1972     unique_id_2 = atom2
1973 gezelter 1610 #endif
1974 gezelter 3506
1975     #ifdef IS_MPI
1976     do i = 1, nSkipsForRowAtom(atom1)
1977     if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then
1978 chrisfen 2229 skip_it = .true.
1979     return
1980     endif
1981     end do
1982 gezelter 3506 #else
1983     do i = 1, nSkipsForLocalAtom(atom1)
1984     if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then
1985     skip_it = .true.
1986     return
1987     endif
1988     end do
1989     #endif
1990 chrisfen 2229
1991     return
1992     end function skipThisPair
1993    
1994 gezelter 3441 function getTopoDistance(atom1, atom2) result(topoDist)
1995     integer, intent(in) :: atom1
1996     integer, intent(in) :: atom2
1997     integer :: topoDist
1998     integer :: unique_id_2
1999     integer :: i
2000    
2001     #ifdef IS_MPI
2002     unique_id_2 = AtomColToGlobal(atom2)
2003     #else
2004     unique_id_2 = atom2
2005     #endif
2006    
2007     ! zero is default for unconnected (i.e. normal) pair interactions
2008    
2009     topoDist = 0
2010    
2011     do i = 1, nTopoPairsForAtom(atom1)
2012     if (toposForAtom(atom1, i) .eq. unique_id_2) then
2013     topoDist = topoDistance(atom1, i)
2014     return
2015     endif
2016     end do
2017    
2018     return
2019     end function getTopoDistance
2020    
2021 chrisfen 2229 function FF_UsesDirectionalAtoms() result(doesit)
2022     logical :: doesit
2023 gezelter 2270 doesit = FF_uses_DirectionalAtoms
2024 chrisfen 2229 end function FF_UsesDirectionalAtoms
2025    
2026     function FF_RequiresPrepairCalc() result(doesit)
2027     logical :: doesit
2028 chuckv 3164 doesit = FF_uses_EAM .or. FF_uses_SC
2029 chrisfen 2229 end function FF_RequiresPrepairCalc
2030    
2031 gezelter 1610 #ifdef PROFILE
2032 chrisfen 2229 function getforcetime() result(totalforcetime)
2033     real(kind=dp) :: totalforcetime
2034     totalforcetime = forcetime
2035     end function getforcetime
2036 gezelter 1610 #endif
2037    
2038 chrisfen 2229 !! This cleans componets of force arrays belonging only to fortran
2039    
2040 gezelter 3126 subroutine add_stress_tensor(dpair, fpair, tau)
2041 chrisfen 2229
2042     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
2043 gezelter 3126 real( kind = dp ), dimension(9), intent(inout) :: tau
2044 chrisfen 2229
2045     ! because the d vector is the rj - ri vector, and
2046     ! because fx, fy, fz are the force on atom i, we need a
2047     ! negative sign here:
2048    
2049 gezelter 3126 tau(1) = tau(1) - dpair(1) * fpair(1)
2050     tau(2) = tau(2) - dpair(1) * fpair(2)
2051     tau(3) = tau(3) - dpair(1) * fpair(3)
2052     tau(4) = tau(4) - dpair(2) * fpair(1)
2053     tau(5) = tau(5) - dpair(2) * fpair(2)
2054     tau(6) = tau(6) - dpair(2) * fpair(3)
2055     tau(7) = tau(7) - dpair(3) * fpair(1)
2056     tau(8) = tau(8) - dpair(3) * fpair(2)
2057     tau(9) = tau(9) - dpair(3) * fpair(3)
2058 chrisfen 2229
2059     end subroutine add_stress_tensor
2060    
2061 gezelter 1610 end module doForces