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 635 by gezelter, Thu Jul 17 20:32:24 2003 UTC vs.
Revision 898 by chuckv, Mon Jan 5 22:49:14 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.11 2003-07-17 20:32:24 gezelter Exp $, $Date: 2003-07-17 20:32:24 $, $Name: not supported by cvs2svn $, $Revision: 1.11 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.16 2004-01-05 22:49:14 chuckv Exp $, $Date: 2004-01-05 22:49:14 $, $Name: not supported by cvs2svn $, $Revision: 1.16 $
13  
14   module sticky_pair
15  
# Line 32 | Line 32 | module sticky_pair
32    real( kind = dp ), save :: SSD_ru = 0.0_dp
33    real( kind = dp ), save :: SSD_rlp = 0.0_dp
34    real( kind = dp ), save :: SSD_rup = 0.0_dp
35 +  real( kind = dp ), save :: SSD_rbig = 0.0_dp
36  
37    public :: check_sticky_FF
38    public :: set_sticky_params
# Line 62 | Line 63 | contains
63      SSD_ru = sticky_ru
64      SSD_rlp = sticky_rlp
65      SSD_rup = sticky_rup
66 +
67 +    if (SSD_ru .gt. SSD_rup) then
68 +       SSD_rbig = SSD_ru
69 +    else
70 +       SSD_rbig = SSD_rup
71 +    endif
72    
73      sticky_initialized = .true.
74      return
# Line 83 | Line 90 | contains
90      real (kind=dp), intent(inout) :: rij, r2
91      real (kind=dp), dimension(3), intent(in) :: d
92      real (kind=dp) :: pot
93 <    real (kind=dp), dimension(9,getNlocal()) :: A
94 <    real (kind=dp), dimension(3,getNlocal()) :: f
95 <    real (kind=dp), dimension(3,getNlocal()) :: t
93 >    real (kind=dp), dimension(9,nLocal) :: A
94 >    real (kind=dp), dimension(3,nLocal) :: f
95 >    real (kind=dp), dimension(3,nLocal) :: t
96      logical, intent(in) :: do_pot, do_stress
97  
98      real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
# Line 104 | Line 111 | contains
111      real (kind=dp) :: rijtest, rjitest
112      real (kind=dp) :: radcomxi, radcomyi, radcomzi
113      real (kind=dp) :: radcomxj, radcomyj, radcomzj
114 <
114 >    integer :: id1, id2
115  
116      if (.not.sticky_initialized) then
117         write(*,*) 'Sticky forces not initialized!'
118         return
119      endif
120  
114    if ( rij .LE. SSD_rup ) then
121  
122 +    if ( rij .LE. SSD_rbig ) then
123 +
124         r3 = r2*rij
125         r5 = r3*r2
126  
# Line 339 | Line 347 | contains
347   #endif
348  
349         if (do_stress) then          
342          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
350  
351 + #ifdef IS_MPI
352 +          id1 = tagRow(atom1)
353 +          id2 = tagColumn(atom2)
354 + #else
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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines