ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
(Generate patch)

Comparing trunk/OOPSE-4/src/UseTheForce/doForces.F90 (file contents):
Revision 3499 by gezelter, Wed Oct 22 20:01:49 2008 UTC vs.
Revision 3500 by gezelter, Wed May 6 20:51:00 2009 UTC

# Line 45 | Line 45
45  
46   !! @author Charles F. Vardeman II
47   !! @author Matthew Meineke
48 < !! @version $Id: doForces.F90,v 1.98 2008-10-22 20:01:48 gezelter Exp $, $Date: 2008-10-22 20:01:48 $, $Name: not supported by cvs2svn $, $Revision: 1.98 $
48 > !! @version $Id: doForces.F90,v 1.99 2009-05-06 20:51:00 gezelter Exp $, $Date: 2009-05-06 20:51:00 $, $Name: not supported by cvs2svn $, $Revision: 1.99 $
49  
50  
51   module doForces
# Line 868 | Line 868 | contains
868      real(kind=dp), dimension(3) :: nChgPos
869      real(kind=dp), dimension(3) :: dipVec
870      real(kind=dp), dimension(3) :: chgVec
871 +    real(kind=dp) :: skch
872  
873      !! initialize box dipole variables
874      if (do_box_dipole) then
# Line 1323 | Line 1324 | contains
1324            iHash = InteractionHash(me_i,me_i)
1325  
1326            if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1327 +
1328 +             ! loop over the excludes to accumulate charge in the
1329 +             ! cutoff sphere that we've left out of the normal pair loop
1330 +             skch = 0.0_dp
1331 +             do i1 = 1, nSkipsForAtom(i)
1332 +                j = skipsForAtom(i, i1)
1333 +                me_j = atid(j)
1334 +                skch = skch + getCharge(me_j)
1335 +             enddo
1336   #ifdef IS_MPI
1337 <             call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1337 >             call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), &
1338                    t, do_pot)
1339   #else
1340 <             call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1340 >             call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), &
1341                    t, do_pot)
1342   #endif
1343            endif

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines