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

Comparing:
branches/mmeineke/OOPSE/libmdtools/calc_LJ_FF.F90 (file contents), Revision 377 by mmeineke, Fri Mar 21 17:42:12 2003 UTC vs.
trunk/OOPSE/libmdtools/calc_LJ_FF.F90 (file contents), Revision 482 by chuckv, Tue Apr 8 22:38: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.1.1.1 2003-03-21 17:42:12 mmeineke Exp $, $Date: 2003-03-21 17:42:12 $, $Name: not supported by cvs2svn $, $Revision: 1.1.1.1 $
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 $
6  
7   module lj
8    use definitions
9    use atype_module
10    use vector_class
11 +  use simulation
12   #ifdef IS_MPI
13    use mpiSimulation
14   #endif
# Line 180 | Line 181 | contains
181      integer, intent(in) ::  atom1, atom2
182      real( kind = dp ), intent(in) :: rij, r2
183      real( kind = dp ) :: pot
184 <    real( kind = dp ), dimension(:,:) :: f    
184 >    real( kind = dp ), dimension(3,getNlocal()) :: f    
185      real( kind = dp ), intent(in), dimension(3) :: d
186      logical, intent(in) :: do_pot, do_stress
187  
# Line 218 | Line 219 | contains
219        
220         dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
221        
222 <       drdx = -d(1) / rij
223 <       drdy = -d(2) / rij
224 <       drdz = -d(3) / rij
222 >       drdx = d(1) / rij
223 >       drdy = d(2) / rij
224 >       drdz = d(3) / rij
225        
226         fx = dudr * drdx
227         fy = dudr * drdy

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines