ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-2.0/src/UseTheForce/DarkSide/sticky.F90
Revision: 1877
Committed: Thu Dec 9 21:15:19 2004 UTC (19 years, 8 months ago) by tim
File size: 15593 byte(s)
Log Message:
sticky module get compiled

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