ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/f_ssdForces.f90
Revision: 11
Committed: Tue Jul 9 18:40:59 2002 UTC (21 years, 11 months ago) by mmeineke
File size: 7118 byte(s)
Log Message:
This commit was generated by cvs2svn to compensate for changes in r10, which
included commits to RCS files with non-trunk default branches.

File Contents

# User Rev Content
1 mmeineke 10 subroutine ssd_forces(natoms, i, j, atom1, atom2, dx, dy, dz, rij, pot, r2, &
2     flx, fly, flz, tlx, tly, tlz, a)
3    
4     implicit none
5    
6     include 'f_ssd.inc'
7    
8     ! This routine does only the sticky portion of the SSD potential
9     ! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
10     ! The Lennard-Jones and dipolar interaction must be handled separately.
11    
12     ! We assume that the rotation matrices have already been calculated
13     ! and placed in the A(9, max_mol) array.
14     !
15     ! i and j are pointers to the two SSD molecules, while atom1 and
16     ! atom2 are the pointers to the location of the atoms in the force
17     ! and position arrays. These are both necessary since the rotation
18     ! matrix is a property of the molecule, while the force is acting on
19     ! the atom. The indexing of atoms and molecules is not necessarily
20     ! the same in simulations of mixtures.
21    
22    
23     integer :: i, j, atom1, atom2, natoms
24     double precision dx, dy, dz, rij, pot, r2
25     double precision, dimension(natoms) ::flx, fly, flz, tlx, tly, tlz
26     double precision, dimension(9,natoms):: a
27    
28    
29     double precision xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
30     double precision r3, r5, r6, s, sp, dsdr, dspdr
31     double precision wi, wj, w, wip, wjp, wp
32     double precision dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
33     double precision dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
34     double precision dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
35     double precision dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
36     double precision zif, zis, zjf, zjs, uglyi, uglyj
37     double precision drdx, drdy, drdz
38     double precision txi, tyi, tzi, txj, tyj, tzj
39     double precision fxii, fyii, fzii, fxjj, fyjj, fzjj
40     double precision fxij, fyij, fzij, fxji, fyji, fzji
41    
42    
43     ! Use molecular positions, since the SSD model has only one atom, and the
44     ! rotation matrix is for the molecule itself:
45    
46     r3 = r2*rij
47     r5 = r3*r2
48    
49     drdx = dx / rij
50     drdy = dy / rij
51     drdz = dz / rij
52    
53     ! rotate the inter-particle separation into the two different
54     ! body-fixed coordinate systems:
55    
56     xi = a(1,i)*dx + a(2,i)*dy + a(3,i)*dz
57     yi = a(4,i)*dx + a(5,i)*dy + a(6,i)*dz
58     zi = a(7,i)*dx + a(8,i)*dy + a(9,i)*dz
59    
60     ! negative sign because this is the vector from j to i:
61    
62     xj = -(a(1,j)*dx + a(2,j)*dy + a(3,j)*dz)
63     yj = -(a(4,j)*dx + a(5,j)*dy + a(6,j)*dz)
64     zj = -(a(7,j)*dx + a(8,j)*dy + a(9,j)*dz)
65    
66     xi2 = xi*xi
67     yi2 = yi*yi
68     zi2 = zi*zi
69    
70     xj2 = xj*xj
71     yj2 = yj*yj
72     zj2 = zj*zj
73    
74     call calc_s(rij, s, sp, dsdr, dspdr)
75    
76     wi = 2.0d0*(xi2-yi2)*zi / r3
77     wj = 2.0d0*(xj2-yj2)*zj / r3
78     w = wi+wj
79    
80     zif = zi/rij - 0.6d0
81     zis = zi/rij + 0.8d0
82    
83     zjf = zj/rij - 0.6d0
84     zjs = zj/rij + 0.8d0
85    
86     wip = zif*zif*zis*zis - SSD_w0
87     wjp = zjf*zjf*zjs*zjs - SSD_w0
88     wp = wip + wjp
89    
90     pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
91    
92     dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
93     dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
94     dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
95    
96     dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
97     dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
98     dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
99    
100     uglyi = zif*zif*zis + zif*zis*zis
101     uglyj = zjf*zjf*zjs + zjf*zjs*zjs
102    
103     dwipdx = -2.0d0*xi*zi*uglyi/r3
104     dwipdy = -2.0d0*yi*zi*uglyi/r3
105     dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
106    
107     dwjpdx = -2.0d0*xj*zj*uglyj/r3
108     dwjpdy = -2.0d0*yj*zj*uglyj/r3
109     dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
110    
111     dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
112     dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
113     dwiduz = - 8.0d0*xi*yi*zi/r3
114    
115     dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
116     dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
117     dwjduz = - 8.0d0*xj*yj*zj/r3
118    
119     dwipdux = 2.0d0*yi*uglyi/rij
120     dwipduy = -2.0d0*xi*uglyi/rij
121     dwipduz = 0.0d0
122    
123     dwjpdux = 2.0d0*yj*uglyj/rij
124     dwjpduy = -2.0d0*xj*uglyj/rij
125     dwjpduz = 0.0d0
126    
127     ! do the torques first since they are easy:
128     ! remember that these are still in the body fixed axes
129    
130     txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
131     tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
132     tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
133    
134     txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
135     tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
136     tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
137    
138     ! go back to lab frame using transpose of rotation matrix:
139    
140     tlx(atom1) = tlx(atom1) + a(1,i)*txi + a(4,i)*tyi + a(7,i)*tzi
141     tly(atom1) = tly(atom1) + a(2,i)*txi + a(5,i)*tyi + a(8,i)*tzi
142     tlz(atom1) = tlz(atom1) + a(3,i)*txi + a(6,i)*tyi + a(9,i)*tzi
143    
144     tlx(atom2) = tlx(atom2) + a(1,j)*txj + a(4,j)*tyj + a(7,j)*tzj
145     tly(atom2) = tly(atom2) + a(2,j)*txj + a(5,j)*tyj + a(8,j)*tzj
146     tlz(atom2) = tlz(atom2) + a(3,j)*txj + a(6,j)*tyj + a(9,j)*tzj
147    
148     ! Now, on to the forces:
149    
150     ! first rotate the i terms back into the lab frame:
151    
152     fxii = a(1,i)*(s*dwidx+sp*dwipdx) + &
153     a(4,i)*(s*dwidy+sp*dwipdy) + &
154     a(7,i)*(s*dwidz+sp*dwipdz)
155     fyii = a(2,i)*(s*dwidx+sp*dwipdx) + &
156     a(5,i)*(s*dwidy+sp*dwipdy) + &
157     a(8,i)*(s*dwidz+sp*dwipdz)
158     fzii = a(3,i)*(s*dwidx+sp*dwipdx) + &
159     a(6,i)*(s*dwidy+sp*dwipdy) + &
160     a(9,i)*(s*dwidz+sp*dwipdz)
161    
162     fxij = -fxii
163     fyij = -fyii
164     fzij = -fzii
165    
166     fxjj = a(1,j)*(s*dwjdx+sp*dwjpdx) + &
167     a(4,j)*(s*dwjdy+sp*dwjpdy) + &
168     a(7,j)*(s*dwjdz+sp*dwjpdz)
169     fyjj = a(2,j)*(s*dwjdx+sp*dwjpdx) + &
170     a(5,j)*(s*dwjdy+sp*dwjpdy) + &
171     a(8,j)*(s*dwjdz+sp*dwjpdz)
172     fzjj = a(3,j)*(s*dwjdx+sp*dwjpdx)+ &
173     a(6,j)*(s*dwjdy+sp*dwjpdy) + &
174     a(9,j)*(s*dwjdz+sp*dwjpdz)
175    
176     fxji = -fxjj
177     fyji = -fyjj
178     fzji = -fzjj
179    
180     ! now assemble these with the radial-only terms:
181    
182     flx(atom1) = flx(atom1) + 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + &
183     fxii + fxji)
184     fly(atom1) = fly(atom1) + 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + &
185     fyii + fyji)
186     flz(atom1) = flz(atom1) + 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + &
187     fzii + fzji)
188    
189     flx(atom2) = flx(atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
190     fxjj + fxij)
191     fly(atom2) = fly(atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
192     fyjj + fyij)
193     flz(atom2) = flz(atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
194     fzjj + fzij)
195    
196     return
197     end subroutine ssd_forces
198    
199     subroutine calc_s(r, s, sp, dsdr, dspdr)
200    
201     ! calculates the switching functions and their derivatives for a given
202     ! value of r
203    
204     double precision r, s, sp, dsdr, dspdr
205     double precision rl, ru, rup
206     ! distances are in angstroms
207     parameter(rl = 2.75d0, ru = 3.35d0, rup = 4.0d0)
208    
209     if (r.lt.rl) then
210     s = 1.0d0
211     sp = 1.0d0
212     dsdr = 0.0d0
213     dspdr = 0.0d0
214     elseif (r.gt.rup) then
215     s = 0.0d0
216     sp = 0.0d0
217     dsdr = 0.0d0
218     dspdr = 0.0d0
219     else
220     sp = ((rup + 2.0d0*r - 3.0d0*rl) * (rup-r)**2)/((rup - rl)**3)
221     dspdr = 6.0d0*(r-rup)*(r-rl)/((rup - rl)**3)
222    
223     if (r.gt.ru) then
224     s = 0.0d0
225     dsdr = 0.0d0
226     else
227     s = ((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2)/((ru - rl)**3)
228     dsdr = 6.0d0*(r-ru)*(r-rl)/((ru - rl)**3)
229     endif
230     endif
231    
232     return
233     end subroutine calc_s
234    
235    
236    
237