ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
Revision: 843
Committed: Wed Oct 29 20:41:39 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 14040 byte(s)
Log Message:
fixed a stdlib.h include error in bass.l

fixed a little bug in the first time step, regarding the setting of ecr and est in fortran

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