ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/UseTheForce/DarkSide/sticky.F90
Revision: 1874
Committed: Thu Dec 9 20:27:59 2004 UTC (19 years, 9 months ago) by gezelter
File size: 15457 byte(s)
Log Message:
sticky module now has the option for multiple sticky types

File Contents

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