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 482 by chuckv, Tue Apr 8 22:38:43 2003 UTC vs.
Revision 483 by gezelter, Wed Apr 9 04:06:43 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.3 2003-04-08 22:38:43 chuckv Exp $, $Date: 2003-04-08 22:38:43 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
5 > !! @version $Id: calc_LJ_FF.F90,v 1.4 2003-04-09 04:06:43 gezelter Exp $, $Date: 2003-04-09 04:06:43 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
6  
7   module lj
8    use definitions
# Line 254 | Line 254 | contains
254   #endif
255        
256         if (do_stress) then
257 <          tau_Temp(1) = tau_Temp(1) + fx * d(1)
258 <          tau_Temp(2) = tau_Temp(2) + fx * d(2)
259 <          tau_Temp(3) = tau_Temp(3) + fx * d(3)
260 <          tau_Temp(4) = tau_Temp(4) + fy * d(1)
261 <          tau_Temp(5) = tau_Temp(5) + fy * d(2)
262 <          tau_Temp(6) = tau_Temp(6) + fy * d(3)
263 <          tau_Temp(7) = tau_Temp(7) + fz * d(1)
264 <          tau_Temp(8) = tau_Temp(8) + fz * d(2)
265 <          tau_Temp(9) = tau_Temp(9) + fz * d(3)
266 <          virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
257 >
258 >          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
259 >             tau_Temp(1) = tau_Temp(1) + fx * d(1)
260 >             tau_Temp(2) = tau_Temp(2) + fx * d(2)
261 >             tau_Temp(3) = tau_Temp(3) + fx * d(3)
262 >             tau_Temp(4) = tau_Temp(4) + fy * d(1)
263 >             tau_Temp(5) = tau_Temp(5) + fy * d(2)
264 >             tau_Temp(6) = tau_Temp(6) + fy * d(3)
265 >             tau_Temp(7) = tau_Temp(7) + fz * d(1)
266 >             tau_Temp(8) = tau_Temp(8) + fz * d(2)
267 >             tau_Temp(9) = tau_Temp(9) + fz * d(3)
268 >             virial_Temp = virial_Temp + &
269 >                  (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
270 >          endif
271 >
272         endif
273        
274      endif

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines