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