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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines