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

Comparing trunk/OOPSE/libmdtools/calc_dipole_dipole.F90 (file contents):
Revision 727 by tim, Wed Aug 27 16:16:01 2003 UTC vs.
Revision 843 by mmeineke, Wed Oct 29 20:41:39 2003 UTC

# Line 26 | Line 26 | contains
26      real(kind=dp), intent(in) :: this_rrf, this_rt
27      rrf = this_rrf
28      rt = this_rt    
29 <     ! pre converts from mu in units of debye to kcal/mol
29 >
30 >      ! pre converts from mu in units of debye to kcal/mol
31      pre = 14.38362_dp
32  
33      dipole_initialized = .true.
# Line 176 | Line 177 | contains
177         t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
178   #endif
179  
180 +       if (do_stress) then          
181 +
182   #ifdef IS_MPI
183            id1 = tagRow(atom1)
184            id2 = tagColumn(atom2)
185   #else
186            id1 = atom1
187            id2 = atom2
188 < #endif      
186 <       if (do_stress) then          
187 <
188 > #endif                
189            if (molMembershipList(id1) .ne. molMembershipList(id2)) then
190  
191               ! because the d vector is the rj - ri vector, and

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines