1 |
chuckv |
298 |
subroutine dipole_dipole(atom1, atom2, atype1, atype2, dx, dy, dz, rij, & |
2 |
|
|
ul1, ul2, rt, rrf, pot) |
3 |
|
|
|
4 |
|
|
include 'sizes.inc' |
5 |
|
|
include 'simulation.inc' |
6 |
|
|
|
7 |
|
|
integer atom1, atom2 |
8 |
|
|
double precision dx, dy, dz, rij |
9 |
|
|
double precision dfact1, dfact2, dip2, r2, r3, r5, pre |
10 |
|
|
double precision dudx, dudy, dudz, dudu1x, dudu1y, dudu1z |
11 |
|
|
double precision dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2 |
12 |
|
|
double precision taper, dtdr, vterm |
13 |
|
|
|
14 |
|
|
real (kind = dp), dimension(3) :: ul1 |
15 |
|
|
real (kind = dp), dimension(3) :: ul2 |
16 |
|
|
type (atype), pointer :: atype1 |
17 |
|
|
type (atype), pointer :: atype2 |
18 |
|
|
|
19 |
|
|
! pre converts from mu in units of debye to kcal/mol |
20 |
|
|
pre = 14.38362d0 |
21 |
|
|
|
22 |
|
|
if (rij.le.rrf) then |
23 |
|
|
|
24 |
|
|
if (rij.lt.rt) then |
25 |
|
|
taper = 1.0d0 |
26 |
|
|
dtdr = 0.0d0 |
27 |
|
|
else |
28 |
|
|
taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3) |
29 |
|
|
dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3) |
30 |
|
|
endif |
31 |
|
|
|
32 |
|
|
r2 = rij*rij |
33 |
|
|
r3 = r2*rij |
34 |
|
|
r5 = r3*r2 |
35 |
|
|
|
36 |
|
|
rdotu1 = dx*ul1(1) + dy*ul1(2) + dz*ul1(3) |
37 |
|
|
rdotu2 = dx*ul2(1) + dy*ul2(2) + dz*ul2(3) |
38 |
|
|
u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3) |
39 |
|
|
|
40 |
|
|
dip2 = pre * atype1%dipole_moment * atype2%dipole_moment |
41 |
|
|
|
42 |
|
|
dfact1 = 3.0d0*dip2 / r2 |
43 |
|
|
dfact2 = 3.0d0*dip2 / r5 |
44 |
|
|
|
45 |
|
|
vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5)) |
46 |
|
|
|
47 |
|
|
pot = pot + vterm*taper |
48 |
|
|
|
49 |
|
|
dudx = (-dfact1 * dx * ((u1dotu2/r3) - & |
50 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
51 |
|
|
dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + & |
52 |
|
|
vterm*dtdr*dx/rij |
53 |
|
|
|
54 |
|
|
dudy = (-dfact1 * dy * ((u1dotu2/r3) - & |
55 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
56 |
|
|
dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + & |
57 |
|
|
vterm*dtdr*dy/rij |
58 |
|
|
|
59 |
|
|
dudz = (-dfact1 * dz * ((u1dotu2/r3) - & |
60 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
61 |
|
|
dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + & |
62 |
|
|
vterm*dtdr*dz/rij |
63 |
|
|
|
64 |
|
|
dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*dx*rdotu2/r5)))*taper |
65 |
|
|
dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*dy*rdotu2/r5)))*taper |
66 |
|
|
dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*dz*rdotu2/r5)))*taper |
67 |
|
|
|
68 |
|
|
dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*dx*rdotu1/r5)))*taper |
69 |
|
|
dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*dy*rdotu1/r5)))*taper |
70 |
|
|
dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*dz*rdotu1/r5)))*taper |
71 |
|
|
|
72 |
|
|
|
73 |
|
|
#ifdef IS_MPI |
74 |
|
|
fRow(1,atom1) = fRow(1,atom1) + dudx |
75 |
|
|
fRow(2,atom1) = fRow(2,atom1) + dudy |
76 |
|
|
fRow(3,atom1) = fRow(3,atom1) + dudz |
77 |
|
|
|
78 |
|
|
fCol(1,atom2) = fCol(1,atom2) - dudx |
79 |
|
|
fCol(2,atom2) = fCol(2,atom2) - dudy |
80 |
|
|
fCol(3,atom2) = fCol(3,atom2) - dudz |
81 |
|
|
|
82 |
|
|
tRow(1,atom1) = tRow(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y |
83 |
|
|
tRow(2,atom1) = tRow(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z |
84 |
|
|
tRow(3,atom1) = tRow(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x |
85 |
|
|
|
86 |
|
|
tCol(1,atom2) = tCol(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y |
87 |
|
|
tCol(2,atom2) = tCol(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z |
88 |
|
|
tCol(3,atom2) = tCol(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x |
89 |
|
|
#else |
90 |
|
|
fTemp(1,atom1) = fTemp(1,atom1) + dudx |
91 |
|
|
fTemp(2,atom1) = fTemp(2,atom1) + dudy |
92 |
|
|
fTemp(3,atom1) = fTemp(3,atom1) + dudz |
93 |
|
|
|
94 |
|
|
fTemp(1,atom2) = fTemp(1,atom2) - dudx |
95 |
|
|
fTemp(2,atom2) = fTemp(2,atom2) - dudy |
96 |
|
|
fTemp(3,atom2) = fTemp(3,atom2) - dudz |
97 |
|
|
|
98 |
|
|
tTemp(1,atom1) = tTemp(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y |
99 |
|
|
tTemp(2,atom1) = tTemp(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z |
100 |
|
|
tTemp(3,atom1) = tTemp(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x |
101 |
|
|
|
102 |
|
|
tTemp(1,atom2) = tTemp(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y |
103 |
|
|
tTemp(2,atom2) = tTemp(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z |
104 |
|
|
tTemp(3,atom2) = tTemp(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x |
105 |
|
|
#endif |
106 |
|
|
|
107 |
|
|
if (do_stress) then |
108 |
|
|
tauTemp(1) = tauTemp(1) + dudx * dx |
109 |
|
|
tauTemp(2) = tauTemp(2) + dudx * dy |
110 |
|
|
tauTemp(3) = tauTemp(3) + dudx * dz |
111 |
|
|
tauTemp(4) = tauTemp(4) + dudy * dx |
112 |
|
|
tauTemp(5) = tauTemp(5) + dudy * dy |
113 |
|
|
tauTemp(6) = tauTemp(6) + dudy * dz |
114 |
|
|
tauTemp(7) = tauTemp(7) + dudz * dx |
115 |
|
|
tauTemp(8) = tauTemp(8) + dudz * dy |
116 |
|
|
tauTemp(9) = tauTemp(9) + dudz * dz |
117 |
|
|
virialTemp = virialTemp + ( tauTemp(1) + tauTemp(5) + tauTemp(9) ) |
118 |
|
|
endif |
119 |
|
|
|
120 |
|
|
endif |
121 |
|
|
|
122 |
|
|
return |
123 |
|
|
end subroutine dipole_dipole |