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

Comparing trunk/OOPSE/libmdtools/calc_sticky_pair.F90 (file contents):
Revision 1160 by gezelter, Tue May 11 21:31:15 2004 UTC vs.
Revision 1192 by gezelter, Mon May 24 21:03:30 2004 UTC

# Line 9 | Line 9
9   !! @author Matthew Meineke
10   !! @author Christopher Fennel
11   !! @author J. Daniel Gezelter
12 < !! @version $Id: calc_sticky_pair.F90,v 1.18 2004-05-11 21:31:14 gezelter Exp $, $Date: 2004-05-11 21:31:14 $, $Name: not supported by cvs2svn $, $Revision: 1.18 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.19 2004-05-24 21:03:25 gezelter Exp $, $Date: 2004-05-24 21:03:25 $, $Name: not supported by cvs2svn $, $Revision: 1.19 $
13  
14   module sticky_pair
15  
# Line 74 | Line 74 | contains
74      return
75    end subroutine set_sticky_params
76  
77 <  subroutine do_sticky_pair(atom1, atom2, d, rij, r2, sw, vpair, pot, A,f, t, &
78 <       do_pot, do_stress)
77 >  subroutine do_sticky_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
78 >       pot, A, f, t, do_pot)
79      
80      !! This routine does only the sticky portion of the SSD potential
81      !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
# Line 89 | Line 89 | contains
89      integer, intent(in) :: atom1, atom2
90      real (kind=dp), intent(inout) :: rij, r2
91      real (kind=dp), dimension(3), intent(in) :: d
92 +    real (kind=dp), dimension(3), intent(inout) :: fpair
93      real (kind=dp) :: pot, vpair, sw
94      real (kind=dp), dimension(9,nLocal) :: A
95      real (kind=dp), dimension(3,nLocal) :: f
96      real (kind=dp), dimension(3,nLocal) :: t
97 <    logical, intent(in) :: do_pot, do_stress
97 >    logical, intent(in) :: do_pot
98  
99      real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
100      real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
# Line 347 | Line 348 | contains
348         f(3,atom2) = f(3,atom2) - fzradial
349   #endif
350  
350       if (do_stress) then          
351
351   #ifdef IS_MPI
352 <          id1 = tagRow(atom1)
353 <          id2 = tagColumn(atom2)
352 >       id1 = tagRow(atom1)
353 >       id2 = tagColumn(atom2)
354   #else
355 <          id1 = atom1
356 <          id2 = atom2
355 >       id1 = atom1
356 >       id2 = atom2
357   #endif
358 <
359 <          if (molMembershipList(id1) .ne. molMembershipList(id2)) then
360 <
361 <             ! because the d vector is the rj - ri vector, and
362 <             ! because fxradial, fyradial, and fzradial are the
363 <             ! (positive) force on atom i (negative on atom j) we need
364 <             ! a negative sign here:
366 <
367 <             tau_Temp(1) = tau_Temp(1) - d(1) * fxradial
368 <             tau_Temp(2) = tau_Temp(2) - d(1) * fyradial
369 <             tau_Temp(3) = tau_Temp(3) - d(1) * fzradial
370 <             tau_Temp(4) = tau_Temp(4) - d(2) * fxradial
371 <             tau_Temp(5) = tau_Temp(5) - d(2) * fyradial
372 <             tau_Temp(6) = tau_Temp(6) - d(2) * fzradial
373 <             tau_Temp(7) = tau_Temp(7) - d(3) * fxradial
374 <             tau_Temp(8) = tau_Temp(8) - d(3) * fyradial
375 <             tau_Temp(9) = tau_Temp(9) - d(3) * fzradial
376 <
377 <             virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
378 <          endif
358 >      
359 >       if (molMembershipList(id1) .ne. molMembershipList(id2)) then
360 >          
361 >          fpair(1) = fpair(1) + fxradial
362 >          fpair(2) = fpair(2) + fyradial
363 >          fpair(3) = fpair(3) + fzradial
364 >          
365         endif
366      endif
381
367    end subroutine do_sticky_pair
368  
369    !! calculates the switching functions and their derivatives for a given

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines