ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 443
Committed: Wed Apr 2 22:19:03 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 12871 byte(s)
Log Message:
dipoles mostly work, but there is a memory leak somewhere.

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.3 2003-04-02 22:19:03 mmeineke Exp $, $Date: 2003-04-02 22:19:03 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
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 real (kind=dp) :: rijtest, rjitest
96
97 if (.not.sticky_initialized) then
98 write(*,*) 'Sticky forces not initialized!'
99 return
100 endif
101
102 r3 = r2*rij
103 r5 = r3*r2
104
105 drdx = d(1) / rij
106 drdy = d(2) / rij
107 drdz = d(3) / rij
108
109 #ifdef IS_MPI
110 ! rotate the inter-particle separation into the two different
111 ! body-fixed coordinate systems:
112
113 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
114 yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
115 zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
116
117 ! negative sign because this is the vector from j to i:
118
119 xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
120 yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
121 zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
122 #else
123 ! rotate the inter-particle separation into the two different
124 ! body-fixed coordinate systems:
125
126 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
127 yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
128 zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
129
130 ! negative sign because this is the vector from j to i:
131
132 xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
133 yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
134 zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
135 #endif
136
137 xi2 = xi*xi
138 yi2 = yi*yi
139 zi2 = zi*zi
140
141 xj2 = xj*xj
142 yj2 = yj*yj
143 zj2 = zj*zj
144
145 call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
146
147 wi = 2.0d0*(xi2-yi2)*zi / r3
148 wj = 2.0d0*(xj2-yj2)*zj / r3
149 w = wi+wj
150
151 zif = zi/rij - 0.6d0
152 zis = zi/rij + 0.8d0
153
154 zjf = zj/rij - 0.6d0
155 zjs = zj/rij + 0.8d0
156
157 wip = zif*zif*zis*zis - SSD_w0
158 wjp = zjf*zjf*zjs*zjs - SSD_w0
159 wp = wip + wjp
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 write(*,'(a4,3es12.3)') 't1= ', t(1:3,atom1)
244 write(*,'(a4,3es12.3)') 't2= ', t(1:3,atom2)
245 write(*,*)
246
247 ! first rotate the i terms back into the lab frame:
248
249 #ifdef IS_MPI
250 fxii = a_Row(1,atom1)*(s*dwidx+sp*dwipdx) + &
251 a_Row(4,atom1)*(s*dwidy+sp*dwipdy) + &
252 a_Row(7,atom1)*(s*dwidz+sp*dwipdz)
253 fyii = a_Row(2,atom1)*(s*dwidx+sp*dwipdx) + &
254 a_Row(5,atom1)*(s*dwidy+sp*dwipdy) + &
255 a_Row(8,atom1)*(s*dwidz+sp*dwipdz)
256 fzii = a_Row(3,atom1)*(s*dwidx+sp*dwipdx) + &
257 a_Row(6,atom1)*(s*dwidy+sp*dwipdy) + &
258 a_Row(9,atom1)*(s*dwidz+sp*dwipdz)
259
260 fxjj = a_Col(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
261 a_Col(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
262 a_Col(7,atom2)*(s*dwjdz+sp*dwjpdz)
263 fyjj = a_Col(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
264 a_Col(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
265 a_Col(8,atom2)*(s*dwjdz+sp*dwjpdz)
266 fzjj = a_Col(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
267 a_Col(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
268 a_Col(9,atom2)*(s*dwjdz+sp*dwjpdz)
269 #else
270 fxii = a(1,atom1)*(s*dwidx+sp*dwipdx) + &
271 a(4,atom1)*(s*dwidy+sp*dwipdy) + &
272 a(7,atom1)*(s*dwidz+sp*dwipdz)
273 fyii = a(2,atom1)*(s*dwidx+sp*dwipdx) + &
274 a(5,atom1)*(s*dwidy+sp*dwipdy) + &
275 a(8,atom1)*(s*dwidz+sp*dwipdz)
276 fzii = a(3,atom1)*(s*dwidx+sp*dwipdx) + &
277 a(6,atom1)*(s*dwidy+sp*dwipdy) + &
278 a(9,atom1)*(s*dwidz+sp*dwipdz)
279
280 fxjj = a(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
281 a(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
282 a(7,atom2)*(s*dwjdz+sp*dwjpdz)
283 fyjj = a(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
284 a(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
285 a(8,atom2)*(s*dwjdz+sp*dwjpdz)
286 fzjj = a(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
287 a(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
288 a(9,atom2)*(s*dwjdz+sp*dwjpdz)
289 #endif
290
291 fxij = -fxii
292 fyij = -fyii
293 fzij = -fzii
294
295 fxji = -fxjj
296 fyji = -fyjj
297 fzji = -fzjj
298
299 ! now assemble these with the radial-only terms:
300
301 fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
302 fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
303 fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
304
305 #ifdef IS_MPI
306 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
307 f_Row(2,atom1) = f_Row(2,atom1) + fyradial
308 f_Row(3,atom1) = f_Row(3,atom1) + fzradial
309
310 f_Col(1,atom2) = f_Col(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - &
311 dspdr*drdx*wp + fxjj + fxij)
312 f_Col(2,atom2) = f_Col(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - &
313 dspdr*drdy*wp + fyjj + fyij)
314 f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - &
315 dspdr*drdz*wp + fzjj + fzij)
316 #else
317 f(1,atom1) = f(1,atom1) + fxradial
318 f(2,atom1) = f(2,atom1) + fyradial
319 f(3,atom1) = f(3,atom1) + fzradial
320
321 f(1,atom2) = f(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
322 fxjj + fxij)
323 f(2,atom2) = f(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
324 fyjj + fyij)
325 f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
326 fzjj + fzij)
327 #endif
328
329 if (do_stress) then
330 tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
331 tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
332 tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
333 tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
334 tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
335 tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
336 tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
337 tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
338 tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
339 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
340 endif
341
342 end subroutine do_sticky_pair
343
344 !! calculates the switching functions and their derivatives for a given
345 subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
346
347 real (kind=dp), intent(in) :: r
348 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
349
350 ! distances must be in angstroms
351
352 if (r.lt.SSD_rl) then
353 s = 1.0d0
354 sp = 1.0d0
355 dsdr = 0.0d0
356 dspdr = 0.0d0
357 elseif (r.gt.SSD_rup) then
358 s = 0.0d0
359 sp = 0.0d0
360 dsdr = 0.0d0
361 dspdr = 0.0d0
362 else
363 sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_rup-r)**2) / &
364 ((SSD_rup - SSD_rl)**3)
365 dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rl)/((SSD_rup - SSD_rl)**3)
366
367 if (r.gt.SSD_ru) then
368 s = 0.0d0
369 dsdr = 0.0d0
370 else
371 s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
372 ((SSD_ru - SSD_rl)**3)
373 dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
374 endif
375 endif
376
377 return
378 end subroutine calc_sw_fnc
379 end module sticky_pair