ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/calc_LJ_FF.F90 (file contents):
Revision 610 by mmeineke, Fri Apr 11 18:46:37 2003 UTC vs.
Revision 611 by gezelter, Tue Jul 15 17:10:50 2003 UTC

# Line 2 | Line 2
2   !! Corresponds to the force field defined in lj_FF.cpp
3   !! @author Charles F. Vardeman II
4   !! @author Matthew Meineke
5 < !! @version $Id: calc_LJ_FF.F90,v 1.6 2003-04-11 18:46:37 mmeineke Exp $, $Date: 2003-04-11 18:46:37 $, $Name: not supported by cvs2svn $, $Revision: 1.6 $
5 > !! @version $Id: calc_LJ_FF.F90,v 1.7 2003-07-15 17:10:50 gezelter Exp $, $Date: 2003-07-15 17:10:50 $, $Name: not supported by cvs2svn $, $Revision: 1.7 $
6  
7   module lj
8    use definitions
# Line 265 | Line 265 | contains
265   #endif
266  
267            if (molMembershipList(id1) .ne. molMembershipList(id2)) then
268 +
269 +             ! because the d vector is the rj - ri vector, and
270 +             ! because fx, fy, fz are the force on atom i, we need a
271 +             ! negative sign here:
272  
273 <             tau_Temp(1) = tau_Temp(1) + fx * d(1)
274 <             tau_Temp(2) = tau_Temp(2) + fx * d(2)
275 <             tau_Temp(3) = tau_Temp(3) + fx * d(3)
276 <             tau_Temp(4) = tau_Temp(4) + fy * d(1)
277 <             tau_Temp(5) = tau_Temp(5) + fy * d(2)
278 <             tau_Temp(6) = tau_Temp(6) + fy * d(3)
279 <             tau_Temp(7) = tau_Temp(7) + fz * d(1)
280 <             tau_Temp(8) = tau_Temp(8) + fz * d(2)
281 <             tau_Temp(9) = tau_Temp(9) + fz * d(3)
273 >             tau_Temp(1) = tau_Temp(1) - d(1) * fx
274 >             tau_Temp(2) = tau_Temp(2) - d(1) * fy
275 >             tau_Temp(3) = tau_Temp(3) - d(1) * fz
276 >             tau_Temp(4) = tau_Temp(4) - d(2) * fx
277 >             tau_Temp(5) = tau_Temp(5) - d(2) * fy
278 >             tau_Temp(6) = tau_Temp(6) - d(2) * fz
279 >             tau_Temp(7) = tau_Temp(7) - d(3) * fx
280 >             tau_Temp(8) = tau_Temp(8) - d(3) * fy
281 >             tau_Temp(9) = tau_Temp(9) - d(3) * fz
282 >            
283               virial_Temp = virial_Temp + &
284                    (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
285  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines