ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 611
Committed: Tue Jul 15 17:10:50 2003 UTC (20 years, 11 months ago) by gezelter
File size: 13219 byte(s)
Log Message:
Fixing  pressure tensor

File Contents

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