ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/sticky.F90
Revision: 1464
Committed: Fri Jul 9 19:29:05 2010 UTC (15 years, 1 month ago) by gezelter
File size: 23547 byte(s)
Log Message:
removing cruft (atom numbers, do_pot, do_stress) from many modules and
force managers

File Contents

# Content
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, fpair, &
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), dimension(3), intent(inout) :: fpair
201 real (kind=dp) :: pot, vpair, sw
202 real (kind=dp), dimension(9) :: A1, A2
203 real (kind=dp), dimension(3) :: f1
204 real (kind=dp), dimension(3) :: t1, t2
205
206 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
207 real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
208 real (kind=dp) :: wi, wj, w, wip, wjp, wp
209 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
210 real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
211 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
212 real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
213 real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
214 real (kind=dp) :: drdx, drdy, drdz
215 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
216 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
217 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
218 real (kind=dp) :: fxradial, fyradial, fzradial
219 real (kind=dp) :: rijtest, rjitest
220 real (kind=dp) :: radcomxi, radcomyi, radcomzi
221 real (kind=dp) :: radcomxj, radcomyj, radcomzj
222 integer :: id1, id2
223
224 real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig, dx
225
226 if (me1.eq.me2) then
227 w0 = StickyMap(me1)%w0
228 v0 = StickyMap(me1)%v0
229 v0p = StickyMap(me1)%v0p
230 rl = StickyMap(me1)%rl
231 ru = StickyMap(me1)%ru
232 rlp = StickyMap(me1)%rlp
233 rup = StickyMap(me1)%rup
234 rbig = StickyMap(me1)%rbig
235 else
236 ! This is silly, but if you want 2 sticky types in your
237 ! simulation, we'll let you do it with the Lorentz-
238 ! Berthelot mixing rules.
239 ! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW)
240 rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl )
241 ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru )
242 rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp )
243 rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup )
244 rbig = max(ru, rup)
245 w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 )
246 v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 )
247 v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p )
248 endif
249
250 if ( rij .LE. rbig ) then
251
252 r3 = r2*rij
253 r5 = r3*r2
254
255 drdx = d(1) / rij
256 drdy = d(2) / rij
257 drdz = d(3) / rij
258
259 ! rotate the inter-particle separation into the two different
260 ! body-fixed coordinate systems:
261
262 xi = A1(1)*d(1) + A1(2)*d(2) + A1(3)*d(3)
263 yi = A1(4)*d(1) + A1(5)*d(2) + A1(6)*d(3)
264 zi = A1(7)*d(1) + A1(8)*d(2) + A1(9)*d(3)
265
266 ! negative sign because this is the vector from j to i:
267
268 xj = -(A2(1)*d(1) + A2(2)*d(2) + A2(3)*d(3))
269 yj = -(A2(4)*d(1) + A2(5)*d(2) + A2(6)*d(3))
270 zj = -(A2(7)*d(1) + A2(8)*d(2) + A2(9)*d(3))
271
272 xi2 = xi*xi
273 yi2 = yi*yi
274 zi2 = zi*zi
275
276 xj2 = xj*xj
277 yj2 = yj*yj
278 zj2 = zj*zj
279
280 ! calculate the switching info. from the splines
281 if (me1.eq.me2) then
282 s = 0.0_dp
283 dsdr = 0.0_dp
284 sp = 0.0_dp
285 dspdr = 0.0_dp
286
287 if (rij.lt.ru) then
288 if (rij.lt.rl) then
289 s = 1.0_dp
290 dsdr = 0.0_dp
291 else
292 ! we are in the switching region
293 dx = rij - rl
294 s = StickyMap(me1)%stickySpline%y(1) + &
295 dx*(dx*(StickyMap(me1)%stickySpline%c(1) + &
296 dx*StickyMap(me1)%stickySpline%d(1)))
297 dsdr = dx*(2.0_dp * StickyMap(me1)%stickySpline%c(1) + &
298 3.0_dp * dx * StickyMap(me1)%stickySpline%d(1))
299 endif
300 endif
301 if (rij.lt.rup) then
302 if (rij.lt.rlp) then
303 sp = 1.0_dp
304 dspdr = 0.0_dp
305 else
306 ! we are in the switching region
307 dx = rij - rlp
308 sp = StickyMap(me1)%stickySplineP%y(1) + &
309 dx*(dx*(StickyMap(me1)%stickySplineP%c(1) + &
310 dx*StickyMap(me1)%stickySplineP%d(1)))
311 dspdr = dx*(2.0_dp * StickyMap(me1)%stickySplineP%c(1) + &
312 3.0_dp * dx * StickyMap(me1)%stickySplineP%d(1))
313 endif
314 endif
315 else
316 ! calculate the switching function explicitly rather than from
317 ! the splines with mixed sticky maps
318 call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
319 endif
320
321 wi = 2.0_dp*(xi2-yi2)*zi / r3
322 wj = 2.0_dp*(xj2-yj2)*zj / r3
323 w = wi+wj
324
325 zif = zi/rij - 0.6_dp
326 zis = zi/rij + 0.8_dp
327
328 zjf = zj/rij - 0.6_dp
329 zjs = zj/rij + 0.8_dp
330
331 wip = zif*zif*zis*zis - w0
332 wjp = zjf*zjf*zjs*zjs - w0
333 wp = wip + wjp
334
335 vpair = vpair + 0.5_dp*(v0*s*w + v0p*sp*wp)
336
337 pot = pot + 0.5_dp*(v0*s*w + v0p*sp*wp)*sw
338
339 dwidx = 4.0_dp*xi*zi/r3 - 6.0_dp*xi*zi*(xi2-yi2)/r5
340 dwidy = - 4.0_dp*yi*zi/r3 - 6.0_dp*yi*zi*(xi2-yi2)/r5
341 dwidz = 2.0_dp*(xi2-yi2)/r3 - 6.0_dp*zi2*(xi2-yi2)/r5
342
343 dwjdx = 4.0_dp*xj*zj/r3 - 6.0_dp*xj*zj*(xj2-yj2)/r5
344 dwjdy = - 4.0_dp*yj*zj/r3 - 6.0_dp*yj*zj*(xj2-yj2)/r5
345 dwjdz = 2.0_dp*(xj2-yj2)/r3 - 6.0_dp*zj2*(xj2-yj2)/r5
346
347 uglyi = zif*zif*zis + zif*zis*zis
348 uglyj = zjf*zjf*zjs + zjf*zjs*zjs
349
350 dwipdx = -2.0_dp*xi*zi*uglyi/r3
351 dwipdy = -2.0_dp*yi*zi*uglyi/r3
352 dwipdz = 2.0_dp*(1.0_dp/rij - zi2/r3)*uglyi
353
354 dwjpdx = -2.0_dp*xj*zj*uglyj/r3
355 dwjpdy = -2.0_dp*yj*zj*uglyj/r3
356 dwjpdz = 2.0_dp*(1.0_dp/rij - zj2/r3)*uglyj
357
358 dwidux = 4.0_dp*(yi*zi2 + 0.5_dp*yi*(xi2-yi2))/r3
359 dwiduy = 4.0_dp*(xi*zi2 - 0.5_dp*xi*(xi2-yi2))/r3
360 dwiduz = - 8.0_dp*xi*yi*zi/r3
361
362 dwjdux = 4.0_dp*(yj*zj2 + 0.5_dp*yj*(xj2-yj2))/r3
363 dwjduy = 4.0_dp*(xj*zj2 - 0.5_dp*xj*(xj2-yj2))/r3
364 dwjduz = - 8.0_dp*xj*yj*zj/r3
365
366 dwipdux = 2.0_dp*yi*uglyi/rij
367 dwipduy = -2.0_dp*xi*uglyi/rij
368 dwipduz = 0.0_dp
369
370 dwjpdux = 2.0_dp*yj*uglyj/rij
371 dwjpduy = -2.0_dp*xj*uglyj/rij
372 dwjpduz = 0.0_dp
373
374 ! do the torques first since they are easy:
375 ! remember that these are still in the body fixed axes
376
377 txi = 0.5_dp*(v0*s*dwidux + v0p*sp*dwipdux)*sw
378 tyi = 0.5_dp*(v0*s*dwiduy + v0p*sp*dwipduy)*sw
379 tzi = 0.5_dp*(v0*s*dwiduz + v0p*sp*dwipduz)*sw
380
381 txj = 0.5_dp*(v0*s*dwjdux + v0p*sp*dwjpdux)*sw
382 tyj = 0.5_dp*(v0*s*dwjduy + v0p*sp*dwjpduy)*sw
383 tzj = 0.5_dp*(v0*s*dwjduz + v0p*sp*dwjpduz)*sw
384
385 ! go back to lab frame using transpose of rotation matrix:
386
387 t1(1) = t1(1) + a1(1)*txi + a1(4)*tyi + a1(7)*tzi
388 t1(2) = t1(2) + a1(2)*txi + a1(5)*tyi + a1(8)*tzi
389 t1(3) = t1(3) + a1(3)*txi + a1(6)*tyi + a1(9)*tzi
390
391 t2(1) = t2(1) + a2(1)*txj + a2(4)*tyj + a2(7)*tzj
392 t2(2) = t2(2) + a2(2)*txj + a2(5)*tyj + a2(8)*tzj
393 t2(3) = t2(3) + a2(3)*txj + a2(6)*tyj + a2(9)*tzj
394
395 ! Now, on to the forces:
396
397 ! first rotate the i terms back into the lab frame:
398
399 radcomxi = (v0*s*dwidx+v0p*sp*dwipdx)*sw
400 radcomyi = (v0*s*dwidy+v0p*sp*dwipdy)*sw
401 radcomzi = (v0*s*dwidz+v0p*sp*dwipdz)*sw
402
403 radcomxj = (v0*s*dwjdx+v0p*sp*dwjpdx)*sw
404 radcomyj = (v0*s*dwjdy+v0p*sp*dwjpdy)*sw
405 radcomzj = (v0*s*dwjdz+v0p*sp*dwjpdz)*sw
406
407 fxii = a1(1)*(radcomxi) + a1(4)*(radcomyi) + a1(7)*(radcomzi)
408 fyii = a1(2)*(radcomxi) + a1(5)*(radcomyi) + a1(8)*(radcomzi)
409 fzii = a1(3)*(radcomxi) + a1(6)*(radcomyi) + a1(9)*(radcomzi)
410
411 fxjj = a2(1)*(radcomxj) + a2(4)*(radcomyj) + a2(7)*(radcomzj)
412 fyjj = a2(2)*(radcomxj) + a2(5)*(radcomyj) + a2(8)*(radcomzj)
413 fzjj = a2(3)*(radcomxj) + a2(6)*(radcomyj) + a2(9)*(radcomzj)
414
415 fxij = -fxii
416 fyij = -fyii
417 fzij = -fzii
418
419 fxji = -fxjj
420 fyji = -fyjj
421 fzji = -fzjj
422
423 ! now assemble these with the radial-only terms:
424
425 fxradial = 0.5_dp*(v0*dsdr*drdx*w + v0p*dspdr*drdx*wp + fxii + fxji)
426 fyradial = 0.5_dp*(v0*dsdr*drdy*w + v0p*dspdr*drdy*wp + fyii + fyji)
427 fzradial = 0.5_dp*(v0*dsdr*drdz*w + v0p*dspdr*drdz*wp + fzii + fzji)
428
429 f1(1) = f1(1) + fxradial
430 f1(2) = f1(2) + fyradial
431 f1(3) = f1(3) + fzradial
432
433 endif
434 end subroutine do_sticky_pair
435
436 !! calculates the switching functions and their derivatives for a given
437 subroutine calc_sw_fnc(r, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
438
439 real (kind=dp), intent(in) :: r, rl, ru, rlp, rup
440 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
441
442 ! distances must be in angstroms
443 s = 0.0_dp
444 dsdr = 0.0_dp
445 sp = 0.0_dp
446 dspdr = 0.0_dp
447
448 if (r.lt.ru) then
449 if (r.lt.rl) then
450 s = 1.0_dp
451 dsdr = 0.0_dp
452 else
453 s = ((ru + 2.0_dp*r - 3.0_dp*rl) * (ru-r)**2) / &
454 ((ru - rl)**3)
455 dsdr = 6.0_dp*(r-ru)*(r-rl)/((ru - rl)**3)
456 endif
457 endif
458
459 if (r.lt.rup) then
460 if (r.lt.rlp) then
461 sp = 1.0_dp
462 dspdr = 0.0_dp
463 else
464 sp = ((rup + 2.0_dp*r - 3.0_dp*rlp) * (rup-r)**2) / &
465 ((rup - rlp)**3)
466 dspdr = 6.0_dp*(r-rup)*(r-rlp)/((rup - rlp)**3)
467 endif
468 endif
469
470 return
471 end subroutine calc_sw_fnc
472
473 subroutine destroyStickyTypes()
474 if(allocated(StickyMap)) deallocate(StickyMap)
475 end subroutine destroyStickyTypes
476
477 subroutine do_sticky_power_pair(me1, me2, d, rij, r2, sw, vpair, fpair, &
478 pot, A1, A2, f1, t1, t2)
479
480 !! i and j are pointers to the two SSD atoms
481
482 real (kind=dp), intent(inout) :: rij, r2
483 real (kind=dp), dimension(3), intent(in) :: d
484 real (kind=dp), dimension(3), intent(inout) :: fpair
485 real (kind=dp) :: pot, vpair, sw
486 real (kind=dp), dimension(9) :: A1, A2
487 real (kind=dp), dimension(3) :: f1
488 real (kind=dp), dimension(3) :: t1, t2
489
490
491 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
492 real (kind=dp) :: xihat, yihat, zihat, xjhat, yjhat, zjhat
493 real (kind=dp) :: rI, rI2, rI3, rI4, rI5, rI6, rI7, s, sp, dsdr, dspdr
494 real (kind=dp) :: wi, wj, w, wi2, wj2, eScale, v0scale
495 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
496 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
497 real (kind=dp) :: drdx, drdy, drdz
498 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
499 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
500 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
501 real (kind=dp) :: fxradial, fyradial, fzradial
502 real (kind=dp) :: rijtest, rjitest
503 real (kind=dp) :: radcomxi, radcomyi, radcomzi
504 real (kind=dp) :: radcomxj, radcomyj, radcomzj
505 integer :: id1, id2
506 integer :: me1, me2
507 real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig
508 real (kind=dp) :: zi3, zi4, zi5, zj3, zj4, zj5
509 real (kind=dp) :: frac1, frac2
510
511 if (.not.allocated(StickyMap)) then
512 call handleError("sticky", "no StickyMap was present before first call of do_sticky_power_pair!")
513 return
514 end if
515
516 if (me1.eq.me2) then
517 w0 = StickyMap(me1)%w0
518 v0 = StickyMap(me1)%v0
519 v0p = StickyMap(me1)%v0p
520 rl = StickyMap(me1)%rl
521 ru = StickyMap(me1)%ru
522 rlp = StickyMap(me1)%rlp
523 rup = StickyMap(me1)%rup
524 rbig = StickyMap(me1)%rbig
525 else
526 ! This is silly, but if you want 2 sticky types in your
527 ! simulation, we'll let you do it with the Lorentz-
528 ! Berthelot mixing rules.
529 ! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW)
530 rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl )
531 ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru )
532 rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp )
533 rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup )
534 rbig = max(ru, rup)
535 w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 )
536 v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 )
537 v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p )
538 endif
539
540 if ( rij .LE. rbig ) then
541
542 rI = 1.0_dp/rij
543 rI2 = rI*rI
544 rI3 = rI2*rI
545 rI4 = rI2*rI2
546 rI5 = rI3*rI2
547 rI6 = rI3*rI3
548 rI7 = rI4*rI3
549
550 drdx = d(1) * rI
551 drdy = d(2) * rI
552 drdz = d(3) * rI
553
554 ! rotate the inter-particle separation into the two different
555 ! body-fixed coordinate systems:
556
557 xi = A1(1)*d(1) + A1(2)*d(2) + A1(3)*d(3)
558 yi = A1(4)*d(1) + A1(5)*d(2) + A1(6)*d(3)
559 zi = A1(7)*d(1) + A1(8)*d(2) + A1(9)*d(3)
560
561 ! negative sign because this is the vector from j to i:
562
563 xj = -(A2(1)*d(1) + A2(2)*d(2) + A2(3)*d(3))
564 yj = -(A2(4)*d(1) + A2(5)*d(2) + A2(6)*d(3))
565 zj = -(A2(7)*d(1) + A2(8)*d(2) + A2(9)*d(3))
566
567 xi2 = xi*xi
568 yi2 = yi*yi
569 zi2 = zi*zi
570 zi3 = zi2*zi
571 zi4 = zi2*zi2
572 zi5 = zi3*zi2
573 xihat = xi*rI
574 yihat = yi*rI
575 zihat = zi*rI
576
577 xj2 = xj*xj
578 yj2 = yj*yj
579 zj2 = zj*zj
580 zj3 = zj2*zj
581 zj4 = zj2*zj2
582 zj5 = zj3*zj2
583 xjhat = xj*rI
584 yjhat = yj*rI
585 zjhat = zj*rI
586
587 call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
588
589 frac1 = 0.25_dp
590 frac2 = 0.75_dp
591
592 wi = 2.0_dp*(xi2-yi2)*zi*rI3
593 wj = 2.0_dp*(xj2-yj2)*zj*rI3
594
595 wi2 = wi*wi
596 wj2 = wj*wj
597
598 w = frac1*wi*wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p
599
600 vpair = vpair + 0.5_dp*(v0*s*w)
601
602 pot = pot + 0.5_dp*(v0*s*w)*sw
603
604 dwidx = ( 4.0_dp*xi*zi*rI3 - 6.0_dp*xi*zi*(xi2-yi2)*rI5 )
605 dwidy = ( -4.0_dp*yi*zi*rI3 - 6.0_dp*yi*zi*(xi2-yi2)*rI5 )
606 dwidz = ( 2.0_dp*(xi2-yi2)*rI3 - 6.0_dp*zi2*(xi2-yi2)*rI5 )
607
608 dwidx = frac1*3.0_dp*wi2*dwidx + frac2*dwidx
609 dwidy = frac1*3.0_dp*wi2*dwidy + frac2*dwidy
610 dwidz = frac1*3.0_dp*wi2*dwidz + frac2*dwidz
611
612 dwjdx = ( 4.0_dp*xj*zj*rI3 - 6.0_dp*xj*zj*(xj2-yj2)*rI5 )
613 dwjdy = ( -4.0_dp*yj*zj*rI3 - 6.0_dp*yj*zj*(xj2-yj2)*rI5 )
614 dwjdz = ( 2.0_dp*(xj2-yj2)*rI3 - 6.0_dp*zj2*(xj2-yj2)*rI5 )
615
616 dwjdx = frac1*3.0_dp*wj2*dwjdx + frac2*dwjdx
617 dwjdy = frac1*3.0_dp*wj2*dwjdy + frac2*dwjdy
618 dwjdz = frac1*3.0_dp*wj2*dwjdz + frac2*dwjdz
619
620 dwidux = ( 4.0_dp*(yi*zi2 + 0.5_dp*yi*(xi2-yi2))*rI3 )
621 dwiduy = ( 4.0_dp*(xi*zi2 - 0.5_dp*xi*(xi2-yi2))*rI3 )
622 dwiduz = ( -8.0_dp*xi*yi*zi*rI3 )
623
624 dwidux = frac1*3.0_dp*wi2*dwidux + frac2*dwidux
625 dwiduy = frac1*3.0_dp*wi2*dwiduy + frac2*dwiduy
626 dwiduz = frac1*3.0_dp*wi2*dwiduz + frac2*dwiduz
627
628 dwjdux = ( 4.0_dp*(yj*zj2 + 0.5_dp*yj*(xj2-yj2))*rI3 )
629 dwjduy = ( 4.0_dp*(xj*zj2 - 0.5_dp*xj*(xj2-yj2))*rI3 )
630 dwjduz = ( -8.0_dp*xj*yj*zj*rI3 )
631
632 dwjdux = frac1*3.0_dp*wj2*dwjdux + frac2*dwjdux
633 dwjduy = frac1*3.0_dp*wj2*dwjduy + frac2*dwjduy
634 dwjduz = frac1*3.0_dp*wj2*dwjduz + frac2*dwjduz
635
636 ! do the torques first since they are easy:
637 ! remember that these are still in the body fixed axes
638
639 txi = 0.5_dp*(v0*s*dwidux)*sw
640 tyi = 0.5_dp*(v0*s*dwiduy)*sw
641 tzi = 0.5_dp*(v0*s*dwiduz)*sw
642
643 txj = 0.5_dp*(v0*s*dwjdux)*sw
644 tyj = 0.5_dp*(v0*s*dwjduy)*sw
645 tzj = 0.5_dp*(v0*s*dwjduz)*sw
646
647 ! go back to lab frame using transpose of rotation matrix:
648
649 t1(1) = t1(1) + a1(1)*txi + a1(4)*tyi + a1(7)*tzi
650 t1(2) = t1(2) + a1(2)*txi + a1(5)*tyi + a1(8)*tzi
651 t1(3) = t1(3) + a1(3)*txi + a1(6)*tyi + a1(9)*tzi
652
653 t2(1) = t2(1) + a2(1)*txj + a2(4)*tyj + a2(7)*tzj
654 t2(2) = t2(2) + a2(2)*txj + a2(5)*tyj + a2(8)*tzj
655 t2(3) = t2(3) + a2(3)*txj + a2(6)*tyj + a2(9)*tzj
656
657 ! Now, on to the forces:
658
659 ! first rotate the i terms back into the lab frame:
660
661 radcomxi = (v0*s*dwidx)*sw
662 radcomyi = (v0*s*dwidy)*sw
663 radcomzi = (v0*s*dwidz)*sw
664
665 radcomxj = (v0*s*dwjdx)*sw
666 radcomyj = (v0*s*dwjdy)*sw
667 radcomzj = (v0*s*dwjdz)*sw
668
669 fxii = a1(1)*(radcomxi) + a1(4)*(radcomyi) + a1(7)*(radcomzi)
670 fyii = a1(2)*(radcomxi) + a1(5)*(radcomyi) + a1(8)*(radcomzi)
671 fzii = a1(3)*(radcomxi) + a1(6)*(radcomyi) + a1(9)*(radcomzi)
672
673 fxjj = a2(1)*(radcomxj) + a2(4)*(radcomyj) + a2(7)*(radcomzj)
674 fyjj = a2(2)*(radcomxj) + a2(5)*(radcomyj) + a2(8)*(radcomzj)
675 fzjj = a2(3)*(radcomxj) + a2(6)*(radcomyj) + a2(9)*(radcomzj)
676
677 fxij = -fxii
678 fyij = -fyii
679 fzij = -fzii
680
681 fxji = -fxjj
682 fyji = -fyjj
683 fzji = -fzjj
684
685 ! now assemble these with the radial-only terms:
686
687 fxradial = 0.5_dp*(v0*dsdr*w*drdx + fxii + fxji)
688 fyradial = 0.5_dp*(v0*dsdr*w*drdy + fyii + fyji)
689 fzradial = 0.5_dp*(v0*dsdr*w*drdz + fzii + fzji)
690
691 f1(1) = f1(1) + fxradial
692 f1(2) = f1(2) + fyradial
693 f1(3) = f1(3) + fzradial
694
695 endif
696 end subroutine do_sticky_power_pair
697
698 end module sticky

Properties

Name Value
svn:keywords Author Id Revision Date