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 443 by mmeineke, Wed Apr 2 22:19:03 2003 UTC vs.
Revision 473 by mmeineke, Mon Apr 7 21:20:38 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-04-02 22:19:03 mmeineke Exp $, $Date: 2003-04-02 22:19:03 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.6 2003-04-07 21:20:38 mmeineke Exp $, $Date: 2003-04-07 21:20:38 $, $Name: not supported by cvs2svn $, $Revision: 1.6 $
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 240 | Line 241 | contains
241   #endif    
242      ! Now, on to the forces:
243      
243    write(*,'(a4,3es12.3)') 't1= ', t(1:3,atom1)
244    write(*,'(a4,3es12.3)') 't2= ', t(1:3,atom2)
245    write(*,*)
246
244      ! first rotate the i terms back into the lab frame:
245  
246   #ifdef IS_MPI    
# Line 307 | 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 - &
313 <         dspdr*drdy*wp + fyjj + fyij)
314 <    f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - &
315 <         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
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 + &
324 <         fyjj + fyij)
325 <    f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
326 <         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          

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines