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 483 by gezelter, Wed Apr 9 04:06:43 2003 UTC vs.
Revision 611 by gezelter, Tue Jul 15 17:10:50 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.8 2003-04-09 04:06:43 gezelter Exp $, $Date: 2003-04-09 04:06:43 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.10 2003-07-15 17:10:50 gezelter Exp $, $Date: 2003-07-15 17:10:50 $, $Name: not supported by cvs2svn $, $Revision: 1.10 $
13  
14   module sticky_pair
15  
# Line 61 | Line 61 | contains
61  
62    subroutine do_sticky_pair(atom1, atom2, d, rij, r2, A, pot, f, t, &
63         do_pot, do_stress)
64 <    
64 >
65      !! This routine does only the sticky portion of the SSD potential
66      !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
67      !! The Lennard-Jones and dipolar interaction must be handled separately.
68 <    
68 >
69      !! We assume that the rotation matrices have already been calculated
70      !! and placed in the A array.
71 <    
71 >
72      !! i and j are pointers to the two SSD atoms
73  
74      integer, intent(in) :: atom1, atom2
# Line 94 | Line 94 | contains
94      real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji      
95      real (kind=dp) :: fxradial, fyradial, fzradial
96      real (kind=dp) :: rijtest, rjitest
97 <      
97 >    real (kind=dp) :: radcomxi, radcomyi, radcomzi
98 >    real (kind=dp) :: radcomxj, radcomyj, radcomzj
99 >
100 >
101      if (.not.sticky_initialized) then
102         write(*,*) 'Sticky forces not initialized!'
103         return
104      endif
105  
106 <    r3 = r2*rij
104 <    r5 = r3*r2
105 <    
106 <    drdx = d(1) / rij
107 <    drdy = d(2) / rij
108 <    drdz = d(3) / rij
109 <    
110 < #ifdef IS_MPI
111 <    ! rotate the inter-particle separation into the two different
112 <    ! body-fixed coordinate systems:
113 <    
114 <    xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
115 <    yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
116 <    zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
117 <    
118 <    ! negative sign because this is the vector from j to i:
119 <    
120 <    xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
121 <    yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
122 <    zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
123 < #else
124 <    ! rotate the inter-particle separation into the two different
125 <    ! body-fixed coordinate systems:
126 <    
127 <    xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
128 <    yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
129 <    zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
130 <    
131 <    ! negative sign because this is the vector from j to i:
132 <    
133 <    xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
134 <    yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
135 <    zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
136 < #endif
137 <    
138 <    xi2 = xi*xi
139 <    yi2 = yi*yi
140 <    zi2 = zi*zi
141 <    
142 <    xj2 = xj*xj
143 <    yj2 = yj*yj
144 <    zj2 = zj*zj
145 <  
146 <    call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
147 <    
148 <    wi = 2.0d0*(xi2-yi2)*zi / r3
149 <    wj = 2.0d0*(xj2-yj2)*zj / r3
150 <    w = wi+wj
151 <    
152 <    zif = zi/rij - 0.6d0
153 <    zis = zi/rij + 0.8d0
154 <    
155 <    zjf = zj/rij - 0.6d0
156 <    zjs = zj/rij + 0.8d0
157 <    
158 <    wip = zif*zif*zis*zis - SSD_w0
159 <    wjp = zjf*zjf*zjs*zjs - SSD_w0
160 <    wp = wip + wjp
106 >    if ( rij .LE. SSD_rup ) then
107  
108 <    if (do_pot) then
109 < #ifdef IS_MPI
110 <       pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
111 <       pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
112 < #else
113 <       pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
114 < #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
173 <    dwidz =   2.0d0*(xi2-yi2)/r3  - 6.0d0*zi2*(xi2-yi2)/r5
174 <    
175 <    dwjdx =   4.0d0*xj*zj/r3  - 6.0d0*xj*zj*(xj2-yj2)/r5
176 <    dwjdy = - 4.0d0*yj*zj/r3  - 6.0d0*yj*zj*(xj2-yj2)/r5
177 <    dwjdz =   2.0d0*(xj2-yj2)/r3  - 6.0d0*zj2*(xj2-yj2)/r5
178 <    
179 <    uglyi = zif*zif*zis + zif*zis*zis
180 <    uglyj = zjf*zjf*zjs + zjf*zjs*zjs
181 <    
182 <    dwipdx = -2.0d0*xi*zi*uglyi/r3
183 <    dwipdy = -2.0d0*yi*zi*uglyi/r3
184 <    dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
185 <    
186 <    dwjpdx = -2.0d0*xj*zj*uglyj/r3
187 <    dwjpdy = -2.0d0*yj*zj*uglyj/r3
188 <    dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
189 <    
190 <    dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
191 <    dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
192 <    dwiduz = - 8.0d0*xi*yi*zi/r3
193 <    
194 <    dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
195 <    dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
196 <    dwjduz = - 8.0d0*xj*yj*zj/r3
197 <    
198 <    dwipdux =  2.0d0*yi*uglyi/rij
199 <    dwipduy = -2.0d0*xi*uglyi/rij
200 <    dwipduz =  0.0d0
201 <    
202 <    dwjpdux =  2.0d0*yj*uglyj/rij
203 <    dwjpduy = -2.0d0*xj*uglyj/rij
204 <    dwjpduz =  0.0d0
205 <    
206 <    ! do the torques first since they are easy:
207 <    ! remember that these are still in the body fixed axes
208 <    
209 <    txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
210 <    tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
211 <    tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
212 <    
213 <    txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
214 <    tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
215 <    tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
216 <    
217 <    ! go back to lab frame using transpose of rotation matrix:
218 <  
108 >       r3 = r2*rij
109 >       r5 = r3*r2
110 >
111 >       drdx = d(1) / rij
112 >       drdy = d(2) / rij
113 >       drdz = d(3) / rij
114 >
115   #ifdef IS_MPI
116 <    t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
117 <         a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
222 <    t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
223 <         a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
224 <    t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
225 <         a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
226 <    
227 <    t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
228 <         a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
229 <    t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
230 <         a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
231 <    t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
232 <         a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
233 < #else
234 <    t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
235 <    t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
236 <    t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
237 <    
238 <    t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
239 <    t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
240 <    t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
241 < #endif    
242 <    ! Now, on to the forces:
243 <    
244 <    ! first rotate the i terms back into the lab frame:
116 >       ! rotate the inter-particle separation into the two different
117 >       ! body-fixed coordinate systems:
118  
119 < #ifdef IS_MPI    
120 <    fxii = a_Row(1,atom1)*(s*dwidx+sp*dwipdx) + &
121 <         a_Row(4,atom1)*(s*dwidy+sp*dwipdy) + &
249 <         a_Row(7,atom1)*(s*dwidz+sp*dwipdz)
250 <    fyii = a_Row(2,atom1)*(s*dwidx+sp*dwipdx) + &
251 <         a_Row(5,atom1)*(s*dwidy+sp*dwipdy) + &
252 <         a_Row(8,atom1)*(s*dwidz+sp*dwipdz)
253 <    fzii = a_Row(3,atom1)*(s*dwidx+sp*dwipdx) + &
254 <         a_Row(6,atom1)*(s*dwidy+sp*dwipdy) + &
255 <         a_Row(9,atom1)*(s*dwidz+sp*dwipdz)
119 >       xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
120 >       yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
121 >       zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
122  
123 <    fxjj = a_Col(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
124 <         a_Col(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
125 <         a_Col(7,atom2)*(s*dwjdz+sp*dwjpdz)
126 <    fyjj = a_Col(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
127 <         a_Col(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
262 <         a_Col(8,atom2)*(s*dwjdz+sp*dwjpdz)
263 <    fzjj = a_Col(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
264 <         a_Col(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
265 <         a_Col(9,atom2)*(s*dwjdz+sp*dwjpdz)
123 >       ! negative sign because this is the vector from j to i:
124 >
125 >       xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
126 >       yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
127 >       zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
128   #else
129 <    fxii = a(1,atom1)*(s*dwidx+sp*dwipdx) + &
130 <         a(4,atom1)*(s*dwidy+sp*dwipdy) + &
269 <         a(7,atom1)*(s*dwidz+sp*dwipdz)
270 <    fyii = a(2,atom1)*(s*dwidx+sp*dwipdx) + &
271 <         a(5,atom1)*(s*dwidy+sp*dwipdy) + &
272 <         a(8,atom1)*(s*dwidz+sp*dwipdz)
273 <    fzii = a(3,atom1)*(s*dwidx+sp*dwipdx) + &
274 <         a(6,atom1)*(s*dwidy+sp*dwipdy) + &
275 <         a(9,atom1)*(s*dwidz+sp*dwipdz)
129 >       ! rotate the inter-particle separation into the two different
130 >       ! body-fixed coordinate systems:
131  
132 <    fxjj = a(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
133 <         a(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
134 <         a(7,atom2)*(s*dwjdz+sp*dwjpdz)
135 <    fyjj = a(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
136 <         a(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
137 <         a(8,atom2)*(s*dwjdz+sp*dwjpdz)
138 <    fzjj = a(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
139 <         a(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
140 <         a(9,atom2)*(s*dwjdz+sp*dwjpdz)
132 >       xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
133 >       yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
134 >       zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
135 >
136 >       ! negative sign because this is the vector from j to i:
137 >
138 >       xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
139 >       yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
140 >       zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
141   #endif
287    
288    fxij = -fxii
289    fyij = -fyii
290    fzij = -fzii
291        
292    fxji = -fxjj
293    fyji = -fyjj
294    fzji = -fzjj
295    
296    ! now assemble these with the radial-only terms:
142  
143 <    fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
144 <    fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
145 <    fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
146 <      
143 >       xi2 = xi*xi
144 >       yi2 = yi*yi
145 >       zi2 = zi*zi
146 >
147 >       xj2 = xj*xj
148 >       yj2 = yj*yj
149 >       zj2 = zj*zj
150 >
151 >       call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
152 >
153 >       wi = 2.0d0*(xi2-yi2)*zi / r3
154 >       wj = 2.0d0*(xj2-yj2)*zj / r3
155 >       w = wi+wj
156 >
157 >       zif = zi/rij - 0.6d0
158 >       zis = zi/rij + 0.8d0
159 >
160 >       zjf = zj/rij - 0.6d0
161 >       zjs = zj/rij + 0.8d0
162 >
163 >       wip = zif*zif*zis*zis - SSD_w0
164 >       wjp = zjf*zjf*zjs*zjs - SSD_w0
165 >       wp = wip + wjp
166 >
167 >       if (do_pot) then
168 > #ifdef IS_MPI
169 >          pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
170 >          pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
171 > #else
172 >          pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
173 > #endif  
174 >       endif
175 >
176 >       dwidx =   4.0d0*xi*zi/r3  - 6.0d0*xi*zi*(xi2-yi2)/r5
177 >       dwidy = - 4.0d0*yi*zi/r3  - 6.0d0*yi*zi*(xi2-yi2)/r5
178 >       dwidz =   2.0d0*(xi2-yi2)/r3  - 6.0d0*zi2*(xi2-yi2)/r5
179 >
180 >       dwjdx =   4.0d0*xj*zj/r3  - 6.0d0*xj*zj*(xj2-yj2)/r5
181 >       dwjdy = - 4.0d0*yj*zj/r3  - 6.0d0*yj*zj*(xj2-yj2)/r5
182 >       dwjdz =   2.0d0*(xj2-yj2)/r3  - 6.0d0*zj2*(xj2-yj2)/r5
183 >
184 >       uglyi = zif*zif*zis + zif*zis*zis
185 >       uglyj = zjf*zjf*zjs + zjf*zjs*zjs
186 >
187 >       dwipdx = -2.0d0*xi*zi*uglyi/r3
188 >       dwipdy = -2.0d0*yi*zi*uglyi/r3
189 >       dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
190 >
191 >       dwjpdx = -2.0d0*xj*zj*uglyj/r3
192 >       dwjpdy = -2.0d0*yj*zj*uglyj/r3
193 >       dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
194 >
195 >       dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
196 >       dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
197 >       dwiduz = - 8.0d0*xi*yi*zi/r3
198 >
199 >       dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
200 >       dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
201 >       dwjduz = - 8.0d0*xj*yj*zj/r3
202 >
203 >       dwipdux =  2.0d0*yi*uglyi/rij
204 >       dwipduy = -2.0d0*xi*uglyi/rij
205 >       dwipduz =  0.0d0
206 >
207 >       dwjpdux =  2.0d0*yj*uglyj/rij
208 >       dwjpduy = -2.0d0*xj*uglyj/rij
209 >       dwjpduz =  0.0d0
210 >
211 >       ! do the torques first since they are easy:
212 >       ! remember that these are still in the body fixed axes
213 >
214 >       txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
215 >       tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
216 >       tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
217 >
218 >       txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
219 >       tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
220 >       tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
221 >
222 >       ! go back to lab frame using transpose of rotation matrix:
223 >
224   #ifdef IS_MPI
225 <    f_Row(1,atom1) = f_Row(1,atom1) + fxradial
226 <    f_Row(2,atom1) = f_Row(2,atom1) + fyradial
227 <    f_Row(3,atom1) = f_Row(3,atom1) + fzradial
228 <    
229 <    f_Col(1,atom2) = f_Col(1,atom2) - fxradial
230 <    f_Col(2,atom2) = f_Col(2,atom2) - fyradial
231 <    f_Col(3,atom2) = f_Col(3,atom2) - fzradial
225 >       t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
226 >            a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
227 >       t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
228 >            a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
229 >       t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
230 >            a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
231 >
232 >       t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
233 >            a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
234 >       t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
235 >            a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
236 >       t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
237 >            a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
238   #else
239 <    f(1,atom1) = f(1,atom1) + fxradial
240 <    f(2,atom1) = f(2,atom1) + fyradial
241 <    f(3,atom1) = f(3,atom1) + fzradial
242 <    
243 <    f(1,atom2) = f(1,atom2) - fxradial
244 <    f(2,atom2) = f(2,atom2) - fyradial
245 <    f(3,atom2) = f(3,atom2) - fzradial
239 >       t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
240 >       t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
241 >       t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
242 >
243 >       t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
244 >       t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
245 >       t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
246 > #endif    
247 >       ! Now, on to the forces:
248 >
249 >       ! first rotate the i terms back into the lab frame:
250 >
251 >       radcomxi = s*dwidx+sp*dwipdx
252 >       radcomyi = s*dwidy+sp*dwipdy
253 >       radcomzi = s*dwidz+sp*dwipdz
254 >
255 >       radcomxj = s*dwjdx+sp*dwjpdx
256 >       radcomyj = s*dwjdy+sp*dwjpdy
257 >       radcomzj = s*dwjdz+sp*dwjpdz
258 >
259 > #ifdef IS_MPI    
260 >       fxii = a_Row(1,atom1)*(radcomxi) + &
261 >            a_Row(4,atom1)*(radcomyi) + &
262 >            a_Row(7,atom1)*(radcomzi)
263 >       fyii = a_Row(2,atom1)*(radcomxi) + &
264 >            a_Row(5,atom1)*(radcomyi) + &
265 >            a_Row(8,atom1)*(radcomzi)
266 >       fzii = a_Row(3,atom1)*(radcomxi) + &
267 >            a_Row(6,atom1)*(radcomyi) + &
268 >            a_Row(9,atom1)*(radcomzi)
269 >
270 >       fxjj = a_Col(1,atom2)*(radcomxj) + &
271 >            a_Col(4,atom2)*(radcomyj) + &
272 >            a_Col(7,atom2)*(radcomzj)
273 >       fyjj = a_Col(2,atom2)*(radcomxj) + &
274 >            a_Col(5,atom2)*(radcomyj) + &
275 >            a_Col(8,atom2)*(radcomzj)
276 >       fzjj = a_Col(3,atom2)*(radcomxj)+ &
277 >            a_Col(6,atom2)*(radcomyj) + &
278 >            a_Col(9,atom2)*(radcomzj)
279 > #else
280 >       fxii = a(1,atom1)*(radcomxi) + &
281 >            a(4,atom1)*(radcomyi) + &
282 >            a(7,atom1)*(radcomzi)
283 >       fyii = a(2,atom1)*(radcomxi) + &
284 >            a(5,atom1)*(radcomyi) + &
285 >            a(8,atom1)*(radcomzi)
286 >       fzii = a(3,atom1)*(radcomxi) + &
287 >            a(6,atom1)*(radcomyi) + &
288 >            a(9,atom1)*(radcomzi)
289 >
290 >       fxjj = a(1,atom2)*(radcomxj) + &
291 >            a(4,atom2)*(radcomyj) + &
292 >            a(7,atom2)*(radcomzj)
293 >       fyjj = a(2,atom2)*(radcomxj) + &
294 >            a(5,atom2)*(radcomyj) + &
295 >            a(8,atom2)*(radcomzj)
296 >       fzjj = a(3,atom2)*(radcomxj)+ &
297 >            a(6,atom2)*(radcomyj) + &
298 >            a(9,atom2)*(radcomzj)
299   #endif
300 <    
301 <    if (do_stress) then          
302 <       if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
303 <          tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
304 <          tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
305 <          tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
306 <          tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
307 <          tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
308 <          tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
309 <          tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
310 <          tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
311 <          tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
312 <          virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
300 >
301 >       fxij = -fxii
302 >       fyij = -fyii
303 >       fzij = -fzii
304 >
305 >       fxji = -fxjj
306 >       fyji = -fyjj
307 >       fzji = -fzjj
308 >
309 >       ! now assemble these with the radial-only terms:
310 >
311 >       fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
312 >       fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
313 >       fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
314 >
315 > #ifdef IS_MPI
316 >       f_Row(1,atom1) = f_Row(1,atom1) + fxradial
317 >       f_Row(2,atom1) = f_Row(2,atom1) + fyradial
318 >       f_Row(3,atom1) = f_Row(3,atom1) + fzradial
319 >
320 >       f_Col(1,atom2) = f_Col(1,atom2) - fxradial
321 >       f_Col(2,atom2) = f_Col(2,atom2) - fyradial
322 >       f_Col(3,atom2) = f_Col(3,atom2) - fzradial
323 > #else
324 >       f(1,atom1) = f(1,atom1) + fxradial
325 >       f(2,atom1) = f(2,atom1) + fyradial
326 >       f(3,atom1) = f(3,atom1) + fzradial
327 >
328 >       f(1,atom2) = f(1,atom2) - fxradial
329 >       f(2,atom2) = f(2,atom2) - fyradial
330 >       f(3,atom2) = f(3,atom2) - fzradial
331 > #endif
332 >
333 >       if (do_stress) then          
334 >          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
335 >
336 >             ! because the d vector is the rj - ri vector, and
337 >             ! because fxradial, fyradial, and fzradial are the
338 >             ! (positive) force on atom i (negative on atom j) we need
339 >             ! a negative sign here:
340 >
341 >             tau_Temp(1) = tau_Temp(1) - d(1) * fxradial
342 >             tau_Temp(2) = tau_Temp(2) - d(1) * fyradial
343 >             tau_Temp(3) = tau_Temp(3) - d(1) * fzradial
344 >             tau_Temp(4) = tau_Temp(4) - d(2) * fxradial
345 >             tau_Temp(5) = tau_Temp(5) - d(2) * fyradial
346 >             tau_Temp(6) = tau_Temp(6) - d(2) * fzradial
347 >             tau_Temp(7) = tau_Temp(7) - d(3) * fxradial
348 >             tau_Temp(8) = tau_Temp(8) - d(3) * fyradial
349 >             tau_Temp(9) = tau_Temp(9) - d(3) * fzradial
350 >
351 >             virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
352 >          endif
353         endif
354      endif
355 <  
355 >
356    end subroutine do_sticky_pair
357  
358    !! calculates the switching functions and their derivatives for a given

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines