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 378 by mmeineke, Fri Mar 21 17:42:12 2003 UTC vs.
Revision 483 by gezelter, Wed Apr 9 04:06:43 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.1.1.1 2003-03-21 17:42:12 mmeineke Exp $, $Date: 2003-03-21 17:42:12 $, $Name: not supported by cvs2svn $, $Revision: 1.1.1.1 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.8 2003-04-09 04:06:43 gezelter Exp $, $Date: 2003-04-09 04:06:43 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $
13  
14   module sticky_pair
15  
16    use force_globals
17    use definitions
18 +  use simulation
19   #ifdef IS_MPI
20    use mpiSimulation
21   #endif
# Line 24 | Line 25 | module sticky_pair
25    PRIVATE
26  
27    logical, save :: sticky_initialized = .false.
28 <  real( kind = dp ), save :: SSD_w0
29 <  real( kind = dp ), save :: SSD_v0
30 <  real( kind = dp ), save :: SSD_rl
31 <  real( kind = dp ), save :: SSD_ru
32 <  real( kind = dp ), save :: SSD_rup
28 >  real( kind = dp ), save :: SSD_w0 = 0.0_dp
29 >  real( kind = dp ), save :: SSD_v0 = 0.0_dp
30 >  real( kind = dp ), save :: SSD_rl = 0.0_dp
31 >  real( kind = dp ), save :: SSD_ru = 0.0_dp
32 >  real( kind = dp ), save :: SSD_rup = 0.0_dp
33  
34    public :: check_sticky_FF
35    public :: set_sticky_params
# Line 74 | Line 75 | contains
75      real (kind=dp), intent(inout) :: rij, r2
76      real (kind=dp), dimension(3), intent(in) :: d
77      real (kind=dp) :: pot
78 <    real (kind=dp), dimension(:,:) :: A
79 <    real (kind=dp), dimension(:,:) :: f
80 <    real (kind=dp), dimension(:,:) :: t
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
# Line 92 | Line 93 | contains
93      real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
94      real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji      
95      real (kind=dp) :: fxradial, fyradial, fzradial
96 +    real (kind=dp) :: rijtest, rjitest
97        
98      if (.not.sticky_initialized) then
99         write(*,*) 'Sticky forces not initialized!'
# Line 140 | Line 142 | contains
142      xj2 = xj*xj
143      yj2 = yj*yj
144      zj2 = zj*zj
145 <    
145 >  
146      call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
147      
148      wi = 2.0d0*(xi2-yi2)*zi / r3
# Line 157 | Line 159 | contains
159      wjp = zjf*zjf*zjs*zjs - SSD_w0
160      wp = wip + wjp
161  
160
162      if (do_pot) then
163   #ifdef IS_MPI
164         pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
# Line 303 | Line 304 | contains
304      f_Row(2,atom1) = f_Row(2,atom1) + fyradial
305      f_Row(3,atom1) = f_Row(3,atom1) + fzradial
306      
307 <    f_Col(1,atom2) = f_Col(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - &
308 <         dspdr*drdx*wp + fxjj + fxij)
309 <    f_Col(2,atom2) = f_Col(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - &
309 <         dspdr*drdy*wp + fyjj + fyij)
310 <    f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - &
311 <         dspdr*drdz*wp + fzjj + fzij)
307 >    f_Col(1,atom2) = f_Col(1,atom2) - fxradial
308 >    f_Col(2,atom2) = f_Col(2,atom2) - fyradial
309 >    f_Col(3,atom2) = f_Col(3,atom2) - fzradial
310   #else
311      f(1,atom1) = f(1,atom1) + fxradial
312      f(2,atom1) = f(2,atom1) + fyradial
313      f(3,atom1) = f(3,atom1) + fzradial
314      
315 <    f(1,atom2) = f(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
316 <         fxjj + fxij)
317 <    f(2,atom2) = f(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
320 <         fyjj + fyij)
321 <    f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
322 <         fzjj + fzij)
315 >    f(1,atom2) = f(1,atom2) - fxradial
316 >    f(2,atom2) = f(2,atom2) - fyradial
317 >    f(3,atom2) = f(3,atom2) - fzradial
318   #endif
319      
320      if (do_stress) then          
321 <       tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
322 <       tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
323 <       tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
324 <       tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
325 <       tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
326 <       tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
327 <       tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
328 <       tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
329 <       tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
330 <       virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
321 >       if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
322 >          tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
323 >          tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
324 >          tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
325 >          tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
326 >          tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
327 >          tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
328 >          tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
329 >          tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
330 >          tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
331 >          virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
332 >       endif
333      endif
334    
335    end subroutine do_sticky_pair

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines