ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 473
Committed: Mon Apr 7 21:20:38 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 12488 byte(s)
Log Message:
fixed a sign error in the radial portion of SSD.

File Contents

# Content
1 !! This Module Calculates forces due to SSD potential and VDW interactions
2 !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
3
4 !! This module contains the Public procedures:
5
6
7 !! Corresponds to the force field defined in ssd_FF.cpp
8 !! @author Charles F. Vardeman II
9 !! @author Matthew Meineke
10 !! @author Christopher Fennel
11 !! @author J. Daniel Gezelter
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
22
23 implicit none
24
25 PRIVATE
26
27 logical, save :: sticky_initialized = .false.
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
36 public :: do_sticky_pair
37
38 contains
39
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
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
58 sticky_initialized = .true.
59 return
60 end subroutine set_sticky_params
61
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)].
67 !! The Lennard-Jones and dipolar interaction must be handled separately.
68
69 !! We assume that the rotation matrices have already been calculated
70 !! and placed in the A array.
71
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, 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) :: 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
88 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
89 real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
90 real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
91 real (kind=dp) :: drdx, drdy, drdz
92 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
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!'
100 return
101 endif
102
103 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
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
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
219 #ifdef IS_MPI
220 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
221 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:
245
246 #ifdef IS_MPI
247 fxii = a_Row(1,atom1)*(s*dwidx+sp*dwipdx) + &
248 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)
256
257 fxjj = a_Col(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
258 a_Col(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
259 a_Col(7,atom2)*(s*dwjdz+sp*dwjpdz)
260 fyjj = a_Col(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
261 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)
266 #else
267 fxii = a(1,atom1)*(s*dwidx+sp*dwipdx) + &
268 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)
276
277 fxjj = a(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
278 a(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
279 a(7,atom2)*(s*dwjdz+sp*dwjpdz)
280 fyjj = a(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
281 a(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
282 a(8,atom2)*(s*dwjdz+sp*dwjpdz)
283 fzjj = a(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
284 a(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
285 a(9,atom2)*(s*dwjdz+sp*dwjpdz)
286 #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:
297
298 fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
299 fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
300 fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
301
302 #ifdef IS_MPI
303 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
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) - 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) + 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))
331 endif
332
333 end subroutine do_sticky_pair
334
335 !! calculates the switching functions and their derivatives for a given
336 subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
337
338 real (kind=dp), intent(in) :: r
339 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
340
341 ! distances must be in angstroms
342
343 if (r.lt.SSD_rl) then
344 s = 1.0d0
345 sp = 1.0d0
346 dsdr = 0.0d0
347 dspdr = 0.0d0
348 elseif (r.gt.SSD_rup) then
349 s = 0.0d0
350 sp = 0.0d0
351 dsdr = 0.0d0
352 dspdr = 0.0d0
353 else
354 sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_rup-r)**2) / &
355 ((SSD_rup - SSD_rl)**3)
356 dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rl)/((SSD_rup - SSD_rl)**3)
357
358 if (r.gt.SSD_ru) then
359 s = 0.0d0
360 dsdr = 0.0d0
361 else
362 s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
363 ((SSD_ru - SSD_rl)**3)
364 dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
365 endif
366 endif
367
368 return
369 end subroutine calc_sw_fnc
370 end module sticky_pair