ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 378
Committed: Fri Mar 21 17:42:12 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 12728 byte(s)
Log Message:
This commit was generated by cvs2svn to compensate for changes in r377,
which included commits to RCS files with non-trunk default branches.

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.1.1.1 2003-03-21 17:42:12 mmeineke Exp $, $Date: 2003-03-21 17:42:12 $, $Name: not supported by cvs2svn $, $Revision: 1.1.1.1 $
13
14 module sticky_pair
15
16 use force_globals
17 use definitions
18 #ifdef IS_MPI
19 use mpiSimulation
20 #endif
21
22 implicit none
23
24 PRIVATE
25
26 logical, save :: sticky_initialized = .false.
27 real( kind = dp ), save :: SSD_w0
28 real( kind = dp ), save :: SSD_v0
29 real( kind = dp ), save :: SSD_rl
30 real( kind = dp ), save :: SSD_ru
31 real( kind = dp ), save :: SSD_rup
32
33 public :: check_sticky_FF
34 public :: set_sticky_params
35 public :: do_sticky_pair
36
37 contains
38
39 subroutine check_sticky_FF(status)
40 integer :: status
41 status = -1
42 if (sticky_initialized) status = 0
43 return
44 end subroutine check_sticky_FF
45
46 subroutine set_sticky_params(sticky_w0, sticky_v0)
47 real( kind = dp ), intent(in) :: sticky_w0, sticky_v0
48
49 ! we could pass all 5 parameters if we felt like it...
50
51 SSD_w0 = sticky_w0
52 SSD_v0 = sticky_v0
53 SSD_rl = 2.75_DP
54 SSD_ru = 3.35_DP
55 SSD_rup = 4.0_DP
56
57 sticky_initialized = .true.
58 return
59 end subroutine set_sticky_params
60
61 subroutine do_sticky_pair(atom1, atom2, d, rij, r2, A, pot, f, t, &
62 do_pot, do_stress)
63
64 !! This routine does only the sticky portion of the SSD potential
65 !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
66 !! The Lennard-Jones and dipolar interaction must be handled separately.
67
68 !! We assume that the rotation matrices have already been calculated
69 !! and placed in the A array.
70
71 !! i and j are pointers to the two SSD atoms
72
73 integer, intent(in) :: atom1, atom2
74 real (kind=dp), intent(inout) :: rij, r2
75 real (kind=dp), dimension(3), intent(in) :: d
76 real (kind=dp) :: pot
77 real (kind=dp), dimension(:,:) :: A
78 real (kind=dp), dimension(:,:) :: f
79 real (kind=dp), dimension(:,:) :: t
80 logical, intent(in) :: do_pot, do_stress
81
82 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
83 real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
84 real (kind=dp) :: wi, wj, w, wip, wjp, wp
85 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
86 real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
87 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
88 real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
89 real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
90 real (kind=dp) :: drdx, drdy, drdz
91 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
92 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
93 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
94 real (kind=dp) :: fxradial, fyradial, fzradial
95
96 if (.not.sticky_initialized) then
97 write(*,*) 'Sticky forces not initialized!'
98 return
99 endif
100
101 r3 = r2*rij
102 r5 = r3*r2
103
104 drdx = d(1) / rij
105 drdy = d(2) / rij
106 drdz = d(3) / rij
107
108 #ifdef IS_MPI
109 ! rotate the inter-particle separation into the two different
110 ! body-fixed coordinate systems:
111
112 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
113 yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
114 zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
115
116 ! negative sign because this is the vector from j to i:
117
118 xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
119 yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
120 zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
121 #else
122 ! rotate the inter-particle separation into the two different
123 ! body-fixed coordinate systems:
124
125 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
126 yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
127 zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
128
129 ! negative sign because this is the vector from j to i:
130
131 xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
132 yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
133 zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
134 #endif
135
136 xi2 = xi*xi
137 yi2 = yi*yi
138 zi2 = zi*zi
139
140 xj2 = xj*xj
141 yj2 = yj*yj
142 zj2 = zj*zj
143
144 call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
145
146 wi = 2.0d0*(xi2-yi2)*zi / r3
147 wj = 2.0d0*(xj2-yj2)*zj / r3
148 w = wi+wj
149
150 zif = zi/rij - 0.6d0
151 zis = zi/rij + 0.8d0
152
153 zjf = zj/rij - 0.6d0
154 zjs = zj/rij + 0.8d0
155
156 wip = zif*zif*zis*zis - SSD_w0
157 wjp = zjf*zjf*zjs*zjs - SSD_w0
158 wp = wip + wjp
159
160
161 if (do_pot) then
162 #ifdef IS_MPI
163 pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
164 pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
165 #else
166 pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
167 #endif
168 endif
169
170 dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
171 dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
172 dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
173
174 dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
175 dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
176 dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
177
178 uglyi = zif*zif*zis + zif*zis*zis
179 uglyj = zjf*zjf*zjs + zjf*zjs*zjs
180
181 dwipdx = -2.0d0*xi*zi*uglyi/r3
182 dwipdy = -2.0d0*yi*zi*uglyi/r3
183 dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
184
185 dwjpdx = -2.0d0*xj*zj*uglyj/r3
186 dwjpdy = -2.0d0*yj*zj*uglyj/r3
187 dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
188
189 dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
190 dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
191 dwiduz = - 8.0d0*xi*yi*zi/r3
192
193 dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
194 dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
195 dwjduz = - 8.0d0*xj*yj*zj/r3
196
197 dwipdux = 2.0d0*yi*uglyi/rij
198 dwipduy = -2.0d0*xi*uglyi/rij
199 dwipduz = 0.0d0
200
201 dwjpdux = 2.0d0*yj*uglyj/rij
202 dwjpduy = -2.0d0*xj*uglyj/rij
203 dwjpduz = 0.0d0
204
205 ! do the torques first since they are easy:
206 ! remember that these are still in the body fixed axes
207
208 txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
209 tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
210 tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
211
212 txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
213 tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
214 tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
215
216 ! go back to lab frame using transpose of rotation matrix:
217
218 #ifdef IS_MPI
219 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
220 a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
221 t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
222 a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
223 t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
224 a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
225
226 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
227 a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
228 t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
229 a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
230 t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
231 a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
232 #else
233 t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
234 t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
235 t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
236
237 t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
238 t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
239 t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
240 #endif
241 ! Now, on to the forces:
242
243 ! first rotate the i terms back into the lab frame:
244
245 #ifdef IS_MPI
246 fxii = a_Row(1,atom1)*(s*dwidx+sp*dwipdx) + &
247 a_Row(4,atom1)*(s*dwidy+sp*dwipdy) + &
248 a_Row(7,atom1)*(s*dwidz+sp*dwipdz)
249 fyii = a_Row(2,atom1)*(s*dwidx+sp*dwipdx) + &
250 a_Row(5,atom1)*(s*dwidy+sp*dwipdy) + &
251 a_Row(8,atom1)*(s*dwidz+sp*dwipdz)
252 fzii = a_Row(3,atom1)*(s*dwidx+sp*dwipdx) + &
253 a_Row(6,atom1)*(s*dwidy+sp*dwipdy) + &
254 a_Row(9,atom1)*(s*dwidz+sp*dwipdz)
255
256 fxjj = a_Col(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
257 a_Col(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
258 a_Col(7,atom2)*(s*dwjdz+sp*dwjpdz)
259 fyjj = a_Col(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
260 a_Col(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
261 a_Col(8,atom2)*(s*dwjdz+sp*dwjpdz)
262 fzjj = a_Col(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
263 a_Col(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
264 a_Col(9,atom2)*(s*dwjdz+sp*dwjpdz)
265 #else
266 fxii = a(1,atom1)*(s*dwidx+sp*dwipdx) + &
267 a(4,atom1)*(s*dwidy+sp*dwipdy) + &
268 a(7,atom1)*(s*dwidz+sp*dwipdz)
269 fyii = a(2,atom1)*(s*dwidx+sp*dwipdx) + &
270 a(5,atom1)*(s*dwidy+sp*dwipdy) + &
271 a(8,atom1)*(s*dwidz+sp*dwipdz)
272 fzii = a(3,atom1)*(s*dwidx+sp*dwipdx) + &
273 a(6,atom1)*(s*dwidy+sp*dwipdy) + &
274 a(9,atom1)*(s*dwidz+sp*dwipdz)
275
276 fxjj = a(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
277 a(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
278 a(7,atom2)*(s*dwjdz+sp*dwjpdz)
279 fyjj = a(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
280 a(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
281 a(8,atom2)*(s*dwjdz+sp*dwjpdz)
282 fzjj = a(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
283 a(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
284 a(9,atom2)*(s*dwjdz+sp*dwjpdz)
285 #endif
286
287 fxij = -fxii
288 fyij = -fyii
289 fzij = -fzii
290
291 fxji = -fxjj
292 fyji = -fyjj
293 fzji = -fzjj
294
295 ! now assemble these with the radial-only terms:
296
297 fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
298 fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
299 fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
300
301 #ifdef IS_MPI
302 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
303 f_Row(2,atom1) = f_Row(2,atom1) + fyradial
304 f_Row(3,atom1) = f_Row(3,atom1) + fzradial
305
306 f_Col(1,atom2) = f_Col(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - &
307 dspdr*drdx*wp + fxjj + fxij)
308 f_Col(2,atom2) = f_Col(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - &
309 dspdr*drdy*wp + fyjj + fyij)
310 f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - &
311 dspdr*drdz*wp + fzjj + fzij)
312 #else
313 f(1,atom1) = f(1,atom1) + fxradial
314 f(2,atom1) = f(2,atom1) + fyradial
315 f(3,atom1) = f(3,atom1) + fzradial
316
317 f(1,atom2) = f(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
318 fxjj + fxij)
319 f(2,atom2) = f(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
320 fyjj + fyij)
321 f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
322 fzjj + fzij)
323 #endif
324
325 if (do_stress) then
326 tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
327 tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
328 tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
329 tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
330 tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
331 tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
332 tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
333 tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
334 tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
335 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
336 endif
337
338 end subroutine do_sticky_pair
339
340 !! calculates the switching functions and their derivatives for a given
341 subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
342
343 real (kind=dp), intent(in) :: r
344 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
345
346 ! distances must be in angstroms
347
348 if (r.lt.SSD_rl) then
349 s = 1.0d0
350 sp = 1.0d0
351 dsdr = 0.0d0
352 dspdr = 0.0d0
353 elseif (r.gt.SSD_rup) then
354 s = 0.0d0
355 sp = 0.0d0
356 dsdr = 0.0d0
357 dspdr = 0.0d0
358 else
359 sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_rup-r)**2) / &
360 ((SSD_rup - SSD_rl)**3)
361 dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rl)/((SSD_rup - SSD_rl)**3)
362
363 if (r.gt.SSD_ru) then
364 s = 0.0d0
365 dsdr = 0.0d0
366 else
367 s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
368 ((SSD_ru - SSD_rl)**3)
369 dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
370 endif
371 endif
372
373 return
374 end subroutine calc_sw_fnc
375 end module sticky_pair