ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1469
Committed: Mon Jul 19 14:07:59 2010 UTC (14 years, 9 months ago) by gezelter
File size: 68141 byte(s)
Log Message:
attempts at c++->fortran linkage.  Still busticated.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date