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.5 2005-03-12 19:05:16 chrisfen Exp $, $Date: 2005-03-12 19:05:16 $, $Name: not supported by cvs2svn $, $Revision: 1.5 $ |
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 |
|
73 |
|
74 |
type :: StickyList |
75 |
integer :: c_ident |
76 |
real( kind = dp ) :: w0 = 0.0_dp |
77 |
real( kind = dp ) :: v0 = 0.0_dp |
78 |
real( kind = dp ) :: v0p = 0.0_dp |
79 |
real( kind = dp ) :: rl = 0.0_dp |
80 |
real( kind = dp ) :: ru = 0.0_dp |
81 |
real( kind = dp ) :: rlp = 0.0_dp |
82 |
real( kind = dp ) :: rup = 0.0_dp |
83 |
real( kind = dp ) :: rbig = 0.0_dp |
84 |
end type StickyList |
85 |
|
86 |
type(StickyList), dimension(:),allocatable :: StickyMap |
87 |
|
88 |
contains |
89 |
|
90 |
subroutine newStickyType(c_ident, w0, v0, v0p, rl, ru, rlp, rup, isError) |
91 |
|
92 |
integer, intent(in) :: c_ident |
93 |
integer, intent(inout) :: isError |
94 |
real( kind = dp ), intent(in) :: w0, v0, v0p |
95 |
real( kind = dp ), intent(in) :: rl, ru |
96 |
real( kind = dp ), intent(in) :: rlp, rup |
97 |
integer :: nATypes, myATID |
98 |
|
99 |
|
100 |
isError = 0 |
101 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
102 |
|
103 |
!! Be simple-minded and assume that we need a StickyMap that |
104 |
!! is the same size as the total number of atom types |
105 |
|
106 |
if (.not.allocated(StickyMap)) then |
107 |
|
108 |
nAtypes = getSize(atypes) |
109 |
|
110 |
if (nAtypes == 0) then |
111 |
isError = -1 |
112 |
return |
113 |
end if |
114 |
|
115 |
if (.not. allocated(StickyMap)) then |
116 |
allocate(StickyMap(nAtypes)) |
117 |
endif |
118 |
|
119 |
end if |
120 |
|
121 |
if (myATID .gt. size(StickyMap)) then |
122 |
isError = -1 |
123 |
return |
124 |
endif |
125 |
|
126 |
! set the values for StickyMap for this atom type: |
127 |
|
128 |
StickyMap(myATID)%c_ident = c_ident |
129 |
|
130 |
! we could pass all 5 parameters if we felt like it... |
131 |
|
132 |
StickyMap(myATID)%w0 = w0 |
133 |
StickyMap(myATID)%v0 = v0 |
134 |
StickyMap(myATID)%v0p = v0p |
135 |
StickyMap(myATID)%rl = rl |
136 |
StickyMap(myATID)%ru = ru |
137 |
StickyMap(myATID)%rlp = rlp |
138 |
StickyMap(myATID)%rup = rup |
139 |
|
140 |
if (StickyMap(myATID)%ru .gt. StickyMap(myATID)%rup) then |
141 |
StickyMap(myATID)%rbig = StickyMap(myATID)%ru |
142 |
else |
143 |
StickyMap(myATID)%rbig = StickyMap(myATID)%rup |
144 |
endif |
145 |
|
146 |
return |
147 |
end subroutine newStickyType |
148 |
|
149 |
subroutine do_sticky_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
150 |
pot, A, f, t, do_pot) |
151 |
|
152 |
!! This routine does only the sticky portion of the SSD potential |
153 |
!! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)]. |
154 |
!! The Lennard-Jones and dipolar interaction must be handled separately. |
155 |
|
156 |
!! We assume that the rotation matrices have already been calculated |
157 |
!! and placed in the A array. |
158 |
|
159 |
!! i and j are pointers to the two SSD atoms |
160 |
|
161 |
integer, intent(in) :: atom1, atom2 |
162 |
real (kind=dp), intent(inout) :: rij, r2 |
163 |
real (kind=dp), dimension(3), intent(in) :: d |
164 |
real (kind=dp), dimension(3), intent(inout) :: fpair |
165 |
real (kind=dp) :: pot, vpair, sw |
166 |
real (kind=dp), dimension(9,nLocal) :: A |
167 |
real (kind=dp), dimension(3,nLocal) :: f |
168 |
real (kind=dp), dimension(3,nLocal) :: t |
169 |
logical, intent(in) :: do_pot |
170 |
|
171 |
real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 |
172 |
real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr |
173 |
real (kind=dp) :: wi, wj, w, wip, wjp, wp |
174 |
real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz |
175 |
real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz |
176 |
real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz |
177 |
real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz |
178 |
real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj |
179 |
real (kind=dp) :: drdx, drdy, drdz |
180 |
real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj |
181 |
real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj |
182 |
real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji |
183 |
real (kind=dp) :: fxradial, fyradial, fzradial |
184 |
real (kind=dp) :: rijtest, rjitest |
185 |
real (kind=dp) :: radcomxi, radcomyi, radcomzi |
186 |
real (kind=dp) :: radcomxj, radcomyj, radcomzj |
187 |
integer :: id1, id2 |
188 |
integer :: me1, me2 |
189 |
real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig |
190 |
|
191 |
if (.not.allocated(StickyMap)) then |
192 |
call handleError("sticky", "no StickyMap was present before first call of do_sticky_pair!") |
193 |
return |
194 |
end if |
195 |
|
196 |
#ifdef IS_MPI |
197 |
me1 = atid_Row(atom1) |
198 |
me2 = atid_Col(atom2) |
199 |
#else |
200 |
me1 = atid(atom1) |
201 |
me2 = atid(atom2) |
202 |
#endif |
203 |
|
204 |
if (me1.eq.me2) then |
205 |
w0 = StickyMap(me1)%w0 |
206 |
v0 = StickyMap(me1)%v0 |
207 |
v0p = StickyMap(me1)%v0p |
208 |
rl = StickyMap(me1)%rl |
209 |
ru = StickyMap(me1)%ru |
210 |
rlp = StickyMap(me1)%rlp |
211 |
rup = StickyMap(me1)%rup |
212 |
rbig = StickyMap(me1)%rbig |
213 |
else |
214 |
! This is silly, but if you want 2 sticky types in your |
215 |
! simulation, we'll let you do it with the Lorentz- |
216 |
! Berthelot mixing rules. |
217 |
! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW) |
218 |
rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl ) |
219 |
ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru ) |
220 |
rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp ) |
221 |
rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup ) |
222 |
rbig = max(ru, rup) |
223 |
w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 ) |
224 |
v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 ) |
225 |
v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p ) |
226 |
endif |
227 |
|
228 |
if ( rij .LE. rbig ) then |
229 |
|
230 |
r3 = r2*rij |
231 |
r5 = r3*r2 |
232 |
|
233 |
drdx = d(1) / rij |
234 |
drdy = d(2) / rij |
235 |
drdz = d(3) / rij |
236 |
|
237 |
#ifdef IS_MPI |
238 |
! rotate the inter-particle separation into the two different |
239 |
! body-fixed coordinate systems: |
240 |
|
241 |
xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3) |
242 |
yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3) |
243 |
zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3) |
244 |
|
245 |
! negative sign because this is the vector from j to i: |
246 |
|
247 |
xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3)) |
248 |
yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3)) |
249 |
zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3)) |
250 |
#else |
251 |
! rotate the inter-particle separation into the two different |
252 |
! body-fixed coordinate systems: |
253 |
|
254 |
xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3) |
255 |
yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3) |
256 |
zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3) |
257 |
|
258 |
! negative sign because this is the vector from j to i: |
259 |
|
260 |
xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3)) |
261 |
yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3)) |
262 |
zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3)) |
263 |
#endif |
264 |
|
265 |
xi2 = xi*xi |
266 |
yi2 = yi*yi |
267 |
zi2 = zi*zi |
268 |
|
269 |
xj2 = xj*xj |
270 |
yj2 = yj*yj |
271 |
zj2 = zj*zj |
272 |
|
273 |
call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr) |
274 |
|
275 |
wi = 2.0d0*(xi2-yi2)*zi / r3 |
276 |
wj = 2.0d0*(xj2-yj2)*zj / r3 |
277 |
w = wi+wj |
278 |
|
279 |
zif = zi/rij - 0.6d0 |
280 |
zis = zi/rij + 0.8d0 |
281 |
|
282 |
zjf = zj/rij - 0.6d0 |
283 |
zjs = zj/rij + 0.8d0 |
284 |
|
285 |
wip = zif*zif*zis*zis - w0 |
286 |
wjp = zjf*zjf*zjs*zjs - w0 |
287 |
wp = wip + wjp |
288 |
|
289 |
vpair = vpair + 0.5d0*(v0*s*w + v0p*sp*wp) |
290 |
if (do_pot) then |
291 |
#ifdef IS_MPI |
292 |
pot_row(atom1) = pot_row(atom1) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw |
293 |
pot_col(atom2) = pot_col(atom2) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw |
294 |
#else |
295 |
pot = pot + 0.5d0*(v0*s*w + v0p*sp*wp)*sw |
296 |
#endif |
297 |
endif |
298 |
|
299 |
dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5 |
300 |
dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5 |
301 |
dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5 |
302 |
|
303 |
dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5 |
304 |
dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5 |
305 |
dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5 |
306 |
|
307 |
uglyi = zif*zif*zis + zif*zis*zis |
308 |
uglyj = zjf*zjf*zjs + zjf*zjs*zjs |
309 |
|
310 |
dwipdx = -2.0d0*xi*zi*uglyi/r3 |
311 |
dwipdy = -2.0d0*yi*zi*uglyi/r3 |
312 |
dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi |
313 |
|
314 |
dwjpdx = -2.0d0*xj*zj*uglyj/r3 |
315 |
dwjpdy = -2.0d0*yj*zj*uglyj/r3 |
316 |
dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj |
317 |
|
318 |
dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3 |
319 |
dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3 |
320 |
dwiduz = - 8.0d0*xi*yi*zi/r3 |
321 |
|
322 |
dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3 |
323 |
dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3 |
324 |
dwjduz = - 8.0d0*xj*yj*zj/r3 |
325 |
|
326 |
dwipdux = 2.0d0*yi*uglyi/rij |
327 |
dwipduy = -2.0d0*xi*uglyi/rij |
328 |
dwipduz = 0.0d0 |
329 |
|
330 |
dwjpdux = 2.0d0*yj*uglyj/rij |
331 |
dwjpduy = -2.0d0*xj*uglyj/rij |
332 |
dwjpduz = 0.0d0 |
333 |
|
334 |
! do the torques first since they are easy: |
335 |
! remember that these are still in the body fixed axes |
336 |
|
337 |
txi = 0.5d0*(v0*s*dwidux + v0p*sp*dwipdux)*sw |
338 |
tyi = 0.5d0*(v0*s*dwiduy + v0p*sp*dwipduy)*sw |
339 |
tzi = 0.5d0*(v0*s*dwiduz + v0p*sp*dwipduz)*sw |
340 |
|
341 |
txj = 0.5d0*(v0*s*dwjdux + v0p*sp*dwjpdux)*sw |
342 |
tyj = 0.5d0*(v0*s*dwjduy + v0p*sp*dwjpduy)*sw |
343 |
tzj = 0.5d0*(v0*s*dwjduz + v0p*sp*dwjpduz)*sw |
344 |
|
345 |
! go back to lab frame using transpose of rotation matrix: |
346 |
|
347 |
#ifdef IS_MPI |
348 |
t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + & |
349 |
a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi |
350 |
t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + & |
351 |
a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi |
352 |
t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + & |
353 |
a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi |
354 |
|
355 |
t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + & |
356 |
a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj |
357 |
t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + & |
358 |
a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj |
359 |
t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + & |
360 |
a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj |
361 |
#else |
362 |
t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi |
363 |
t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi |
364 |
t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi |
365 |
|
366 |
t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj |
367 |
t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj |
368 |
t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj |
369 |
#endif |
370 |
! Now, on to the forces: |
371 |
|
372 |
! first rotate the i terms back into the lab frame: |
373 |
|
374 |
radcomxi = (v0*s*dwidx+v0p*sp*dwipdx)*sw |
375 |
radcomyi = (v0*s*dwidy+v0p*sp*dwipdy)*sw |
376 |
radcomzi = (v0*s*dwidz+v0p*sp*dwipdz)*sw |
377 |
|
378 |
radcomxj = (v0*s*dwjdx+v0p*sp*dwjpdx)*sw |
379 |
radcomyj = (v0*s*dwjdy+v0p*sp*dwjpdy)*sw |
380 |
radcomzj = (v0*s*dwjdz+v0p*sp*dwjpdz)*sw |
381 |
|
382 |
#ifdef IS_MPI |
383 |
fxii = a_Row(1,atom1)*(radcomxi) + & |
384 |
a_Row(4,atom1)*(radcomyi) + & |
385 |
a_Row(7,atom1)*(radcomzi) |
386 |
fyii = a_Row(2,atom1)*(radcomxi) + & |
387 |
a_Row(5,atom1)*(radcomyi) + & |
388 |
a_Row(8,atom1)*(radcomzi) |
389 |
fzii = a_Row(3,atom1)*(radcomxi) + & |
390 |
a_Row(6,atom1)*(radcomyi) + & |
391 |
a_Row(9,atom1)*(radcomzi) |
392 |
|
393 |
fxjj = a_Col(1,atom2)*(radcomxj) + & |
394 |
a_Col(4,atom2)*(radcomyj) + & |
395 |
a_Col(7,atom2)*(radcomzj) |
396 |
fyjj = a_Col(2,atom2)*(radcomxj) + & |
397 |
a_Col(5,atom2)*(radcomyj) + & |
398 |
a_Col(8,atom2)*(radcomzj) |
399 |
fzjj = a_Col(3,atom2)*(radcomxj)+ & |
400 |
a_Col(6,atom2)*(radcomyj) + & |
401 |
a_Col(9,atom2)*(radcomzj) |
402 |
#else |
403 |
fxii = a(1,atom1)*(radcomxi) + & |
404 |
a(4,atom1)*(radcomyi) + & |
405 |
a(7,atom1)*(radcomzi) |
406 |
fyii = a(2,atom1)*(radcomxi) + & |
407 |
a(5,atom1)*(radcomyi) + & |
408 |
a(8,atom1)*(radcomzi) |
409 |
fzii = a(3,atom1)*(radcomxi) + & |
410 |
a(6,atom1)*(radcomyi) + & |
411 |
a(9,atom1)*(radcomzi) |
412 |
|
413 |
fxjj = a(1,atom2)*(radcomxj) + & |
414 |
a(4,atom2)*(radcomyj) + & |
415 |
a(7,atom2)*(radcomzj) |
416 |
fyjj = a(2,atom2)*(radcomxj) + & |
417 |
a(5,atom2)*(radcomyj) + & |
418 |
a(8,atom2)*(radcomzj) |
419 |
fzjj = a(3,atom2)*(radcomxj)+ & |
420 |
a(6,atom2)*(radcomyj) + & |
421 |
a(9,atom2)*(radcomzj) |
422 |
#endif |
423 |
|
424 |
fxij = -fxii |
425 |
fyij = -fyii |
426 |
fzij = -fzii |
427 |
|
428 |
fxji = -fxjj |
429 |
fyji = -fyjj |
430 |
fzji = -fzjj |
431 |
|
432 |
! now assemble these with the radial-only terms: |
433 |
|
434 |
fxradial = 0.5d0*(v0*dsdr*drdx*w + v0p*dspdr*drdx*wp + fxii + fxji) |
435 |
fyradial = 0.5d0*(v0*dsdr*drdy*w + v0p*dspdr*drdy*wp + fyii + fyji) |
436 |
fzradial = 0.5d0*(v0*dsdr*drdz*w + v0p*dspdr*drdz*wp + fzii + fzji) |
437 |
|
438 |
#ifdef IS_MPI |
439 |
f_Row(1,atom1) = f_Row(1,atom1) + fxradial |
440 |
f_Row(2,atom1) = f_Row(2,atom1) + fyradial |
441 |
f_Row(3,atom1) = f_Row(3,atom1) + fzradial |
442 |
|
443 |
f_Col(1,atom2) = f_Col(1,atom2) - fxradial |
444 |
f_Col(2,atom2) = f_Col(2,atom2) - fyradial |
445 |
f_Col(3,atom2) = f_Col(3,atom2) - fzradial |
446 |
#else |
447 |
f(1,atom1) = f(1,atom1) + fxradial |
448 |
f(2,atom1) = f(2,atom1) + fyradial |
449 |
f(3,atom1) = f(3,atom1) + fzradial |
450 |
|
451 |
f(1,atom2) = f(1,atom2) - fxradial |
452 |
f(2,atom2) = f(2,atom2) - fyradial |
453 |
f(3,atom2) = f(3,atom2) - fzradial |
454 |
#endif |
455 |
|
456 |
#ifdef IS_MPI |
457 |
id1 = AtomRowToGlobal(atom1) |
458 |
id2 = AtomColToGlobal(atom2) |
459 |
#else |
460 |
id1 = atom1 |
461 |
id2 = atom2 |
462 |
#endif |
463 |
|
464 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
465 |
|
466 |
fpair(1) = fpair(1) + fxradial |
467 |
fpair(2) = fpair(2) + fyradial |
468 |
fpair(3) = fpair(3) + fzradial |
469 |
|
470 |
endif |
471 |
endif |
472 |
end subroutine do_sticky_pair |
473 |
|
474 |
!! calculates the switching functions and their derivatives for a given |
475 |
subroutine calc_sw_fnc(r, rl, ru, rlp, rup, s, sp, dsdr, dspdr) |
476 |
|
477 |
real (kind=dp), intent(in) :: r, rl, ru, rlp, rup |
478 |
real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr |
479 |
|
480 |
! distances must be in angstroms |
481 |
|
482 |
if (r.lt.rl) then |
483 |
s = 1.0d0 |
484 |
dsdr = 0.0d0 |
485 |
elseif (r.gt.ru) then |
486 |
s = 0.0d0 |
487 |
dsdr = 0.0d0 |
488 |
else |
489 |
s = ((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) / & |
490 |
((ru - rl)**3) |
491 |
dsdr = 6.0d0*(r-ru)*(r-rl)/((ru - rl)**3) |
492 |
endif |
493 |
|
494 |
if (r.lt.rlp) then |
495 |
sp = 1.0d0 |
496 |
dspdr = 0.0d0 |
497 |
elseif (r.gt.rup) then |
498 |
sp = 0.0d0 |
499 |
dspdr = 0.0d0 |
500 |
else |
501 |
sp = ((rup + 2.0d0*r - 3.0d0*rlp) * (rup-r)**2) / & |
502 |
((rup - rlp)**3) |
503 |
dspdr = 6.0d0*(r-rup)*(r-rlp)/((rup - rlp)**3) |
504 |
endif |
505 |
|
506 |
return |
507 |
end subroutine calc_sw_fnc |
508 |
end module sticky |