1 |
chuckv |
307 |
module dipole_dipole |
2 |
|
|
|
3 |
|
|
use simulation |
4 |
|
|
use definitions |
5 |
|
|
use forceGlobals |
6 |
|
|
use atype_typedefs |
7 |
gezelter |
309 |
#ifdef IS_MPI |
8 |
|
|
use mpiSimulation |
9 |
|
|
#endif |
10 |
|
|
implicit none |
11 |
chuckv |
307 |
|
12 |
|
|
contains |
13 |
|
|
|
14 |
gezelter |
309 |
subroutine do_dipole_pair(atom1, atom2, atype1, atype2, dx, dy, dz, rij, & |
15 |
|
|
rt, rrf, pot, u_l, f, t) |
16 |
chuckv |
307 |
|
17 |
|
|
integer atom1, atom2 |
18 |
|
|
double precision dx, dy, dz, rij |
19 |
|
|
double precision dfact1, dfact2, dip2, r2, r3, r5, pre |
20 |
|
|
double precision dudx, dudy, dudz, dudu1x, dudu1y, dudu1z |
21 |
|
|
double precision dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2 |
22 |
|
|
double precision taper, dtdr, vterm |
23 |
gezelter |
309 |
|
24 |
|
|
real( kind = dp ) :: pot, rt, rrf |
25 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
26 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: f |
27 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: t |
28 |
chuckv |
307 |
|
29 |
|
|
real (kind = dp), dimension(3) :: ul1 |
30 |
|
|
real (kind = dp), dimension(3) :: ul2 |
31 |
|
|
type (atype), pointer :: atype1 |
32 |
|
|
type (atype), pointer :: atype2 |
33 |
|
|
|
34 |
gezelter |
309 |
#ifdef IS_MPI |
35 |
|
|
ul1(1) = u_l_Row(1,atom1) |
36 |
|
|
ul1(2) = u_l_Row(2,atom1) |
37 |
|
|
ul1(3) = u_l_Row(3,atom1) |
38 |
|
|
|
39 |
|
|
ul2(1) = u_l_Col(1,atom2) |
40 |
|
|
ul2(2) = u_l_Col(2,atom2) |
41 |
|
|
ul2(3) = u_l_Col(3,atom2) |
42 |
|
|
#else |
43 |
|
|
ul1(1) = u_l(1,atom1) |
44 |
|
|
ul1(2) = u_l(2,atom1) |
45 |
|
|
ul1(3) = u_l(3,atom1) |
46 |
|
|
|
47 |
|
|
ul1(1) = u_l(1,atom2) |
48 |
|
|
ul1(2) = u_l(2,atom2) |
49 |
|
|
ul1(3) = u_l(3,atom2) |
50 |
|
|
#endif |
51 |
|
|
|
52 |
chuckv |
307 |
! pre converts from mu in units of debye to kcal/mol |
53 |
|
|
pre = 14.38362d0 |
54 |
|
|
|
55 |
|
|
if (rij.le.rrf) then |
56 |
|
|
|
57 |
|
|
if (rij.lt.rt) then |
58 |
|
|
taper = 1.0d0 |
59 |
|
|
dtdr = 0.0d0 |
60 |
|
|
else |
61 |
|
|
taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3) |
62 |
|
|
dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3) |
63 |
|
|
endif |
64 |
|
|
|
65 |
|
|
r2 = rij*rij |
66 |
|
|
r3 = r2*rij |
67 |
|
|
r5 = r3*r2 |
68 |
|
|
|
69 |
|
|
rdotu1 = dx*ul1(1) + dy*ul1(2) + dz*ul1(3) |
70 |
|
|
rdotu2 = dx*ul2(1) + dy*ul2(2) + dz*ul2(3) |
71 |
|
|
u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3) |
72 |
|
|
|
73 |
|
|
dip2 = pre * atype1%dipole_moment * atype2%dipole_moment |
74 |
|
|
|
75 |
|
|
dfact1 = 3.0d0*dip2 / r2 |
76 |
|
|
dfact2 = 3.0d0*dip2 / r5 |
77 |
|
|
|
78 |
|
|
vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5)) |
79 |
|
|
|
80 |
gezelter |
309 |
#ifdef IS_MPI |
81 |
|
|
pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper |
82 |
|
|
pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper |
83 |
|
|
#else |
84 |
chuckv |
307 |
pot = pot + vterm*taper |
85 |
gezelter |
309 |
#endif |
86 |
chuckv |
307 |
|
87 |
|
|
dudx = (-dfact1 * dx * ((u1dotu2/r3) - & |
88 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
89 |
|
|
dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + & |
90 |
|
|
vterm*dtdr*dx/rij |
91 |
|
|
|
92 |
|
|
dudy = (-dfact1 * dy * ((u1dotu2/r3) - & |
93 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
94 |
|
|
dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + & |
95 |
|
|
vterm*dtdr*dy/rij |
96 |
|
|
|
97 |
|
|
dudz = (-dfact1 * dz * ((u1dotu2/r3) - & |
98 |
|
|
(5.0d0*(rdotu1*rdotu2)/r5)) - & |
99 |
|
|
dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + & |
100 |
|
|
vterm*dtdr*dz/rij |
101 |
|
|
|
102 |
|
|
dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*dx*rdotu2/r5)))*taper |
103 |
|
|
dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*dy*rdotu2/r5)))*taper |
104 |
|
|
dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*dz*rdotu2/r5)))*taper |
105 |
|
|
|
106 |
|
|
dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*dx*rdotu1/r5)))*taper |
107 |
|
|
dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*dy*rdotu1/r5)))*taper |
108 |
|
|
dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*dz*rdotu1/r5)))*taper |
109 |
|
|
|
110 |
|
|
|
111 |
|
|
#ifdef IS_MPI |
112 |
gezelter |
309 |
f_Row(1,atom1) = f_Row(1,atom1) + dudx |
113 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + dudy |
114 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + dudz |
115 |
chuckv |
307 |
|
116 |
gezelter |
309 |
f_Col(1,atom2) = f_Col(1,atom2) - dudx |
117 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - dudy |
118 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - dudz |
119 |
chuckv |
307 |
|
120 |
gezelter |
309 |
t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y |
121 |
|
|
t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z |
122 |
|
|
t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x |
123 |
chuckv |
307 |
|
124 |
gezelter |
309 |
t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y |
125 |
|
|
t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z |
126 |
|
|
t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x |
127 |
chuckv |
307 |
#else |
128 |
gezelter |
309 |
f(1,atom1) = f(1,atom1) + dudx |
129 |
|
|
f(2,atom1) = f(2,atom1) + dudy |
130 |
|
|
f(3,atom1) = f(3,atom1) + dudz |
131 |
chuckv |
307 |
|
132 |
gezelter |
309 |
f(1,atom2) = f(1,atom2) - dudx |
133 |
|
|
f(2,atom2) = f(2,atom2) - dudy |
134 |
|
|
f(3,atom2) = f(3,atom2) - dudz |
135 |
chuckv |
307 |
|
136 |
gezelter |
309 |
t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y |
137 |
|
|
t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z |
138 |
|
|
t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x |
139 |
chuckv |
307 |
|
140 |
gezelter |
309 |
t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y |
141 |
|
|
t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z |
142 |
|
|
t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x |
143 |
chuckv |
307 |
#endif |
144 |
|
|
|
145 |
gezelter |
309 |
if (doStress()) then |
146 |
|
|
tau_Temp(1) = tau_Temp(1) + dudx * dx |
147 |
|
|
tau_Temp(2) = tau_Temp(2) + dudx * dy |
148 |
|
|
tau_Temp(3) = tau_Temp(3) + dudx * dz |
149 |
|
|
tau_Temp(4) = tau_Temp(4) + dudy * dx |
150 |
|
|
tau_Temp(5) = tau_Temp(5) + dudy * dy |
151 |
|
|
tau_Temp(6) = tau_Temp(6) + dudy * dz |
152 |
|
|
tau_Temp(7) = tau_Temp(7) + dudz * dx |
153 |
|
|
tau_Temp(8) = tau_Temp(8) + dudz * dy |
154 |
|
|
tau_Temp(9) = tau_Temp(9) + dudz * dz |
155 |
|
|
virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) |
156 |
chuckv |
307 |
endif |
157 |
|
|
|
158 |
|
|
endif |
159 |
|
|
|
160 |
|
|
return |
161 |
gezelter |
309 |
end subroutine do_dipole_pair |
162 |
|
|
|
163 |
chuckv |
307 |
end module dipole_dipole |