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 322 by gezelter, Wed Mar 12 15:17:52 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.2 2003-03-12 15:17:52 gezelter Exp $, $Date: 2003-03-12 15:17:52 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $
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 25 | Line 25 | module sticky_pair
25    PRIVATE
26  
27    logical, save :: sticky_initialized = .false.
28 <  real( kind = dp ), save :: SSD_w0
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
33  
34 <  public :: initialize_sticky
34 >  public :: check_sticky_FF
35 >  public :: set_sticky_params
36    public :: do_sticky_pair
37  
38   contains
39  
40 <  subroutine initialize_sticky(sticky_w0, sticky_v0)
40 >  subroutine check_sticky_FF(status)
41 >    integer :: status
42 >    status = -1
43 >    if (sticky_initialized) status = 0
44 >    return
45 >  end subroutine check_sticky_FF
46 >
47 >  subroutine set_sticky_params(sticky_w0, sticky_v0)
48      real( kind = dp ), intent(in) :: sticky_w0, sticky_v0
49      
50      ! we could pass all 5 parameters if we felt like it...
51 <
51 >    
52      SSD_w0 = sticky_w0
53      SSD_v0 = sticky_v0
54      SSD_rl = 2.75_DP
55      SSD_ru = 3.35_DP
56      SSD_rup = 4.0_DP
57 <
57 >  
58      sticky_initialized = .true.
59      return
60 <  end subroutine initialize_sticky
60 >  end subroutine set_sticky_params
61  
62 <  subroutine do_sticky_pair(atom1, atom2, d, rij, 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 63 | Line 72 | contains
72      !! i and j are pointers to the two SSD atoms
73  
74      integer, intent(in) :: atom1, atom2
75 <    real (kind=dp), intent(inout) :: rij    
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(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) :: r2, r3, r5, r6, s, sp, dsdr, dspdr
84 >    real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
85      real (kind=dp) :: wi, wj, w, wip, wjp, wp
86      real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
87      real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
# Line 88 | Line 99 | contains
99         return
100      endif
101  
91    r2 = rij*rij
102      r3 = r2*rij
103      r5 = r3*r2
104      
# Line 148 | Line 158 | contains
158      wjp = zjf*zjf*zjs*zjs - SSD_w0
159      wp = wip + wjp
160  
161 < #ifdef IS_MPI
162 <    pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
163 <    pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
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)
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 310 | 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