ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/sticky.F90
Revision: 2251
Committed: Sun May 29 21:16:25 2005 UTC (19 years, 3 months ago) by chrisfen
File size: 28876 byte(s)
Log Message:
Removed balance from the Darkside (files)

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. 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 !! 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: sticky.F90,v 1.13 2005-05-29 21:16:25 chrisfen Exp $, $Date: 2005-05-29 21:16:25 $, $Name: not supported by cvs2svn $, $Revision: 1.13 $
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 #ifdef IS_MPI
64 use mpiSimulation
65 #endif
66 implicit none
67
68 PRIVATE
69
70 public :: newStickyType
71 public :: do_sticky_pair
72 public :: destroyStickyTypes
73 public :: do_sticky_power_pair
74
75
76 type :: StickyList
77 integer :: c_ident
78 real( kind = dp ) :: w0 = 0.0_dp
79 real( kind = dp ) :: v0 = 0.0_dp
80 real( kind = dp ) :: v0p = 0.0_dp
81 real( kind = dp ) :: rl = 0.0_dp
82 real( kind = dp ) :: ru = 0.0_dp
83 real( kind = dp ) :: rlp = 0.0_dp
84 real( kind = dp ) :: rup = 0.0_dp
85 real( kind = dp ) :: rbig = 0.0_dp
86 end type StickyList
87
88 type(StickyList), dimension(:),allocatable :: StickyMap
89
90 contains
91
92 subroutine newStickyType(c_ident, w0, v0, v0p, rl, ru, rlp, rup, isError)
93
94 integer, intent(in) :: c_ident
95 integer, intent(inout) :: isError
96 real( kind = dp ), intent(in) :: w0, v0, v0p
97 real( kind = dp ), intent(in) :: rl, ru
98 real( kind = dp ), intent(in) :: rlp, rup
99 integer :: nATypes, myATID
100
101
102 isError = 0
103 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
104
105 !! Be simple-minded and assume that we need a StickyMap that
106 !! is the same size as the total number of atom types
107
108 if (.not.allocated(StickyMap)) then
109
110 nAtypes = getSize(atypes)
111
112 if (nAtypes == 0) then
113 isError = -1
114 return
115 end if
116
117 if (.not. allocated(StickyMap)) then
118 allocate(StickyMap(nAtypes))
119 endif
120
121 end if
122
123 if (myATID .gt. size(StickyMap)) then
124 isError = -1
125 return
126 endif
127
128 ! set the values for StickyMap for this atom type:
129
130 StickyMap(myATID)%c_ident = c_ident
131
132 ! we could pass all 5 parameters if we felt like it...
133
134 StickyMap(myATID)%w0 = w0
135 StickyMap(myATID)%v0 = v0
136 StickyMap(myATID)%v0p = v0p
137 StickyMap(myATID)%rl = rl
138 StickyMap(myATID)%ru = ru
139 StickyMap(myATID)%rlp = rlp
140 StickyMap(myATID)%rup = rup
141
142 if (StickyMap(myATID)%ru .gt. StickyMap(myATID)%rup) then
143 StickyMap(myATID)%rbig = StickyMap(myATID)%ru
144 else
145 StickyMap(myATID)%rbig = StickyMap(myATID)%rup
146 endif
147
148 return
149 end subroutine newStickyType
150
151 subroutine do_sticky_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
152 pot, A, f, t, do_pot)
153
154 !! This routine does only the sticky portion of the SSD potential
155 !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
156 !! The Lennard-Jones and dipolar interaction must be handled separately.
157
158 !! We assume that the rotation matrices have already been calculated
159 !! and placed in the A array.
160
161 !! i and j are pointers to the two SSD atoms
162
163 integer, intent(in) :: atom1, atom2
164 real (kind=dp), intent(inout) :: rij, r2
165 real (kind=dp), dimension(3), intent(in) :: d
166 real (kind=dp), dimension(3), intent(inout) :: fpair
167 real (kind=dp) :: pot, vpair, sw
168 real (kind=dp), dimension(9,nLocal) :: A
169 real (kind=dp), dimension(3,nLocal) :: f
170 real (kind=dp), dimension(3,nLocal) :: t
171 logical, intent(in) :: do_pot
172
173 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
174 real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
175 real (kind=dp) :: wi, wj, w, wip, wjp, wp
176 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
177 real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
178 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
179 real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
180 real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
181 real (kind=dp) :: drdx, drdy, drdz
182 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
183 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
184 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
185 real (kind=dp) :: fxradial, fyradial, fzradial
186 real (kind=dp) :: rijtest, rjitest
187 real (kind=dp) :: radcomxi, radcomyi, radcomzi
188 real (kind=dp) :: radcomxj, radcomyj, radcomzj
189 integer :: id1, id2
190 integer :: me1, me2
191 real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig
192
193 if (.not.allocated(StickyMap)) then
194 call handleError("sticky", "no StickyMap was present before first call of do_sticky_pair!")
195 return
196 end if
197
198 #ifdef IS_MPI
199 me1 = atid_Row(atom1)
200 me2 = atid_Col(atom2)
201 #else
202 me1 = atid(atom1)
203 me2 = atid(atom2)
204 #endif
205
206 if (me1.eq.me2) then
207 w0 = StickyMap(me1)%w0
208 v0 = StickyMap(me1)%v0
209 v0p = StickyMap(me1)%v0p
210 rl = StickyMap(me1)%rl
211 ru = StickyMap(me1)%ru
212 rlp = StickyMap(me1)%rlp
213 rup = StickyMap(me1)%rup
214 rbig = StickyMap(me1)%rbig
215 else
216 ! This is silly, but if you want 2 sticky types in your
217 ! simulation, we'll let you do it with the Lorentz-
218 ! Berthelot mixing rules.
219 ! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW)
220 rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl )
221 ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru )
222 rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp )
223 rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup )
224 rbig = max(ru, rup)
225 w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 )
226 v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 )
227 v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p )
228 endif
229
230 if ( rij .LE. rbig ) then
231
232 r3 = r2*rij
233 r5 = r3*r2
234
235 drdx = d(1) / rij
236 drdy = d(2) / rij
237 drdz = d(3) / rij
238
239 #ifdef IS_MPI
240 ! rotate the inter-particle separation into the two different
241 ! body-fixed coordinate systems:
242
243 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
244 yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
245 zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
246
247 ! negative sign because this is the vector from j to i:
248
249 xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
250 yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
251 zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
252 #else
253 ! rotate the inter-particle separation into the two different
254 ! body-fixed coordinate systems:
255
256 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
257 yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
258 zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
259
260 ! negative sign because this is the vector from j to i:
261
262 xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
263 yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
264 zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
265 #endif
266
267 xi2 = xi*xi
268 yi2 = yi*yi
269 zi2 = zi*zi
270
271 xj2 = xj*xj
272 yj2 = yj*yj
273 zj2 = zj*zj
274
275 call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
276
277 wi = 2.0d0*(xi2-yi2)*zi / r3
278 wj = 2.0d0*(xj2-yj2)*zj / r3
279 w = wi+wj
280
281 zif = zi/rij - 0.6d0
282 zis = zi/rij + 0.8d0
283
284 zjf = zj/rij - 0.6d0
285 zjs = zj/rij + 0.8d0
286
287 wip = zif*zif*zis*zis - w0
288 wjp = zjf*zjf*zjs*zjs - w0
289 wp = wip + wjp
290
291 vpair = vpair + 0.5d0*(v0*s*w + v0p*sp*wp)
292 if (do_pot) then
293 #ifdef IS_MPI
294 pot_row(atom1) = pot_row(atom1) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw
295 pot_col(atom2) = pot_col(atom2) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw
296 #else
297 pot = pot + 0.5d0*(v0*s*w + v0p*sp*wp)*sw
298 #endif
299 endif
300
301 dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
302 dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
303 dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
304
305 dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
306 dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
307 dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
308
309 uglyi = zif*zif*zis + zif*zis*zis
310 uglyj = zjf*zjf*zjs + zjf*zjs*zjs
311
312 dwipdx = -2.0d0*xi*zi*uglyi/r3
313 dwipdy = -2.0d0*yi*zi*uglyi/r3
314 dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
315
316 dwjpdx = -2.0d0*xj*zj*uglyj/r3
317 dwjpdy = -2.0d0*yj*zj*uglyj/r3
318 dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
319
320 dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
321 dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
322 dwiduz = - 8.0d0*xi*yi*zi/r3
323
324 dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
325 dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
326 dwjduz = - 8.0d0*xj*yj*zj/r3
327
328 dwipdux = 2.0d0*yi*uglyi/rij
329 dwipduy = -2.0d0*xi*uglyi/rij
330 dwipduz = 0.0d0
331
332 dwjpdux = 2.0d0*yj*uglyj/rij
333 dwjpduy = -2.0d0*xj*uglyj/rij
334 dwjpduz = 0.0d0
335
336 ! do the torques first since they are easy:
337 ! remember that these are still in the body fixed axes
338
339 txi = 0.5d0*(v0*s*dwidux + v0p*sp*dwipdux)*sw
340 tyi = 0.5d0*(v0*s*dwiduy + v0p*sp*dwipduy)*sw
341 tzi = 0.5d0*(v0*s*dwiduz + v0p*sp*dwipduz)*sw
342
343 txj = 0.5d0*(v0*s*dwjdux + v0p*sp*dwjpdux)*sw
344 tyj = 0.5d0*(v0*s*dwjduy + v0p*sp*dwjpduy)*sw
345 tzj = 0.5d0*(v0*s*dwjduz + v0p*sp*dwjpduz)*sw
346
347 ! go back to lab frame using transpose of rotation matrix:
348
349 #ifdef IS_MPI
350 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
351 a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
352 t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
353 a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
354 t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
355 a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
356
357 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
358 a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
359 t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
360 a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
361 t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
362 a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
363 #else
364 t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
365 t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
366 t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
367
368 t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
369 t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
370 t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
371 #endif
372 ! Now, on to the forces:
373
374 ! first rotate the i terms back into the lab frame:
375
376 radcomxi = (v0*s*dwidx+v0p*sp*dwipdx)*sw
377 radcomyi = (v0*s*dwidy+v0p*sp*dwipdy)*sw
378 radcomzi = (v0*s*dwidz+v0p*sp*dwipdz)*sw
379
380 radcomxj = (v0*s*dwjdx+v0p*sp*dwjpdx)*sw
381 radcomyj = (v0*s*dwjdy+v0p*sp*dwjpdy)*sw
382 radcomzj = (v0*s*dwjdz+v0p*sp*dwjpdz)*sw
383
384 #ifdef IS_MPI
385 fxii = a_Row(1,atom1)*(radcomxi) + &
386 a_Row(4,atom1)*(radcomyi) + &
387 a_Row(7,atom1)*(radcomzi)
388 fyii = a_Row(2,atom1)*(radcomxi) + &
389 a_Row(5,atom1)*(radcomyi) + &
390 a_Row(8,atom1)*(radcomzi)
391 fzii = a_Row(3,atom1)*(radcomxi) + &
392 a_Row(6,atom1)*(radcomyi) + &
393 a_Row(9,atom1)*(radcomzi)
394
395 fxjj = a_Col(1,atom2)*(radcomxj) + &
396 a_Col(4,atom2)*(radcomyj) + &
397 a_Col(7,atom2)*(radcomzj)
398 fyjj = a_Col(2,atom2)*(radcomxj) + &
399 a_Col(5,atom2)*(radcomyj) + &
400 a_Col(8,atom2)*(radcomzj)
401 fzjj = a_Col(3,atom2)*(radcomxj)+ &
402 a_Col(6,atom2)*(radcomyj) + &
403 a_Col(9,atom2)*(radcomzj)
404 #else
405 fxii = a(1,atom1)*(radcomxi) + &
406 a(4,atom1)*(radcomyi) + &
407 a(7,atom1)*(radcomzi)
408 fyii = a(2,atom1)*(radcomxi) + &
409 a(5,atom1)*(radcomyi) + &
410 a(8,atom1)*(radcomzi)
411 fzii = a(3,atom1)*(radcomxi) + &
412 a(6,atom1)*(radcomyi) + &
413 a(9,atom1)*(radcomzi)
414
415 fxjj = a(1,atom2)*(radcomxj) + &
416 a(4,atom2)*(radcomyj) + &
417 a(7,atom2)*(radcomzj)
418 fyjj = a(2,atom2)*(radcomxj) + &
419 a(5,atom2)*(radcomyj) + &
420 a(8,atom2)*(radcomzj)
421 fzjj = a(3,atom2)*(radcomxj)+ &
422 a(6,atom2)*(radcomyj) + &
423 a(9,atom2)*(radcomzj)
424 #endif
425
426 fxij = -fxii
427 fyij = -fyii
428 fzij = -fzii
429
430 fxji = -fxjj
431 fyji = -fyjj
432 fzji = -fzjj
433
434 ! now assemble these with the radial-only terms:
435
436 fxradial = 0.5d0*(v0*dsdr*drdx*w + v0p*dspdr*drdx*wp + fxii + fxji)
437 fyradial = 0.5d0*(v0*dsdr*drdy*w + v0p*dspdr*drdy*wp + fyii + fyji)
438 fzradial = 0.5d0*(v0*dsdr*drdz*w + v0p*dspdr*drdz*wp + fzii + fzji)
439
440 #ifdef IS_MPI
441 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
442 f_Row(2,atom1) = f_Row(2,atom1) + fyradial
443 f_Row(3,atom1) = f_Row(3,atom1) + fzradial
444
445 f_Col(1,atom2) = f_Col(1,atom2) - fxradial
446 f_Col(2,atom2) = f_Col(2,atom2) - fyradial
447 f_Col(3,atom2) = f_Col(3,atom2) - fzradial
448 #else
449 f(1,atom1) = f(1,atom1) + fxradial
450 f(2,atom1) = f(2,atom1) + fyradial
451 f(3,atom1) = f(3,atom1) + fzradial
452
453 f(1,atom2) = f(1,atom2) - fxradial
454 f(2,atom2) = f(2,atom2) - fyradial
455 f(3,atom2) = f(3,atom2) - fzradial
456 #endif
457
458 #ifdef IS_MPI
459 id1 = AtomRowToGlobal(atom1)
460 id2 = AtomColToGlobal(atom2)
461 #else
462 id1 = atom1
463 id2 = atom2
464 #endif
465
466 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
467
468 fpair(1) = fpair(1) + fxradial
469 fpair(2) = fpair(2) + fyradial
470 fpair(3) = fpair(3) + fzradial
471
472 endif
473 endif
474 end subroutine do_sticky_pair
475
476 !! calculates the switching functions and their derivatives for a given
477 subroutine calc_sw_fnc(r, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
478
479 real (kind=dp), intent(in) :: r, rl, ru, rlp, rup
480 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
481
482 ! distances must be in angstroms
483
484 if (r.lt.rl) then
485 s = 1.0d0
486 dsdr = 0.0d0
487 elseif (r.gt.ru) then
488 s = 0.0d0
489 dsdr = 0.0d0
490 else
491 s = ((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) / &
492 ((ru - rl)**3)
493 dsdr = 6.0d0*(r-ru)*(r-rl)/((ru - rl)**3)
494 endif
495
496 if (r.lt.rlp) then
497 sp = 1.0d0
498 dspdr = 0.0d0
499 elseif (r.gt.rup) then
500 sp = 0.0d0
501 dspdr = 0.0d0
502 else
503 sp = ((rup + 2.0d0*r - 3.0d0*rlp) * (rup-r)**2) / &
504 ((rup - rlp)**3)
505 dspdr = 6.0d0*(r-rup)*(r-rlp)/((rup - rlp)**3)
506 endif
507
508 return
509 end subroutine calc_sw_fnc
510
511 subroutine destroyStickyTypes()
512 if(allocated(StickyMap)) deallocate(StickyMap)
513 end subroutine destroyStickyTypes
514
515 subroutine do_sticky_power_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
516 pot, A, f, t, do_pot)
517 !! We assume that the rotation matrices have already been calculated
518 !! and placed in the A array.
519
520 !! i and j are pointers to the two SSD atoms
521
522 integer, intent(in) :: atom1, atom2
523 real (kind=dp), intent(inout) :: rij, r2
524 real (kind=dp), dimension(3), intent(in) :: d
525 real (kind=dp), dimension(3), intent(inout) :: fpair
526 real (kind=dp) :: pot, vpair, sw
527 real (kind=dp), dimension(9,nLocal) :: A
528 real (kind=dp), dimension(3,nLocal) :: f
529 real (kind=dp), dimension(3,nLocal) :: t
530 logical, intent(in) :: do_pot
531
532 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
533 real (kind=dp) :: xihat, yihat, zihat, xjhat, yjhat, zjhat
534 real (kind=dp) :: rI, rI2, rI3, rI4, rI5, rI6, rI7, s, sp, dsdr, dspdr
535 real (kind=dp) :: wi, wj, w, wi2, wj2, eScale, v0scale
536 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
537 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
538 real (kind=dp) :: drdx, drdy, drdz
539 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
540 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
541 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
542 real (kind=dp) :: fxradial, fyradial, fzradial
543 real (kind=dp) :: rijtest, rjitest
544 real (kind=dp) :: radcomxi, radcomyi, radcomzi
545 real (kind=dp) :: radcomxj, radcomyj, radcomzj
546 integer :: id1, id2
547 integer :: me1, me2
548 real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig
549 real (kind=dp) :: zi3, zi4, zi5, zj3, zj4, zj5
550 real (kind=dp) :: frac1, frac2
551
552 if (.not.allocated(StickyMap)) then
553 call handleError("sticky", "no StickyMap was present before first call of do_sticky_power_pair!")
554 return
555 end if
556
557 #ifdef IS_MPI
558 me1 = atid_Row(atom1)
559 me2 = atid_Col(atom2)
560 #else
561 me1 = atid(atom1)
562 me2 = atid(atom2)
563 #endif
564
565 if (me1.eq.me2) then
566 w0 = StickyMap(me1)%w0
567 v0 = StickyMap(me1)%v0
568 v0p = StickyMap(me1)%v0p
569 rl = StickyMap(me1)%rl
570 ru = StickyMap(me1)%ru
571 rlp = StickyMap(me1)%rlp
572 rup = StickyMap(me1)%rup
573 rbig = StickyMap(me1)%rbig
574 else
575 ! This is silly, but if you want 2 sticky types in your
576 ! simulation, we'll let you do it with the Lorentz-
577 ! Berthelot mixing rules.
578 ! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW)
579 rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl )
580 ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru )
581 rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp )
582 rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup )
583 rbig = max(ru, rup)
584 w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 )
585 v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 )
586 v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p )
587 endif
588
589 if ( rij .LE. rbig ) then
590
591 rI = 1.0d0/rij
592 rI2 = rI*rI
593 rI3 = rI2*rI
594 rI4 = rI2*rI2
595 rI5 = rI3*rI2
596 rI6 = rI3*rI3
597 rI7 = rI4*rI3
598
599 drdx = d(1) * rI
600 drdy = d(2) * rI
601 drdz = d(3) * rI
602
603 #ifdef IS_MPI
604 ! rotate the inter-particle separation into the two different
605 ! body-fixed coordinate systems:
606
607 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
608 yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
609 zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
610
611 ! negative sign because this is the vector from j to i:
612
613 xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
614 yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
615 zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
616 #else
617 ! rotate the inter-particle separation into the two different
618 ! body-fixed coordinate systems:
619
620 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
621 yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
622 zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
623
624 ! negative sign because this is the vector from j to i:
625
626 xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
627 yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
628 zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
629 #endif
630
631 xi2 = xi*xi
632 yi2 = yi*yi
633 zi2 = zi*zi
634 zi3 = zi2*zi
635 zi4 = zi2*zi2
636 zi5 = zi3*zi2
637 xihat = xi*rI
638 yihat = yi*rI
639 zihat = zi*rI
640
641 xj2 = xj*xj
642 yj2 = yj*yj
643 zj2 = zj*zj
644 zj3 = zj2*zj
645 zj4 = zj2*zj2
646 zj5 = zj3*zj2
647 xjhat = xj*rI
648 yjhat = yj*rI
649 zjhat = zj*rI
650
651 call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
652
653 frac1 = 0.25d0
654 frac2 = 0.75d0
655
656 wi = 2.0d0*(xi2-yi2)*zi*rI3
657 wj = 2.0d0*(xj2-yj2)*zj*rI3
658
659 wi2 = wi*wi
660 wj2 = wj*wj
661
662 w = frac1*wi*wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p
663
664 vpair = vpair + 0.5d0*(v0*s*w)
665
666 if (do_pot) then
667 #ifdef IS_MPI
668 pot_row(atom1) = pot_row(atom1) + 0.25d0*(v0*s*w)*sw
669 pot_col(atom2) = pot_col(atom2) + 0.25d0*(v0*s*w)*sw
670 #else
671 pot = pot + 0.5d0*(v0*s*w)*sw
672 #endif
673 endif
674
675 dwidx = ( 4.0d0*xi*zi*rI3 - 6.0d0*xi*zi*(xi2-yi2)*rI5 )
676 dwidy = ( -4.0d0*yi*zi*rI3 - 6.0d0*yi*zi*(xi2-yi2)*rI5 )
677 dwidz = ( 2.0d0*(xi2-yi2)*rI3 - 6.0d0*zi2*(xi2-yi2)*rI5 )
678
679 dwidx = frac1*3.0d0*wi2*dwidx + frac2*dwidx
680 dwidy = frac1*3.0d0*wi2*dwidy + frac2*dwidy
681 dwidz = frac1*3.0d0*wi2*dwidz + frac2*dwidz
682
683 dwjdx = ( 4.0d0*xj*zj*rI3 - 6.0d0*xj*zj*(xj2-yj2)*rI5 )
684 dwjdy = ( -4.0d0*yj*zj*rI3 - 6.0d0*yj*zj*(xj2-yj2)*rI5 )
685 dwjdz = ( 2.0d0*(xj2-yj2)*rI3 - 6.0d0*zj2*(xj2-yj2)*rI5 )
686
687 dwjdx = frac1*3.0d0*wj2*dwjdx + frac2*dwjdx
688 dwjdy = frac1*3.0d0*wj2*dwjdy + frac2*dwjdy
689 dwjdz = frac1*3.0d0*wj2*dwjdz + frac2*dwjdz
690
691 dwidux = ( 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))*rI3 )
692 dwiduy = ( 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))*rI3 )
693 dwiduz = ( -8.0d0*xi*yi*zi*rI3 )
694
695 dwidux = frac1*3.0d0*wi2*dwidux + frac2*dwidux
696 dwiduy = frac1*3.0d0*wi2*dwiduy + frac2*dwiduy
697 dwiduz = frac1*3.0d0*wi2*dwiduz + frac2*dwiduz
698
699 dwjdux = ( 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))*rI3 )
700 dwjduy = ( 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))*rI3 )
701 dwjduz = ( -8.0d0*xj*yj*zj*rI3 )
702
703 dwjdux = frac1*3.0d0*wj2*dwjdux + frac2*dwjdux
704 dwjduy = frac1*3.0d0*wj2*dwjduy + frac2*dwjduy
705 dwjduz = frac1*3.0d0*wj2*dwjduz + frac2*dwjduz
706
707 ! do the torques first since they are easy:
708 ! remember that these are still in the body fixed axes
709
710 txi = 0.5d0*(v0*s*dwidux)*sw
711 tyi = 0.5d0*(v0*s*dwiduy)*sw
712 tzi = 0.5d0*(v0*s*dwiduz)*sw
713
714 txj = 0.5d0*(v0*s*dwjdux)*sw
715 tyj = 0.5d0*(v0*s*dwjduy)*sw
716 tzj = 0.5d0*(v0*s*dwjduz)*sw
717
718 ! go back to lab frame using transpose of rotation matrix:
719
720 #ifdef IS_MPI
721 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
722 a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
723 t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
724 a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
725 t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
726 a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
727
728 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
729 a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
730 t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
731 a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
732 t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
733 a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
734 #else
735 t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
736 t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
737 t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
738
739 t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
740 t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
741 t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
742 #endif
743 ! Now, on to the forces:
744
745 ! first rotate the i terms back into the lab frame:
746
747 radcomxi = (v0*s*dwidx)*sw
748 radcomyi = (v0*s*dwidy)*sw
749 radcomzi = (v0*s*dwidz)*sw
750
751 radcomxj = (v0*s*dwjdx)*sw
752 radcomyj = (v0*s*dwjdy)*sw
753 radcomzj = (v0*s*dwjdz)*sw
754
755 #ifdef IS_MPI
756 fxii = a_Row(1,atom1)*(radcomxi) + &
757 a_Row(4,atom1)*(radcomyi) + &
758 a_Row(7,atom1)*(radcomzi)
759 fyii = a_Row(2,atom1)*(radcomxi) + &
760 a_Row(5,atom1)*(radcomyi) + &
761 a_Row(8,atom1)*(radcomzi)
762 fzii = a_Row(3,atom1)*(radcomxi) + &
763 a_Row(6,atom1)*(radcomyi) + &
764 a_Row(9,atom1)*(radcomzi)
765
766 fxjj = a_Col(1,atom2)*(radcomxj) + &
767 a_Col(4,atom2)*(radcomyj) + &
768 a_Col(7,atom2)*(radcomzj)
769 fyjj = a_Col(2,atom2)*(radcomxj) + &
770 a_Col(5,atom2)*(radcomyj) + &
771 a_Col(8,atom2)*(radcomzj)
772 fzjj = a_Col(3,atom2)*(radcomxj)+ &
773 a_Col(6,atom2)*(radcomyj) + &
774 a_Col(9,atom2)*(radcomzj)
775 #else
776 fxii = a(1,atom1)*(radcomxi) + &
777 a(4,atom1)*(radcomyi) + &
778 a(7,atom1)*(radcomzi)
779 fyii = a(2,atom1)*(radcomxi) + &
780 a(5,atom1)*(radcomyi) + &
781 a(8,atom1)*(radcomzi)
782 fzii = a(3,atom1)*(radcomxi) + &
783 a(6,atom1)*(radcomyi) + &
784 a(9,atom1)*(radcomzi)
785
786 fxjj = a(1,atom2)*(radcomxj) + &
787 a(4,atom2)*(radcomyj) + &
788 a(7,atom2)*(radcomzj)
789 fyjj = a(2,atom2)*(radcomxj) + &
790 a(5,atom2)*(radcomyj) + &
791 a(8,atom2)*(radcomzj)
792 fzjj = a(3,atom2)*(radcomxj)+ &
793 a(6,atom2)*(radcomyj) + &
794 a(9,atom2)*(radcomzj)
795 #endif
796
797 fxij = -fxii
798 fyij = -fyii
799 fzij = -fzii
800
801 fxji = -fxjj
802 fyji = -fyjj
803 fzji = -fzjj
804
805 ! now assemble these with the radial-only terms:
806
807 fxradial = 0.5d0*(v0*dsdr*w*drdx + fxii + fxji)
808 fyradial = 0.5d0*(v0*dsdr*w*drdy + fyii + fyji)
809 fzradial = 0.5d0*(v0*dsdr*w*drdz + fzii + fzji)
810
811 #ifdef IS_MPI
812 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
813 f_Row(2,atom1) = f_Row(2,atom1) + fyradial
814 f_Row(3,atom1) = f_Row(3,atom1) + fzradial
815
816 f_Col(1,atom2) = f_Col(1,atom2) - fxradial
817 f_Col(2,atom2) = f_Col(2,atom2) - fyradial
818 f_Col(3,atom2) = f_Col(3,atom2) - fzradial
819 #else
820 f(1,atom1) = f(1,atom1) + fxradial
821 f(2,atom1) = f(2,atom1) + fyradial
822 f(3,atom1) = f(3,atom1) + fzradial
823
824 f(1,atom2) = f(1,atom2) - fxradial
825 f(2,atom2) = f(2,atom2) - fyradial
826 f(3,atom2) = f(3,atom2) - fzradial
827 #endif
828
829 #ifdef IS_MPI
830 id1 = AtomRowToGlobal(atom1)
831 id2 = AtomColToGlobal(atom2)
832 #else
833 id1 = atom1
834 id2 = atom2
835 #endif
836
837 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
838
839 fpair(1) = fpair(1) + fxradial
840 fpair(2) = fpair(2) + fyradial
841 fpair(3) = fpair(3) + fzradial
842
843 endif
844 endif
845 end subroutine do_sticky_power_pair
846
847 end module sticky