ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_sticky_pair.F90
(Generate patch)

Comparing trunk/OOPSE_old/src/mdtools/libmdCode/calc_sticky_pair.F90 (file contents):
Revision 324 by gezelter, Wed Mar 12 16:38:17 2003 UTC vs.
Revision 326 by gezelter, Wed Mar 12 19:31:55 2003 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.4 2003-03-12 16:38:17 gezelter Exp $, $Date: 2003-03-12 16:38:17 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.5 2003-03-12 19:31:55 gezelter Exp $, $Date: 2003-03-12 19:31:55 $, $Name: not supported by cvs2svn $, $Revision: 1.5 $
13  
14   module sticky_pair
15  
# Line 59 | Line 59 | contains
59      return
60    end subroutine set_sticky_params
61  
62 <  subroutine do_sticky_pair(atom1, atom2, d, rij, r2, A, pot, f, t)
62 >  subroutine do_sticky_pair(atom1, atom2, d, rij, r2, A, pot, f, t, &
63 >       do_pot, do_stress)
64      
65      !! This routine does only the sticky portion of the SSD potential
66      !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
# Line 77 | Line 78 | contains
78      real (kind=dp), dimension(9,getNlocal()) :: A
79      real (kind=dp), dimension(3,getNlocal()) :: f
80      real (kind=dp), dimension(3,getNlocal()) :: t
81 +    logical, intent(in) :: do_pot, do_stress
82  
83      real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
84      real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
# Line 156 | Line 158 | contains
158      wjp = zjf*zjf*zjs*zjs - SSD_w0
159      wp = wip + wjp
160  
161 +
162 +    if (do_pot) then
163   #ifdef IS_MPI
164 <    pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
165 <    pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
164 >       pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
165 >       pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
166   #else
167 <    pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
168 < #endif    
167 >       pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
168 > #endif  
169 >    endif
170      
171      dwidx =   4.0d0*xi*zi/r3  - 6.0d0*xi*zi*(xi2-yi2)/r5
172      dwidy = - 4.0d0*yi*zi/r3  - 6.0d0*yi*zi*(xi2-yi2)/r5
# Line 318 | Line 323 | contains
323           fzjj + fzij)
324   #endif
325      
326 <    if (doStress()) then          
326 >    if (do_stress) then          
327         tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
328         tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
329         tau_Temp(3) = tau_Temp(3) + fxradial * d(3)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines