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 623 by chuckv, Wed Jul 16 16:40:03 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.9 2003-07-16 16:40:03 chuckv Exp $, $Date: 2003-07-16 16:40:03 $, $Name: not supported by cvs2svn $, $Revision: 1.9 $
6  
7   module lj
8    use definitions
# Line 21 | Line 21 | module lj
21   #include "fForceField.h"
22  
23    integer, save :: LJ_Mixing_Policy
24 <  integer, save :: LJ_rcut
24 >  real(kind=DP), save :: LJ_rcut
25    
26    !! Logical has lj force field module been initialized?
27    
# Line 85 | Line 85 | contains
85    subroutine LJ_new_rcut(rcut, status)
86      integer :: status, myStatus
87      real(kind=dp) :: rcut
88 <
89 <    LJ_rcut = rcut
88 >    
89      status = 0
91    call createMixingList(myStatus)
92    if (myStatus /= 0) then
93       status = -1
94       return
95    end if
90  
91 +    if ( rcut .ne. LJ_rcut ) then
92 +       LJ_rcut = rcut
93 +       call createMixingList(myStatus)
94 +       if (myStatus /= 0) then
95 +          status = -1
96 +          return
97 +       end if
98 +    endif
99 +    
100 +    
101      return
102    end subroutine LJ_new_rcut
103    
# Line 115 | Line 119 | contains
119          
120      if (.not. associated(ljMixed)) then
121         allocate(ljMixed(nAtypes, nAtypes))
122 <    else
119 <       status = -1
120 <       return
121 <    end if
122 >    endif
123  
124      rcut6 = LJ_rcut**6
125  
# Line 195 | Line 196 | contains
196      real( kind = dp ) :: t6
197      real( kind = dp ) :: t12
198      real( kind = dp ) :: delta
199 +    integer :: id1, id2
200  
201  
202      if (rij.lt.LJ_rcut)  then
# Line 254 | Line 256 | contains
256   #endif
257        
258         if (do_stress) 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 + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
259 >
260 > #ifdef IS_MPI
261 >          id1 = tagRow(atom1)
262 >          id2 = tagColumn(atom2)
263 > #else
264 >          id1 = atom1
265 >          id2 = atom2
266 > #endif
267 >
268 >          if (molMembershipList(id1) .ne. molMembershipList(id2)) then
269 >
270 >             ! because the d vector is the rj - ri vector, and
271 >             ! because fx, fy, fz are the force on atom i, we need a
272 >             ! negative sign here:
273 >
274 >             tau_Temp(1) = tau_Temp(1) - d(1) * fx
275 >             tau_Temp(2) = tau_Temp(2) - d(1) * fy
276 >             tau_Temp(3) = tau_Temp(3) - d(1) * fz
277 >             tau_Temp(4) = tau_Temp(4) - d(2) * fx
278 >             tau_Temp(5) = tau_Temp(5) - d(2) * fy
279 >             tau_Temp(6) = tau_Temp(6) - d(2) * fz
280 >             tau_Temp(7) = tau_Temp(7) - d(3) * fx
281 >             tau_Temp(8) = tau_Temp(8) - d(3) * fy
282 >             tau_Temp(9) = tau_Temp(9) - d(3) * fz
283 >            
284 >             virial_Temp = virial_Temp + &
285 >                  (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
286 >
287 >          endif
288         endif
289        
290      endif

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines