ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/sticky.F90
Revision: 2756
Committed: Wed May 17 15:37:15 2006 UTC (18 years, 2 months ago) by gezelter
File size: 31629 byte(s)
Log Message:
Getting fortran side prepped for single precision...

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