1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
|
42 |
!! This Module Calculates forces due to SSD potential and VDW interactions |
43 |
!! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)]. |
44 |
|
45 |
!! This module contains the Public procedures: |
46 |
|
47 |
|
48 |
!! Corresponds to the force field defined in ssd_FF.cpp |
49 |
!! @author Charles F. Vardeman II |
50 |
!! @author Matthew Meineke |
51 |
!! @author Christopher Fennell |
52 |
!! @author J. Daniel Gezelter |
53 |
!! @version $Id$, $Date$, $Name: not supported by cvs2svn $, $Revision$ |
54 |
|
55 |
module sticky |
56 |
|
57 |
use force_globals |
58 |
use definitions |
59 |
use atype_module |
60 |
use vector_class |
61 |
use simulation |
62 |
use status |
63 |
use interpolation |
64 |
implicit none |
65 |
|
66 |
PRIVATE |
67 |
#define __FORTRAN90 |
68 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
69 |
|
70 |
public :: newStickyType |
71 |
public :: do_sticky_pair |
72 |
public :: destroyStickyTypes |
73 |
public :: do_sticky_power_pair |
74 |
public :: getStickyCut |
75 |
public :: getStickyPowerCut |
76 |
|
77 |
type :: StickyList |
78 |
integer :: c_ident |
79 |
real( kind = dp ) :: w0 = 0.0_dp |
80 |
real( kind = dp ) :: v0 = 0.0_dp |
81 |
real( kind = dp ) :: v0p = 0.0_dp |
82 |
real( kind = dp ) :: rl = 0.0_dp |
83 |
real( kind = dp ) :: ru = 0.0_dp |
84 |
real( kind = dp ) :: rlp = 0.0_dp |
85 |
real( kind = dp ) :: rup = 0.0_dp |
86 |
real( kind = dp ) :: rbig = 0.0_dp |
87 |
type(cubicSpline) :: stickySpline |
88 |
type(cubicSpline) :: stickySplineP |
89 |
end type StickyList |
90 |
|
91 |
type(StickyList), dimension(:),allocatable :: StickyMap |
92 |
logical, save :: hasStickyMap = .false. |
93 |
|
94 |
contains |
95 |
|
96 |
subroutine newStickyType(c_ident, w0, v0, v0p, rl, ru, rlp, rup, isError) |
97 |
|
98 |
integer, intent(in) :: c_ident |
99 |
integer, intent(inout) :: isError |
100 |
real( kind = dp ), intent(in) :: w0, v0, v0p |
101 |
real( kind = dp ), intent(in) :: rl, ru |
102 |
real( kind = dp ), intent(in) :: rlp, rup |
103 |
real( kind = dp ), dimension(2) :: rCubVals, sCubVals, rpCubVals, spCubVals |
104 |
integer :: nATypes, myATID |
105 |
|
106 |
|
107 |
isError = 0 |
108 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
109 |
|
110 |
!! Be simple-minded and assume that we need a StickyMap that |
111 |
!! is the same size as the total number of atom types |
112 |
|
113 |
if (.not.allocated(StickyMap)) then |
114 |
|
115 |
nAtypes = getSize(atypes) |
116 |
|
117 |
if (nAtypes == 0) then |
118 |
isError = -1 |
119 |
return |
120 |
end if |
121 |
|
122 |
if (.not. allocated(StickyMap)) then |
123 |
allocate(StickyMap(nAtypes)) |
124 |
endif |
125 |
|
126 |
end if |
127 |
|
128 |
if (myATID .gt. size(StickyMap)) then |
129 |
isError = -1 |
130 |
return |
131 |
endif |
132 |
|
133 |
! set the values for StickyMap for this atom type: |
134 |
|
135 |
StickyMap(myATID)%c_ident = c_ident |
136 |
|
137 |
! we could pass all 5 parameters if we felt like it... |
138 |
|
139 |
StickyMap(myATID)%w0 = w0 |
140 |
StickyMap(myATID)%v0 = v0 |
141 |
StickyMap(myATID)%v0p = v0p |
142 |
StickyMap(myATID)%rl = rl |
143 |
StickyMap(myATID)%ru = ru |
144 |
StickyMap(myATID)%rlp = rlp |
145 |
StickyMap(myATID)%rup = rup |
146 |
|
147 |
if (StickyMap(myATID)%ru .gt. StickyMap(myATID)%rup) then |
148 |
StickyMap(myATID)%rbig = StickyMap(myATID)%ru |
149 |
else |
150 |
StickyMap(myATID)%rbig = StickyMap(myATID)%rup |
151 |
endif |
152 |
|
153 |
! build the 2 cubic splines for the sticky switching functions |
154 |
|
155 |
rCubVals(1) = rl |
156 |
rCubVals(2) = ru |
157 |
sCubVals(1) = 1.0_dp |
158 |
sCubVals(2) = 0.0_dp |
159 |
call newSpline(StickyMap(myATID)%stickySpline, rCubVals, sCubVals, .true.) |
160 |
rpCubVals(1) = rlp |
161 |
rpCubVals(2) = rup |
162 |
spCubVals(1) = 1.0_dp |
163 |
spCubVals(2) = 0.0_dp |
164 |
call newSpline(StickyMap(myATID)%stickySplineP,rpCubVals,spCubVals,.true.) |
165 |
|
166 |
hasStickyMap = .true. |
167 |
|
168 |
return |
169 |
end subroutine newStickyType |
170 |
|
171 |
function getStickyCut(atomID) result(cutValue) |
172 |
integer, intent(in) :: atomID |
173 |
real(kind=dp) :: cutValue |
174 |
|
175 |
cutValue = StickyMap(atomID)%rbig |
176 |
end function getStickyCut |
177 |
|
178 |
function getStickyPowerCut(atomID) result(cutValue) |
179 |
integer, intent(in) :: atomID |
180 |
real(kind=dp) :: cutValue |
181 |
|
182 |
cutValue = StickyMap(atomID)%rbig |
183 |
end function getStickyPowerCut |
184 |
|
185 |
subroutine do_sticky_pair(me1, me2, d, rij, r2, sw, vpair, & |
186 |
pot, A1, A2, f1, t1, t2) |
187 |
|
188 |
!! This routine does only the sticky portion of the SSD potential |
189 |
!! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)]. |
190 |
!! The Lennard-Jones and dipolar interaction must be handled separately. |
191 |
|
192 |
!! We assume that the rotation matrices have already been calculated |
193 |
!! and placed in the A array. |
194 |
|
195 |
!! i and j are pointers to the two SSD atoms |
196 |
|
197 |
integer, intent(in) :: me1, me2 |
198 |
real (kind=dp), intent(inout) :: rij, r2 |
199 |
real (kind=dp), dimension(3), intent(in) :: d |
200 |
real (kind=dp) :: pot, vpair, sw |
201 |
real (kind=dp), dimension(9) :: A1, A2 |
202 |
real (kind=dp), dimension(3) :: f1 |
203 |
real (kind=dp), dimension(3) :: t1, t2 |
204 |
|
205 |
real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 |
206 |
real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr |
207 |
real (kind=dp) :: wi, wj, w, wip, wjp, wp |
208 |
real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz |
209 |
real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz |
210 |
real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz |
211 |
real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz |
212 |
real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj |
213 |
real (kind=dp) :: drdx, drdy, drdz |
214 |
real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj |
215 |
real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj |
216 |
real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji |
217 |
real (kind=dp) :: fxradial, fyradial, fzradial |
218 |
real (kind=dp) :: rijtest, rjitest |
219 |
real (kind=dp) :: radcomxi, radcomyi, radcomzi |
220 |
real (kind=dp) :: radcomxj, radcomyj, radcomzj |
221 |
integer :: id1, id2 |
222 |
|
223 |
real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig, dx |
224 |
|
225 |
if (me1.eq.me2) then |
226 |
w0 = StickyMap(me1)%w0 |
227 |
v0 = StickyMap(me1)%v0 |
228 |
v0p = StickyMap(me1)%v0p |
229 |
rl = StickyMap(me1)%rl |
230 |
ru = StickyMap(me1)%ru |
231 |
rlp = StickyMap(me1)%rlp |
232 |
rup = StickyMap(me1)%rup |
233 |
rbig = StickyMap(me1)%rbig |
234 |
else |
235 |
! This is silly, but if you want 2 sticky types in your |
236 |
! simulation, we'll let you do it with the Lorentz- |
237 |
! Berthelot mixing rules. |
238 |
! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW) |
239 |
rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl ) |
240 |
ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru ) |
241 |
rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp ) |
242 |
rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup ) |
243 |
rbig = max(ru, rup) |
244 |
w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 ) |
245 |
v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 ) |
246 |
v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p ) |
247 |
endif |
248 |
|
249 |
if ( rij .LE. rbig ) then |
250 |
|
251 |
r3 = r2*rij |
252 |
r5 = r3*r2 |
253 |
|
254 |
drdx = d(1) / rij |
255 |
drdy = d(2) / rij |
256 |
drdz = d(3) / rij |
257 |
|
258 |
! rotate the inter-particle separation into the two different |
259 |
! body-fixed coordinate systems: |
260 |
|
261 |
xi = A1(1)*d(1) + A1(2)*d(2) + A1(3)*d(3) |
262 |
yi = A1(4)*d(1) + A1(5)*d(2) + A1(6)*d(3) |
263 |
zi = A1(7)*d(1) + A1(8)*d(2) + A1(9)*d(3) |
264 |
|
265 |
! negative sign because this is the vector from j to i: |
266 |
|
267 |
xj = -(A2(1)*d(1) + A2(2)*d(2) + A2(3)*d(3)) |
268 |
yj = -(A2(4)*d(1) + A2(5)*d(2) + A2(6)*d(3)) |
269 |
zj = -(A2(7)*d(1) + A2(8)*d(2) + A2(9)*d(3)) |
270 |
|
271 |
xi2 = xi*xi |
272 |
yi2 = yi*yi |
273 |
zi2 = zi*zi |
274 |
|
275 |
xj2 = xj*xj |
276 |
yj2 = yj*yj |
277 |
zj2 = zj*zj |
278 |
|
279 |
! calculate the switching info. from the splines |
280 |
if (me1.eq.me2) then |
281 |
s = 0.0_dp |
282 |
dsdr = 0.0_dp |
283 |
sp = 0.0_dp |
284 |
dspdr = 0.0_dp |
285 |
|
286 |
if (rij.lt.ru) then |
287 |
if (rij.lt.rl) then |
288 |
s = 1.0_dp |
289 |
dsdr = 0.0_dp |
290 |
else |
291 |
! we are in the switching region |
292 |
dx = rij - rl |
293 |
s = StickyMap(me1)%stickySpline%y(1) + & |
294 |
dx*(dx*(StickyMap(me1)%stickySpline%c(1) + & |
295 |
dx*StickyMap(me1)%stickySpline%d(1))) |
296 |
dsdr = dx*(2.0_dp * StickyMap(me1)%stickySpline%c(1) + & |
297 |
3.0_dp * dx * StickyMap(me1)%stickySpline%d(1)) |
298 |
endif |
299 |
endif |
300 |
if (rij.lt.rup) then |
301 |
if (rij.lt.rlp) then |
302 |
sp = 1.0_dp |
303 |
dspdr = 0.0_dp |
304 |
else |
305 |
! we are in the switching region |
306 |
dx = rij - rlp |
307 |
sp = StickyMap(me1)%stickySplineP%y(1) + & |
308 |
dx*(dx*(StickyMap(me1)%stickySplineP%c(1) + & |
309 |
dx*StickyMap(me1)%stickySplineP%d(1))) |
310 |
dspdr = dx*(2.0_dp * StickyMap(me1)%stickySplineP%c(1) + & |
311 |
3.0_dp * dx * StickyMap(me1)%stickySplineP%d(1)) |
312 |
endif |
313 |
endif |
314 |
else |
315 |
! calculate the switching function explicitly rather than from |
316 |
! the splines with mixed sticky maps |
317 |
call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr) |
318 |
endif |
319 |
|
320 |
wi = 2.0_dp*(xi2-yi2)*zi / r3 |
321 |
wj = 2.0_dp*(xj2-yj2)*zj / r3 |
322 |
w = wi+wj |
323 |
|
324 |
zif = zi/rij - 0.6_dp |
325 |
zis = zi/rij + 0.8_dp |
326 |
|
327 |
zjf = zj/rij - 0.6_dp |
328 |
zjs = zj/rij + 0.8_dp |
329 |
|
330 |
wip = zif*zif*zis*zis - w0 |
331 |
wjp = zjf*zjf*zjs*zjs - w0 |
332 |
wp = wip + wjp |
333 |
|
334 |
vpair = vpair + 0.5_dp*(v0*s*w + v0p*sp*wp) |
335 |
|
336 |
pot = pot + 0.5_dp*(v0*s*w + v0p*sp*wp)*sw |
337 |
|
338 |
dwidx = 4.0_dp*xi*zi/r3 - 6.0_dp*xi*zi*(xi2-yi2)/r5 |
339 |
dwidy = - 4.0_dp*yi*zi/r3 - 6.0_dp*yi*zi*(xi2-yi2)/r5 |
340 |
dwidz = 2.0_dp*(xi2-yi2)/r3 - 6.0_dp*zi2*(xi2-yi2)/r5 |
341 |
|
342 |
dwjdx = 4.0_dp*xj*zj/r3 - 6.0_dp*xj*zj*(xj2-yj2)/r5 |
343 |
dwjdy = - 4.0_dp*yj*zj/r3 - 6.0_dp*yj*zj*(xj2-yj2)/r5 |
344 |
dwjdz = 2.0_dp*(xj2-yj2)/r3 - 6.0_dp*zj2*(xj2-yj2)/r5 |
345 |
|
346 |
uglyi = zif*zif*zis + zif*zis*zis |
347 |
uglyj = zjf*zjf*zjs + zjf*zjs*zjs |
348 |
|
349 |
dwipdx = -2.0_dp*xi*zi*uglyi/r3 |
350 |
dwipdy = -2.0_dp*yi*zi*uglyi/r3 |
351 |
dwipdz = 2.0_dp*(1.0_dp/rij - zi2/r3)*uglyi |
352 |
|
353 |
dwjpdx = -2.0_dp*xj*zj*uglyj/r3 |
354 |
dwjpdy = -2.0_dp*yj*zj*uglyj/r3 |
355 |
dwjpdz = 2.0_dp*(1.0_dp/rij - zj2/r3)*uglyj |
356 |
|
357 |
dwidux = 4.0_dp*(yi*zi2 + 0.5_dp*yi*(xi2-yi2))/r3 |
358 |
dwiduy = 4.0_dp*(xi*zi2 - 0.5_dp*xi*(xi2-yi2))/r3 |
359 |
dwiduz = - 8.0_dp*xi*yi*zi/r3 |
360 |
|
361 |
dwjdux = 4.0_dp*(yj*zj2 + 0.5_dp*yj*(xj2-yj2))/r3 |
362 |
dwjduy = 4.0_dp*(xj*zj2 - 0.5_dp*xj*(xj2-yj2))/r3 |
363 |
dwjduz = - 8.0_dp*xj*yj*zj/r3 |
364 |
|
365 |
dwipdux = 2.0_dp*yi*uglyi/rij |
366 |
dwipduy = -2.0_dp*xi*uglyi/rij |
367 |
dwipduz = 0.0_dp |
368 |
|
369 |
dwjpdux = 2.0_dp*yj*uglyj/rij |
370 |
dwjpduy = -2.0_dp*xj*uglyj/rij |
371 |
dwjpduz = 0.0_dp |
372 |
|
373 |
! do the torques first since they are easy: |
374 |
! remember that these are still in the body fixed axes |
375 |
|
376 |
txi = 0.5_dp*(v0*s*dwidux + v0p*sp*dwipdux)*sw |
377 |
tyi = 0.5_dp*(v0*s*dwiduy + v0p*sp*dwipduy)*sw |
378 |
tzi = 0.5_dp*(v0*s*dwiduz + v0p*sp*dwipduz)*sw |
379 |
|
380 |
txj = 0.5_dp*(v0*s*dwjdux + v0p*sp*dwjpdux)*sw |
381 |
tyj = 0.5_dp*(v0*s*dwjduy + v0p*sp*dwjpduy)*sw |
382 |
tzj = 0.5_dp*(v0*s*dwjduz + v0p*sp*dwjpduz)*sw |
383 |
|
384 |
! go back to lab frame using transpose of rotation matrix: |
385 |
|
386 |
t1(1) = t1(1) + a1(1)*txi + a1(4)*tyi + a1(7)*tzi |
387 |
t1(2) = t1(2) + a1(2)*txi + a1(5)*tyi + a1(8)*tzi |
388 |
t1(3) = t1(3) + a1(3)*txi + a1(6)*tyi + a1(9)*tzi |
389 |
|
390 |
t2(1) = t2(1) + a2(1)*txj + a2(4)*tyj + a2(7)*tzj |
391 |
t2(2) = t2(2) + a2(2)*txj + a2(5)*tyj + a2(8)*tzj |
392 |
t2(3) = t2(3) + a2(3)*txj + a2(6)*tyj + a2(9)*tzj |
393 |
|
394 |
! Now, on to the forces: |
395 |
|
396 |
! first rotate the i terms back into the lab frame: |
397 |
|
398 |
radcomxi = (v0*s*dwidx+v0p*sp*dwipdx)*sw |
399 |
radcomyi = (v0*s*dwidy+v0p*sp*dwipdy)*sw |
400 |
radcomzi = (v0*s*dwidz+v0p*sp*dwipdz)*sw |
401 |
|
402 |
radcomxj = (v0*s*dwjdx+v0p*sp*dwjpdx)*sw |
403 |
radcomyj = (v0*s*dwjdy+v0p*sp*dwjpdy)*sw |
404 |
radcomzj = (v0*s*dwjdz+v0p*sp*dwjpdz)*sw |
405 |
|
406 |
fxii = a1(1)*(radcomxi) + a1(4)*(radcomyi) + a1(7)*(radcomzi) |
407 |
fyii = a1(2)*(radcomxi) + a1(5)*(radcomyi) + a1(8)*(radcomzi) |
408 |
fzii = a1(3)*(radcomxi) + a1(6)*(radcomyi) + a1(9)*(radcomzi) |
409 |
|
410 |
fxjj = a2(1)*(radcomxj) + a2(4)*(radcomyj) + a2(7)*(radcomzj) |
411 |
fyjj = a2(2)*(radcomxj) + a2(5)*(radcomyj) + a2(8)*(radcomzj) |
412 |
fzjj = a2(3)*(radcomxj) + a2(6)*(radcomyj) + a2(9)*(radcomzj) |
413 |
|
414 |
fxij = -fxii |
415 |
fyij = -fyii |
416 |
fzij = -fzii |
417 |
|
418 |
fxji = -fxjj |
419 |
fyji = -fyjj |
420 |
fzji = -fzjj |
421 |
|
422 |
! now assemble these with the radial-only terms: |
423 |
|
424 |
fxradial = 0.5_dp*(v0*dsdr*drdx*w + v0p*dspdr*drdx*wp + fxii + fxji) |
425 |
fyradial = 0.5_dp*(v0*dsdr*drdy*w + v0p*dspdr*drdy*wp + fyii + fyji) |
426 |
fzradial = 0.5_dp*(v0*dsdr*drdz*w + v0p*dspdr*drdz*wp + fzii + fzji) |
427 |
|
428 |
f1(1) = f1(1) + fxradial |
429 |
f1(2) = f1(2) + fyradial |
430 |
f1(3) = f1(3) + fzradial |
431 |
|
432 |
endif |
433 |
end subroutine do_sticky_pair |
434 |
|
435 |
!! calculates the switching functions and their derivatives for a given |
436 |
subroutine calc_sw_fnc(r, rl, ru, rlp, rup, s, sp, dsdr, dspdr) |
437 |
|
438 |
real (kind=dp), intent(in) :: r, rl, ru, rlp, rup |
439 |
real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr |
440 |
|
441 |
! distances must be in angstroms |
442 |
s = 0.0_dp |
443 |
dsdr = 0.0_dp |
444 |
sp = 0.0_dp |
445 |
dspdr = 0.0_dp |
446 |
|
447 |
if (r.lt.ru) then |
448 |
if (r.lt.rl) then |
449 |
s = 1.0_dp |
450 |
dsdr = 0.0_dp |
451 |
else |
452 |
s = ((ru + 2.0_dp*r - 3.0_dp*rl) * (ru-r)**2) / & |
453 |
((ru - rl)**3) |
454 |
dsdr = 6.0_dp*(r-ru)*(r-rl)/((ru - rl)**3) |
455 |
endif |
456 |
endif |
457 |
|
458 |
if (r.lt.rup) then |
459 |
if (r.lt.rlp) then |
460 |
sp = 1.0_dp |
461 |
dspdr = 0.0_dp |
462 |
else |
463 |
sp = ((rup + 2.0_dp*r - 3.0_dp*rlp) * (rup-r)**2) / & |
464 |
((rup - rlp)**3) |
465 |
dspdr = 6.0_dp*(r-rup)*(r-rlp)/((rup - rlp)**3) |
466 |
endif |
467 |
endif |
468 |
|
469 |
return |
470 |
end subroutine calc_sw_fnc |
471 |
|
472 |
subroutine destroyStickyTypes() |
473 |
if(allocated(StickyMap)) deallocate(StickyMap) |
474 |
end subroutine destroyStickyTypes |
475 |
|
476 |
subroutine do_sticky_power_pair(me1, me2, d, rij, r2, sw, vpair, & |
477 |
pot, A1, A2, f1, t1, t2) |
478 |
|
479 |
!! i and j are pointers to the two SSD atoms |
480 |
|
481 |
real (kind=dp), intent(inout) :: rij, r2 |
482 |
real (kind=dp), dimension(3), intent(in) :: d |
483 |
real (kind=dp) :: pot, vpair, sw |
484 |
real (kind=dp), dimension(9) :: A1, A2 |
485 |
real (kind=dp), dimension(3) :: f1 |
486 |
real (kind=dp), dimension(3) :: t1, t2 |
487 |
|
488 |
|
489 |
real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 |
490 |
real (kind=dp) :: xihat, yihat, zihat, xjhat, yjhat, zjhat |
491 |
real (kind=dp) :: rI, rI2, rI3, rI4, rI5, rI6, rI7, s, sp, dsdr, dspdr |
492 |
real (kind=dp) :: wi, wj, w, wi2, wj2, eScale, v0scale |
493 |
real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz |
494 |
real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz |
495 |
real (kind=dp) :: drdx, drdy, drdz |
496 |
real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj |
497 |
real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj |
498 |
real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji |
499 |
real (kind=dp) :: fxradial, fyradial, fzradial |
500 |
real (kind=dp) :: rijtest, rjitest |
501 |
real (kind=dp) :: radcomxi, radcomyi, radcomzi |
502 |
real (kind=dp) :: radcomxj, radcomyj, radcomzj |
503 |
integer :: id1, id2 |
504 |
integer :: me1, me2 |
505 |
real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig |
506 |
real (kind=dp) :: zi3, zi4, zi5, zj3, zj4, zj5 |
507 |
real (kind=dp) :: frac1, frac2 |
508 |
|
509 |
if (.not.allocated(StickyMap)) then |
510 |
call handleError("sticky", "no StickyMap was present before first call of do_sticky_power_pair!") |
511 |
return |
512 |
end if |
513 |
|
514 |
if (me1.eq.me2) then |
515 |
w0 = StickyMap(me1)%w0 |
516 |
v0 = StickyMap(me1)%v0 |
517 |
v0p = StickyMap(me1)%v0p |
518 |
rl = StickyMap(me1)%rl |
519 |
ru = StickyMap(me1)%ru |
520 |
rlp = StickyMap(me1)%rlp |
521 |
rup = StickyMap(me1)%rup |
522 |
rbig = StickyMap(me1)%rbig |
523 |
else |
524 |
! This is silly, but if you want 2 sticky types in your |
525 |
! simulation, we'll let you do it with the Lorentz- |
526 |
! Berthelot mixing rules. |
527 |
! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW) |
528 |
rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl ) |
529 |
ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru ) |
530 |
rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp ) |
531 |
rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup ) |
532 |
rbig = max(ru, rup) |
533 |
w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 ) |
534 |
v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 ) |
535 |
v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p ) |
536 |
endif |
537 |
|
538 |
if ( rij .LE. rbig ) then |
539 |
|
540 |
rI = 1.0_dp/rij |
541 |
rI2 = rI*rI |
542 |
rI3 = rI2*rI |
543 |
rI4 = rI2*rI2 |
544 |
rI5 = rI3*rI2 |
545 |
rI6 = rI3*rI3 |
546 |
rI7 = rI4*rI3 |
547 |
|
548 |
drdx = d(1) * rI |
549 |
drdy = d(2) * rI |
550 |
drdz = d(3) * rI |
551 |
|
552 |
! rotate the inter-particle separation into the two different |
553 |
! body-fixed coordinate systems: |
554 |
|
555 |
xi = A1(1)*d(1) + A1(2)*d(2) + A1(3)*d(3) |
556 |
yi = A1(4)*d(1) + A1(5)*d(2) + A1(6)*d(3) |
557 |
zi = A1(7)*d(1) + A1(8)*d(2) + A1(9)*d(3) |
558 |
|
559 |
! negative sign because this is the vector from j to i: |
560 |
|
561 |
xj = -(A2(1)*d(1) + A2(2)*d(2) + A2(3)*d(3)) |
562 |
yj = -(A2(4)*d(1) + A2(5)*d(2) + A2(6)*d(3)) |
563 |
zj = -(A2(7)*d(1) + A2(8)*d(2) + A2(9)*d(3)) |
564 |
|
565 |
xi2 = xi*xi |
566 |
yi2 = yi*yi |
567 |
zi2 = zi*zi |
568 |
zi3 = zi2*zi |
569 |
zi4 = zi2*zi2 |
570 |
zi5 = zi3*zi2 |
571 |
xihat = xi*rI |
572 |
yihat = yi*rI |
573 |
zihat = zi*rI |
574 |
|
575 |
xj2 = xj*xj |
576 |
yj2 = yj*yj |
577 |
zj2 = zj*zj |
578 |
zj3 = zj2*zj |
579 |
zj4 = zj2*zj2 |
580 |
zj5 = zj3*zj2 |
581 |
xjhat = xj*rI |
582 |
yjhat = yj*rI |
583 |
zjhat = zj*rI |
584 |
|
585 |
call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr) |
586 |
|
587 |
frac1 = 0.25_dp |
588 |
frac2 = 0.75_dp |
589 |
|
590 |
wi = 2.0_dp*(xi2-yi2)*zi*rI3 |
591 |
wj = 2.0_dp*(xj2-yj2)*zj*rI3 |
592 |
|
593 |
wi2 = wi*wi |
594 |
wj2 = wj*wj |
595 |
|
596 |
w = frac1*wi*wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p |
597 |
|
598 |
vpair = vpair + 0.5_dp*(v0*s*w) |
599 |
|
600 |
pot = pot + 0.5_dp*(v0*s*w)*sw |
601 |
|
602 |
dwidx = ( 4.0_dp*xi*zi*rI3 - 6.0_dp*xi*zi*(xi2-yi2)*rI5 ) |
603 |
dwidy = ( -4.0_dp*yi*zi*rI3 - 6.0_dp*yi*zi*(xi2-yi2)*rI5 ) |
604 |
dwidz = ( 2.0_dp*(xi2-yi2)*rI3 - 6.0_dp*zi2*(xi2-yi2)*rI5 ) |
605 |
|
606 |
dwidx = frac1*3.0_dp*wi2*dwidx + frac2*dwidx |
607 |
dwidy = frac1*3.0_dp*wi2*dwidy + frac2*dwidy |
608 |
dwidz = frac1*3.0_dp*wi2*dwidz + frac2*dwidz |
609 |
|
610 |
dwjdx = ( 4.0_dp*xj*zj*rI3 - 6.0_dp*xj*zj*(xj2-yj2)*rI5 ) |
611 |
dwjdy = ( -4.0_dp*yj*zj*rI3 - 6.0_dp*yj*zj*(xj2-yj2)*rI5 ) |
612 |
dwjdz = ( 2.0_dp*(xj2-yj2)*rI3 - 6.0_dp*zj2*(xj2-yj2)*rI5 ) |
613 |
|
614 |
dwjdx = frac1*3.0_dp*wj2*dwjdx + frac2*dwjdx |
615 |
dwjdy = frac1*3.0_dp*wj2*dwjdy + frac2*dwjdy |
616 |
dwjdz = frac1*3.0_dp*wj2*dwjdz + frac2*dwjdz |
617 |
|
618 |
dwidux = ( 4.0_dp*(yi*zi2 + 0.5_dp*yi*(xi2-yi2))*rI3 ) |
619 |
dwiduy = ( 4.0_dp*(xi*zi2 - 0.5_dp*xi*(xi2-yi2))*rI3 ) |
620 |
dwiduz = ( -8.0_dp*xi*yi*zi*rI3 ) |
621 |
|
622 |
dwidux = frac1*3.0_dp*wi2*dwidux + frac2*dwidux |
623 |
dwiduy = frac1*3.0_dp*wi2*dwiduy + frac2*dwiduy |
624 |
dwiduz = frac1*3.0_dp*wi2*dwiduz + frac2*dwiduz |
625 |
|
626 |
dwjdux = ( 4.0_dp*(yj*zj2 + 0.5_dp*yj*(xj2-yj2))*rI3 ) |
627 |
dwjduy = ( 4.0_dp*(xj*zj2 - 0.5_dp*xj*(xj2-yj2))*rI3 ) |
628 |
dwjduz = ( -8.0_dp*xj*yj*zj*rI3 ) |
629 |
|
630 |
dwjdux = frac1*3.0_dp*wj2*dwjdux + frac2*dwjdux |
631 |
dwjduy = frac1*3.0_dp*wj2*dwjduy + frac2*dwjduy |
632 |
dwjduz = frac1*3.0_dp*wj2*dwjduz + frac2*dwjduz |
633 |
|
634 |
! do the torques first since they are easy: |
635 |
! remember that these are still in the body fixed axes |
636 |
|
637 |
txi = 0.5_dp*(v0*s*dwidux)*sw |
638 |
tyi = 0.5_dp*(v0*s*dwiduy)*sw |
639 |
tzi = 0.5_dp*(v0*s*dwiduz)*sw |
640 |
|
641 |
txj = 0.5_dp*(v0*s*dwjdux)*sw |
642 |
tyj = 0.5_dp*(v0*s*dwjduy)*sw |
643 |
tzj = 0.5_dp*(v0*s*dwjduz)*sw |
644 |
|
645 |
! go back to lab frame using transpose of rotation matrix: |
646 |
|
647 |
t1(1) = t1(1) + a1(1)*txi + a1(4)*tyi + a1(7)*tzi |
648 |
t1(2) = t1(2) + a1(2)*txi + a1(5)*tyi + a1(8)*tzi |
649 |
t1(3) = t1(3) + a1(3)*txi + a1(6)*tyi + a1(9)*tzi |
650 |
|
651 |
t2(1) = t2(1) + a2(1)*txj + a2(4)*tyj + a2(7)*tzj |
652 |
t2(2) = t2(2) + a2(2)*txj + a2(5)*tyj + a2(8)*tzj |
653 |
t2(3) = t2(3) + a2(3)*txj + a2(6)*tyj + a2(9)*tzj |
654 |
|
655 |
! Now, on to the forces: |
656 |
|
657 |
! first rotate the i terms back into the lab frame: |
658 |
|
659 |
radcomxi = (v0*s*dwidx)*sw |
660 |
radcomyi = (v0*s*dwidy)*sw |
661 |
radcomzi = (v0*s*dwidz)*sw |
662 |
|
663 |
radcomxj = (v0*s*dwjdx)*sw |
664 |
radcomyj = (v0*s*dwjdy)*sw |
665 |
radcomzj = (v0*s*dwjdz)*sw |
666 |
|
667 |
fxii = a1(1)*(radcomxi) + a1(4)*(radcomyi) + a1(7)*(radcomzi) |
668 |
fyii = a1(2)*(radcomxi) + a1(5)*(radcomyi) + a1(8)*(radcomzi) |
669 |
fzii = a1(3)*(radcomxi) + a1(6)*(radcomyi) + a1(9)*(radcomzi) |
670 |
|
671 |
fxjj = a2(1)*(radcomxj) + a2(4)*(radcomyj) + a2(7)*(radcomzj) |
672 |
fyjj = a2(2)*(radcomxj) + a2(5)*(radcomyj) + a2(8)*(radcomzj) |
673 |
fzjj = a2(3)*(radcomxj) + a2(6)*(radcomyj) + a2(9)*(radcomzj) |
674 |
|
675 |
fxij = -fxii |
676 |
fyij = -fyii |
677 |
fzij = -fzii |
678 |
|
679 |
fxji = -fxjj |
680 |
fyji = -fyjj |
681 |
fzji = -fzjj |
682 |
|
683 |
! now assemble these with the radial-only terms: |
684 |
|
685 |
fxradial = 0.5_dp*(v0*dsdr*w*drdx + fxii + fxji) |
686 |
fyradial = 0.5_dp*(v0*dsdr*w*drdy + fyii + fyji) |
687 |
fzradial = 0.5_dp*(v0*dsdr*w*drdz + fzii + fzji) |
688 |
|
689 |
f1(1) = f1(1) + fxradial |
690 |
f1(2) = f1(2) + fyradial |
691 |
f1(3) = f1(3) + fzradial |
692 |
|
693 |
endif |
694 |
end subroutine do_sticky_power_pair |
695 |
|
696 |
end module sticky |