ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 1160
Committed: Tue May 11 21:31:15 2004 UTC (20 years, 1 month ago) by gezelter
File size: 14169 byte(s)
Log Message:
Fortran-side changes for group-based cutoffs

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