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