ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/sticky.F90
Revision: 2277
Committed: Fri Aug 26 21:30:41 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 29302 byte(s)
Log Message:
added some probably nonfunctional get*cut routines

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