ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_sticky_pair.F90
Revision: 322
Committed: Wed Mar 12 15:17:52 2003 UTC (21 years, 5 months ago) by gezelter
File size: 12450 byte(s)
Log Message:
Massive sticky rewrite

File Contents

# User Rev Content
1 gezelter 319 !! 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 gezelter 322 !! @author Christopher Fennel
11     !! @author J. Daniel Gezelter
12     !! @version $Id: calc_sticky_pair.F90,v 1.2 2003-03-12 15:17:52 gezelter Exp $, $Date: 2003-03-12 15:17:52 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $
13 gezelter 319
14 gezelter 322 module sticky_pair
15 gezelter 319
16     use simulation
17     use definitions
18 gezelter 322 use forceGlobals
19 gezelter 319 #ifdef IS_MPI
20     use mpiSimulation
21     #endif
22 gezelter 322
23 gezelter 319 implicit none
24 gezelter 322
25 gezelter 319 PRIVATE
26    
27 gezelter 322 logical, save :: sticky_initialized = .false.
28     real( kind = dp ), save :: SSD_w0
29     real( kind = dp ), save :: SSD_v0
30     real( kind = dp ), save :: SSD_rl
31     real( kind = dp ), save :: SSD_ru
32     real( kind = dp ), save :: SSD_rup
33 gezelter 319
34 gezelter 322 public :: initialize_sticky
35     public :: do_sticky_pair
36 gezelter 319
37     contains
38    
39 gezelter 322 subroutine initialize_sticky(sticky_w0, sticky_v0)
40     real( kind = dp ), intent(in) :: sticky_w0, sticky_v0
41    
42     ! we could pass all 5 parameters if we felt like it...
43 gezelter 319
44 gezelter 322 SSD_w0 = sticky_w0
45     SSD_v0 = sticky_v0
46     SSD_rl = 2.75_DP
47     SSD_ru = 3.35_DP
48     SSD_rup = 4.0_DP
49 gezelter 319
50 gezelter 322 sticky_initialized = .true.
51     return
52     end subroutine initialize_sticky
53 gezelter 319
54 gezelter 322 subroutine do_sticky_pair(atom1, atom2, d, rij, A, pot, f, t)
55    
56     !! This routine does only the sticky portion of the SSD potential
57     !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
58     !! The Lennard-Jones and dipolar interaction must be handled separately.
59    
60     !! We assume that the rotation matrices have already been calculated
61     !! and placed in the A array.
62    
63     !! i and j are pointers to the two SSD atoms
64 gezelter 319
65 gezelter 322 integer, intent(in) :: atom1, atom2
66     real (kind=dp), intent(inout) :: rij
67     real (kind=dp), dimension(3), intent(in) :: d
68     real (kind=dp), dimension(9,getNlocal()) :: A
69     real (kind=dp), dimension(3,getNlocal()) :: f
70     real (kind=dp), dimension(3,getNlocal()) :: t
71 gezelter 319
72 gezelter 322 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
73     real (kind=dp) :: r2, r3, r5, r6, s, sp, dsdr, dspdr
74     real (kind=dp) :: wi, wj, w, wip, wjp, wp
75     real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
76     real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
77     real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
78     real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
79     real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
80     real (kind=dp) :: drdx, drdy, drdz
81     real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
82     real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
83     real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
84     real (kind=dp) :: fxradial, fyradial, fzradial
85    
86     if (.not.sticky_initialized) then
87     write(*,*) 'Sticky forces not initialized!'
88 gezelter 319 return
89     endif
90    
91 gezelter 322 r2 = rij*rij
92     r3 = r2*rij
93     r5 = r3*r2
94    
95     drdx = d(1) / rij
96     drdy = d(2) / rij
97     drdz = d(3) / rij
98    
99 gezelter 319 #ifdef IS_MPI
100 gezelter 322 ! rotate the inter-particle separation into the two different
101     ! body-fixed coordinate systems:
102 gezelter 319
103 gezelter 322 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
104     yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
105     zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
106 gezelter 319
107 gezelter 322 ! negative sign because this is the vector from j to i:
108    
109     xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
110     yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
111     zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
112 gezelter 319 #else
113 gezelter 322 ! rotate the inter-particle separation into the two different
114     ! body-fixed coordinate systems:
115 gezelter 319
116 gezelter 322 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
117     yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
118     zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
119    
120     ! negative sign because this is the vector from j to i:
121    
122     xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
123     yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
124     zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
125 gezelter 319 #endif
126    
127 gezelter 322 xi2 = xi*xi
128     yi2 = yi*yi
129     zi2 = zi*zi
130 gezelter 319
131 gezelter 322 xj2 = xj*xj
132     yj2 = yj*yj
133     zj2 = zj*zj
134 gezelter 319
135 gezelter 322 call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
136    
137     wi = 2.0d0*(xi2-yi2)*zi / r3
138     wj = 2.0d0*(xj2-yj2)*zj / r3
139     w = wi+wj
140    
141     zif = zi/rij - 0.6d0
142     zis = zi/rij + 0.8d0
143    
144     zjf = zj/rij - 0.6d0
145     zjs = zj/rij + 0.8d0
146    
147     wip = zif*zif*zis*zis - SSD_w0
148     wjp = zjf*zjf*zjs*zjs - SSD_w0
149     wp = wip + wjp
150 gezelter 319
151 gezelter 322 #ifdef IS_MPI
152     pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
153     pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
154     #else
155     pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
156     #endif
157 gezelter 319
158 gezelter 322 dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
159     dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
160     dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
161    
162     dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
163     dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
164     dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
165    
166     uglyi = zif*zif*zis + zif*zis*zis
167     uglyj = zjf*zjf*zjs + zjf*zjs*zjs
168    
169     dwipdx = -2.0d0*xi*zi*uglyi/r3
170     dwipdy = -2.0d0*yi*zi*uglyi/r3
171     dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
172    
173     dwjpdx = -2.0d0*xj*zj*uglyj/r3
174     dwjpdy = -2.0d0*yj*zj*uglyj/r3
175     dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
176    
177     dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
178     dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
179     dwiduz = - 8.0d0*xi*yi*zi/r3
180    
181     dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
182     dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
183     dwjduz = - 8.0d0*xj*yj*zj/r3
184    
185     dwipdux = 2.0d0*yi*uglyi/rij
186     dwipduy = -2.0d0*xi*uglyi/rij
187     dwipduz = 0.0d0
188    
189     dwjpdux = 2.0d0*yj*uglyj/rij
190     dwjpduy = -2.0d0*xj*uglyj/rij
191     dwjpduz = 0.0d0
192    
193     ! do the torques first since they are easy:
194     ! remember that these are still in the body fixed axes
195    
196     txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
197     tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
198     tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
199    
200     txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
201     tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
202     tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
203    
204     ! go back to lab frame using transpose of rotation matrix:
205 gezelter 319
206     #ifdef IS_MPI
207 gezelter 322 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
208     a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
209     t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
210     a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
211     t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
212     a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
213 gezelter 319
214 gezelter 322 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
215     a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
216     t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
217     a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
218     t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
219     a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
220 gezelter 319 #else
221 gezelter 322 t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
222     t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
223     t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
224    
225     t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
226     t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
227     t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
228     #endif
229     ! Now, on to the forces:
230    
231     ! first rotate the i terms back into the lab frame:
232 gezelter 319
233 gezelter 322 #ifdef IS_MPI
234     fxii = a_Row(1,atom1)*(s*dwidx+sp*dwipdx) + &
235     a_Row(4,atom1)*(s*dwidy+sp*dwipdy) + &
236     a_Row(7,atom1)*(s*dwidz+sp*dwipdz)
237     fyii = a_Row(2,atom1)*(s*dwidx+sp*dwipdx) + &
238     a_Row(5,atom1)*(s*dwidy+sp*dwipdy) + &
239     a_Row(8,atom1)*(s*dwidz+sp*dwipdz)
240     fzii = a_Row(3,atom1)*(s*dwidx+sp*dwipdx) + &
241     a_Row(6,atom1)*(s*dwidy+sp*dwipdy) + &
242     a_Row(9,atom1)*(s*dwidz+sp*dwipdz)
243 gezelter 319
244 gezelter 322 fxjj = a_Col(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
245     a_Col(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
246     a_Col(7,atom2)*(s*dwjdz+sp*dwjpdz)
247     fyjj = a_Col(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
248     a_Col(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
249     a_Col(8,atom2)*(s*dwjdz+sp*dwjpdz)
250     fzjj = a_Col(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
251     a_Col(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
252     a_Col(9,atom2)*(s*dwjdz+sp*dwjpdz)
253 gezelter 319 #else
254 gezelter 322 fxii = a(1,atom1)*(s*dwidx+sp*dwipdx) + &
255     a(4,atom1)*(s*dwidy+sp*dwipdy) + &
256     a(7,atom1)*(s*dwidz+sp*dwipdz)
257     fyii = a(2,atom1)*(s*dwidx+sp*dwipdx) + &
258     a(5,atom1)*(s*dwidy+sp*dwipdy) + &
259     a(8,atom1)*(s*dwidz+sp*dwipdz)
260     fzii = a(3,atom1)*(s*dwidx+sp*dwipdx) + &
261     a(6,atom1)*(s*dwidy+sp*dwipdy) + &
262     a(9,atom1)*(s*dwidz+sp*dwipdz)
263 gezelter 319
264 gezelter 322 fxjj = a(1,atom2)*(s*dwjdx+sp*dwjpdx) + &
265     a(4,atom2)*(s*dwjdy+sp*dwjpdy) + &
266     a(7,atom2)*(s*dwjdz+sp*dwjpdz)
267     fyjj = a(2,atom2)*(s*dwjdx+sp*dwjpdx) + &
268     a(5,atom2)*(s*dwjdy+sp*dwjpdy) + &
269     a(8,atom2)*(s*dwjdz+sp*dwjpdz)
270     fzjj = a(3,atom2)*(s*dwjdx+sp*dwjpdx)+ &
271     a(6,atom2)*(s*dwjdy+sp*dwjpdy) + &
272     a(9,atom2)*(s*dwjdz+sp*dwjpdz)
273 gezelter 319 #endif
274 gezelter 322
275     fxij = -fxii
276     fyij = -fyii
277     fzij = -fzii
278    
279     fxji = -fxjj
280     fyji = -fyjj
281     fzji = -fzjj
282    
283     ! now assemble these with the radial-only terms:
284 gezelter 319
285 gezelter 322 fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
286     fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
287     fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
288 gezelter 319
289     #ifdef IS_MPI
290 gezelter 322 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
291     f_Row(2,atom1) = f_Row(2,atom1) + fyradial
292     f_Row(3,atom1) = f_Row(3,atom1) + fzradial
293    
294     f_Col(1,atom2) = f_Col(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - &
295     dspdr*drdx*wp + fxjj + fxij)
296     f_Col(2,atom2) = f_Col(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - &
297     dspdr*drdy*wp + fyjj + fyij)
298     f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - &
299     dspdr*drdz*wp + fzjj + fzij)
300 gezelter 319 #else
301 gezelter 322 f(1,atom1) = f(1,atom1) + fxradial
302     f(2,atom1) = f(2,atom1) + fyradial
303     f(3,atom1) = f(3,atom1) + fzradial
304    
305     f(1,atom2) = f(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
306     fxjj + fxij)
307     f(2,atom2) = f(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
308     fyjj + fyij)
309     f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
310     fzjj + fzij)
311 gezelter 319 #endif
312 gezelter 322
313     if (doStress()) then
314     tau_Temp(1) = tau_Temp(1) + fxradial * d(1)
315     tau_Temp(2) = tau_Temp(2) + fxradial * d(2)
316     tau_Temp(3) = tau_Temp(3) + fxradial * d(3)
317     tau_Temp(4) = tau_Temp(4) + fyradial * d(1)
318     tau_Temp(5) = tau_Temp(5) + fyradial * d(2)
319     tau_Temp(6) = tau_Temp(6) + fyradial * d(3)
320     tau_Temp(7) = tau_Temp(7) + fzradial * d(1)
321     tau_Temp(8) = tau_Temp(8) + fzradial * d(2)
322     tau_Temp(9) = tau_Temp(9) + fzradial * d(3)
323     virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
324     endif
325    
326     end subroutine do_sticky_pair
327 gezelter 319
328 gezelter 322 !! calculates the switching functions and their derivatives for a given
329     subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
330 gezelter 319
331 gezelter 322 real (kind=dp), intent(in) :: r
332     real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
333 gezelter 319
334 gezelter 322 ! distances must be in angstroms
335 gezelter 319
336 gezelter 322 if (r.lt.SSD_rl) then
337     s = 1.0d0
338     sp = 1.0d0
339     dsdr = 0.0d0
340     dspdr = 0.0d0
341     elseif (r.gt.SSD_rup) then
342     s = 0.0d0
343     sp = 0.0d0
344     dsdr = 0.0d0
345     dspdr = 0.0d0
346     else
347     sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_rup-r)**2) / &
348     ((SSD_rup - SSD_rl)**3)
349     dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rl)/((SSD_rup - SSD_rl)**3)
350 gezelter 319
351 gezelter 322 if (r.gt.SSD_ru) then
352     s = 0.0d0
353     dsdr = 0.0d0
354 gezelter 319 else
355 gezelter 322 s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
356     ((SSD_ru - SSD_rl)**3)
357     dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
358 gezelter 319 endif
359 gezelter 322 endif
360 gezelter 319
361 gezelter 322 return
362     end subroutine calc_sw_fnc
363     end module sticky_pair