1 |
module gb_pair |
2 |
use force_globals |
3 |
use definitions |
4 |
#ifdef IS_MPI |
5 |
use mpiSimulation |
6 |
#endif |
7 |
|
8 |
implicit none |
9 |
|
10 |
PRIVATE |
11 |
|
12 |
logical, save :: gb_pair_initialized = .false. |
13 |
real(kind=dp), save :: gb_sigma |
14 |
real(kind=dp), save :: gb_l2b_ratio |
15 |
real(kind=dp), save :: gb_eps |
16 |
real(kind=dp), save :: gb_eps_ratio |
17 |
real(kind=dp), save :: gb_mu |
18 |
real(kind=dp), save :: gb_nu |
19 |
|
20 |
public :: check_gb_pair_FF |
21 |
public :: set_gb_pair_params |
22 |
public :: do_gb_pair |
23 |
|
24 |
contains |
25 |
|
26 |
subroutine check_gb_pair_FF(status) |
27 |
integer :: status |
28 |
status = -1 |
29 |
if (gb_pair_initialized) status = 0 |
30 |
return |
31 |
end subroutine check_gb_pair_FF |
32 |
|
33 |
subroutine set_gb_pair_params(sigma, l2b_ratio, eps, eps_ratio, mu, nu) |
34 |
real( kind = dp ), intent(in) :: sigma, l2b_ratio, eps, eps_ratio |
35 |
real( kind = dp ), intent(in) :: mu, nu |
36 |
|
37 |
gb_sigma = sigma |
38 |
gb_l2b_ratio = l2b_ratio |
39 |
gb_eps = eps |
40 |
gb_eps_ratio = eps_ratio |
41 |
gb_mu = mu |
42 |
gb_nu = nu |
43 |
|
44 |
gb_pair_initialized = .true. |
45 |
return |
46 |
end subroutine set_gb_pair_params |
47 |
|
48 |
|
49 |
subroutine do_gb_pair(atom1, atom2, d, r, r2, u_l, pot, f, t, & |
50 |
do_pot, do_stress) |
51 |
|
52 |
integer, intent(in) :: atom1, atom2 |
53 |
real (kind=dp), intent(inout) :: r, r2 |
54 |
real (kind=dp), dimension(3), intent(in) :: d |
55 |
real (kind=dp) :: pot |
56 |
real (kind=dp), dimension(:,:) :: u_l |
57 |
real (kind=dp), dimension(:,:) :: f |
58 |
real (kind=dp), dimension(:,:) :: t |
59 |
logical, intent(in) :: do_pot, do_stress |
60 |
real (kind = dp), dimension(3) :: ul1 |
61 |
real (kind = dp), dimension(3) :: ul2 |
62 |
|
63 |
real(kind=dp) :: chi, chiprime, emu, s2 |
64 |
real(kind=dp) :: r4, rdotu1, rdotu2, u1dotu2, g, gp, gpi, gmu, gmum |
65 |
real(kind=dp) :: curlyE, enu, enum, eps, dotsum, dotdiff, ds2, dd2 |
66 |
real(kind=dp) :: opXdot, omXdot, opXpdot, omXpdot, pref, gfact |
67 |
real(kind=dp) :: BigR, Ri, Ri2, Ri6, Ri7, Ri12, Ri13, R126, R137 |
68 |
real(kind=dp) :: dru1dx, dru1dy, dru1dz |
69 |
real(kind=dp) :: dru2dx, dru2dy, dru2dz |
70 |
real(kind=dp) :: dBigRdx, dBigRdy, dBigRdz |
71 |
real(kind=dp) :: dBigRdu1x, dBigRdu1y, dBigRdu1z |
72 |
real(kind=dp) :: dBigRdu2x, dBigRdu2y, dBigRdu2z |
73 |
real(kind=dp) :: dUdx, dUdy, dUdz |
74 |
real(kind=dp) :: dUdu1x, dUdu1y, dUdu1z, dUdu2x, dUdu2y, dUdu2z |
75 |
real(kind=dp) :: dcE, dcEdu1x, dcEdu1y, dcEdu1z, dcEdu2x, dcEdu2y, dcEdu2z |
76 |
real(kind=dp) :: depsdu1x, depsdu1y, depsdu1z, depsdu2x, depsdu2y, depsdu2z |
77 |
real(kind=dp) :: drdx, drdy, drdz |
78 |
real(kind=dp) :: dgdx, dgdy, dgdz |
79 |
real(kind=dp) :: dgdu1x, dgdu1y, dgdu1z, dgdu2x, dgdu2y, dgdu2z |
80 |
real(kind=dp) :: dgpdx, dgpdy, dgpdz |
81 |
real(kind=dp) :: dgpdu1x, dgpdu1y, dgpdu1z, dgpdu2x, dgpdu2y, dgpdu2z |
82 |
real(kind=dp) :: line1a, line1bx, line1by, line1bz |
83 |
real(kind=dp) :: line2a, line2bx, line2by, line2bz |
84 |
real(kind=dp) :: line3a, line3b, line3, line3x, line3y, line3z |
85 |
real(kind=dp) :: term1x, term1y, term1z, term1u1x, term1u1y, term1u1z |
86 |
real(kind=dp) :: term1u2x, term1u2y, term1u2z |
87 |
real(kind=dp) :: term2a, term2b, term2u1x, term2u1y, term2u1z |
88 |
real(kind=dp) :: term2u2x, term2u2y, term2u2z |
89 |
real(kind=dp) :: yick1, yick2, mess1, mess2 |
90 |
|
91 |
s2 = (gb_l2b_ratio)**2 |
92 |
emu = (gb_eps_ratio)**(1.0d0/gb_mu) |
93 |
|
94 |
chi = (s2 - 1.0d0)/(s2 + 1.0d0) |
95 |
chiprime = (1.0d0 - emu)/(1.0d0 + emu) |
96 |
|
97 |
r4 = r2*r2 |
98 |
|
99 |
#ifdef IS_MPI |
100 |
ul1(1) = u_l_Row(1,atom1) |
101 |
ul1(2) = u_l_Row(2,atom1) |
102 |
ul1(3) = u_l_Row(3,atom1) |
103 |
|
104 |
ul2(1) = u_l_Col(1,atom2) |
105 |
ul2(2) = u_l_Col(2,atom2) |
106 |
ul2(3) = u_l_Col(3,atom2) |
107 |
#else |
108 |
ul1(1) = u_l(1,atom1) |
109 |
ul1(2) = u_l(2,atom1) |
110 |
ul1(3) = u_l(3,atom1) |
111 |
|
112 |
ul2(1) = u_l(1,atom2) |
113 |
ul2(2) = u_l(2,atom2) |
114 |
ul2(3) = u_l(3,atom2) |
115 |
#endif |
116 |
|
117 |
dru1dx = ul1(1) |
118 |
dru2dx = ul2(1) |
119 |
dru1dy = ul1(2) |
120 |
dru2dy = ul2(2) |
121 |
dru1dz = ul1(3) |
122 |
dru2dz = ul2(3) |
123 |
|
124 |
drdx = d(1) / r |
125 |
drdy = d(2) / r |
126 |
drdz = d(3) / r |
127 |
|
128 |
! do some dot products: |
129 |
! NB the r in these dot products is the actual intermolecular vector, |
130 |
! and is not the unit vector in that direction. |
131 |
|
132 |
rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3) |
133 |
rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3) |
134 |
u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3) |
135 |
|
136 |
! This stuff is all for the calculation of g(Chi) and dgdx |
137 |
! Line numbers roughly follow the lines in equation A25 of Luckhurst |
138 |
! et al. Liquid Crystals 8, 451-464 (1990). |
139 |
! We note however, that there are some major typos in that Appendix |
140 |
! of the Luckhurst paper, particularly in equations A23, A29 and A31 |
141 |
! We have attempted to correct them below. |
142 |
|
143 |
dotsum = rdotu1+rdotu2 |
144 |
dotdiff = rdotu1-rdotu2 |
145 |
ds2 = dotsum*dotsum |
146 |
dd2 = dotdiff*dotdiff |
147 |
|
148 |
opXdot = 1.0d0 + Chi*u1dotu2 |
149 |
omXdot = 1.0d0 - Chi*u1dotu2 |
150 |
opXpdot = 1.0d0 + ChiPrime*u1dotu2 |
151 |
omXpdot = 1.0d0 - ChiPrime*u1dotu2 |
152 |
|
153 |
line1a = dotsum/opXdot |
154 |
line1bx = dru1dx + dru2dx |
155 |
line1by = dru1dy + dru2dy |
156 |
line1bz = dru1dz + dru2dz |
157 |
|
158 |
line2a = dotdiff/omXdot |
159 |
line2bx = dru1dx - dru2dx |
160 |
line2by = dru1dy - dru2dy |
161 |
line2bz = dru1dz - dru2dz |
162 |
|
163 |
term1x = -Chi*(line1a*line1bx + line2a*line2bx)/r2 |
164 |
term1y = -Chi*(line1a*line1by + line2a*line2by)/r2 |
165 |
term1z = -Chi*(line1a*line1bz + line2a*line2bz)/r2 |
166 |
|
167 |
line3a = ds2/opXdot |
168 |
line3b = dd2/omXdot |
169 |
line3 = Chi*(line3a + line3b)/r4 |
170 |
line3x = d(1)*line3 |
171 |
line3y = d(2)*line3 |
172 |
line3z = d(3)*line3 |
173 |
|
174 |
dgdx = term1x + line3x |
175 |
dgdy = term1y + line3y |
176 |
dgdz = term1z + line3z |
177 |
|
178 |
term1u1x = 2.0d0*(line1a+line2a)*d(1) |
179 |
term1u1y = 2.0d0*(line1a+line2a)*d(2) |
180 |
term1u1z = 2.0d0*(line1a+line2a)*d(3) |
181 |
term1u2x = 2.0d0*(line1a-line2a)*d(1) |
182 |
term1u2y = 2.0d0*(line1a-line2a)*d(2) |
183 |
term1u2z = 2.0d0*(line1a-line2a)*d(3) |
184 |
|
185 |
term2a = -line3a/opXdot |
186 |
term2b = line3b/omXdot |
187 |
|
188 |
term2u1x = Chi*ul2(1)*(term2a + term2b) |
189 |
term2u1y = Chi*ul2(2)*(term2a + term2b) |
190 |
term2u1z = Chi*ul2(3)*(term2a + term2b) |
191 |
term2u2x = Chi*ul1(1)*(term2a + term2b) |
192 |
term2u2y = Chi*ul1(2)*(term2a + term2b) |
193 |
term2u2z = Chi*ul1(3)*(term2a + term2b) |
194 |
|
195 |
pref = -Chi*0.5d0/r2 |
196 |
|
197 |
dgdu1x = pref*(term1u1x+term2u1x) |
198 |
dgdu1y = pref*(term1u1y+term2u1y) |
199 |
dgdu1z = pref*(term1u1z+term2u1z) |
200 |
dgdu2x = pref*(term1u2x+term2u2x) |
201 |
dgdu2y = pref*(term1u2y+term2u2y) |
202 |
dgdu2z = pref*(term1u2z+term2u2z) |
203 |
|
204 |
g = 1.0d0 - Chi*(line3a + line3b)/(2.0d0*r2) |
205 |
|
206 |
BigR = (r - gb_sigma*(g**(-0.5d0)) + gb_sigma)/gb_sigma |
207 |
Ri = 1.0d0/BigR |
208 |
Ri2 = Ri*Ri |
209 |
Ri6 = Ri2*Ri2*Ri2 |
210 |
Ri7 = Ri6*Ri |
211 |
Ri12 = Ri6*Ri6 |
212 |
Ri13 = Ri6*Ri7 |
213 |
|
214 |
gfact = (g**(-1.5d0))*0.5d0 |
215 |
|
216 |
dBigRdx = drdx/gb_sigma + dgdx*gfact |
217 |
dBigRdy = drdy/gb_sigma + dgdy*gfact |
218 |
dBigRdz = drdz/gb_sigma + dgdz*gfact |
219 |
dBigRdu1x = dgdu1x*gfact |
220 |
dBigRdu1y = dgdu1y*gfact |
221 |
dBigRdu1z = dgdu1z*gfact |
222 |
dBigRdu2x = dgdu2x*gfact |
223 |
dBigRdu2y = dgdu2y*gfact |
224 |
dBigRdu2z = dgdu2z*gfact |
225 |
|
226 |
! Now, we must do it again for g(ChiPrime) and dgpdx |
227 |
|
228 |
line1a = dotsum/opXpdot |
229 |
line2a = dotdiff/omXpdot |
230 |
term1x = -ChiPrime*(line1a*line1bx + line2a*line2bx)/r2 |
231 |
term1y = -ChiPrime*(line1a*line1by + line2a*line2by)/r2 |
232 |
term1z = -ChiPrime*(line1a*line1bz + line2a*line2bz)/r2 |
233 |
line3a = ds2/opXpdot |
234 |
line3b = dd2/omXpdot |
235 |
line3 = ChiPrime*(line3a + line3b)/r4 |
236 |
line3x = d(1)*line3 |
237 |
line3y = d(2)*line3 |
238 |
line3z = d(3)*line3 |
239 |
|
240 |
dgpdx = term1x + line3x |
241 |
dgpdy = term1y + line3y |
242 |
dgpdz = term1z + line3z |
243 |
|
244 |
term1u1x = 2.0d0*(line1a+line2a)*d(1) |
245 |
term1u1y = 2.0d0*(line1a+line2a)*d(2) |
246 |
term1u1z = 2.0d0*(line1a+line2a)*d(3) |
247 |
term1u2x = 2.0d0*(line1a-line2a)*d(1) |
248 |
term1u2y = 2.0d0*(line1a-line2a)*d(2) |
249 |
term1u2z = 2.0d0*(line1a-line2a)*d(3) |
250 |
|
251 |
term2a = -line3a/opXpdot |
252 |
term2b = line3b/omXpdot |
253 |
|
254 |
term2u1x = ChiPrime*ul2(1)*(term2a + term2b) |
255 |
term2u1y = ChiPrime*ul2(2)*(term2a + term2b) |
256 |
term2u1z = ChiPrime*ul2(3)*(term2a + term2b) |
257 |
term2u2x = ChiPrime*ul1(1)*(term2a + term2b) |
258 |
term2u2y = ChiPrime*ul1(2)*(term2a + term2b) |
259 |
term2u2z = ChiPrime*ul1(3)*(term2a + term2b) |
260 |
|
261 |
pref = -ChiPrime*0.5d0/r2 |
262 |
|
263 |
dgpdu1x = pref*(term1u1x+term2u1x) |
264 |
dgpdu1y = pref*(term1u1y+term2u1y) |
265 |
dgpdu1z = pref*(term1u1z+term2u1z) |
266 |
dgpdu2x = pref*(term1u2x+term2u2x) |
267 |
dgpdu2y = pref*(term1u2y+term2u2y) |
268 |
dgpdu2z = pref*(term1u2z+term2u2z) |
269 |
|
270 |
gp = 1.0d0 - ChiPrime*(line3a + line3b)/(2.0d0*r2) |
271 |
gmu = gp**gb_mu |
272 |
gpi = 1.0d0 / gp |
273 |
gmum = gmu*gpi |
274 |
|
275 |
! write(*,*) atom1, atom2, Chi, u1dotu2 |
276 |
curlyE = 1.0d0/dsqrt(1.0d0 - Chi*Chi*u1dotu2*u1dotu2) |
277 |
|
278 |
dcE = (curlyE**3)*Chi*Chi*u1dotu2 |
279 |
|
280 |
dcEdu1x = dcE*ul2(1) |
281 |
dcEdu1y = dcE*ul2(2) |
282 |
dcEdu1z = dcE*ul2(3) |
283 |
dcEdu2x = dcE*ul1(1) |
284 |
dcEdu2y = dcE*ul1(2) |
285 |
dcEdu2z = dcE*ul1(3) |
286 |
|
287 |
enu = curlyE**gb_nu |
288 |
enum = enu/curlyE |
289 |
|
290 |
eps = gb_eps*enu*gmu |
291 |
|
292 |
yick1 = gb_eps*enu*gb_mu*gmum |
293 |
yick2 = gb_eps*gmu*gb_nu*enum |
294 |
|
295 |
depsdu1x = yick1*dgpdu1x + yick2*dcEdu1x |
296 |
depsdu1y = yick1*dgpdu1y + yick2*dcEdu1y |
297 |
depsdu1z = yick1*dgpdu1z + yick2*dcEdu1z |
298 |
depsdu2x = yick1*dgpdu2x + yick2*dcEdu2x |
299 |
depsdu2y = yick1*dgpdu2y + yick2*dcEdu2y |
300 |
depsdu2z = yick1*dgpdu2z + yick2*dcEdu2z |
301 |
|
302 |
R126 = Ri12 - Ri6 |
303 |
R137 = 6.0d0*Ri7 - 12.0d0*Ri13 |
304 |
|
305 |
mess1 = gmu*R137 |
306 |
mess2 = R126*gb_mu*gmum |
307 |
|
308 |
dUdx = 4.0d0*gb_eps*enu*(mess1*dBigRdx + mess2*dgpdx) |
309 |
dUdy = 4.0d0*gb_eps*enu*(mess1*dBigRdy + mess2*dgpdy) |
310 |
dUdz = 4.0d0*gb_eps*enu*(mess1*dBigRdz + mess2*dgpdz) |
311 |
|
312 |
dUdu1x = 4.0d0*(R126*depsdu1x + eps*R137*dBigRdu1x) |
313 |
dUdu1y = 4.0d0*(R126*depsdu1y + eps*R137*dBigRdu1y) |
314 |
dUdu1z = 4.0d0*(R126*depsdu1z + eps*R137*dBigRdu1z) |
315 |
dUdu2x = 4.0d0*(R126*depsdu2x + eps*R137*dBigRdu2x) |
316 |
dUdu2y = 4.0d0*(R126*depsdu2y + eps*R137*dBigRdu2y) |
317 |
dUdu2z = 4.0d0*(R126*depsdu2z + eps*R137*dBigRdu2z) |
318 |
|
319 |
#ifdef IS_MPI |
320 |
f_Row(1,atom1) = f_Row(1,atom1) + dUdx |
321 |
f_Row(2,atom1) = f_Row(2,atom1) + dUdy |
322 |
f_Row(3,atom1) = f_Row(3,atom1) + dUdz |
323 |
|
324 |
f_Col(1,atom2) = f_Col(1,atom2) - dUdx |
325 |
f_Col(2,atom2) = f_Col(2,atom2) - dUdy |
326 |
f_Col(3,atom2) = f_Col(3,atom2) - dUdz |
327 |
|
328 |
t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y |
329 |
t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z |
330 |
t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x |
331 |
|
332 |
t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y |
333 |
t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z |
334 |
t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x |
335 |
#else |
336 |
f(1,atom1) = f(1,atom1) + dUdx |
337 |
f(2,atom1) = f(2,atom1) + dUdy |
338 |
f(3,atom1) = f(3,atom1) + dUdz |
339 |
|
340 |
f(1,atom2) = f(1,atom2) - dUdx |
341 |
f(2,atom2) = f(2,atom2) - dUdy |
342 |
f(3,atom2) = f(3,atom2) - dUdz |
343 |
|
344 |
t(1,atom1) = t(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y |
345 |
t(2,atom1) = t(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z |
346 |
t(3,atom1) = t(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x |
347 |
|
348 |
t(1,atom2) = t(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y |
349 |
t(2,atom2) = t(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z |
350 |
t(3,atom2) = t(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x |
351 |
#endif |
352 |
|
353 |
if (do_stress) then |
354 |
tau_Temp(1) = tau_Temp(1) + dUdx * d(1) |
355 |
tau_Temp(2) = tau_Temp(2) + dUdx * d(2) |
356 |
tau_Temp(3) = tau_Temp(3) + dUdx * d(3) |
357 |
tau_Temp(4) = tau_Temp(4) + dUdy * d(1) |
358 |
tau_Temp(5) = tau_Temp(5) + dUdy * d(2) |
359 |
tau_Temp(6) = tau_Temp(6) + dUdy * d(3) |
360 |
tau_Temp(7) = tau_Temp(7) + dUdz * d(1) |
361 |
tau_Temp(8) = tau_Temp(8) + dUdz * d(2) |
362 |
tau_Temp(9) = tau_Temp(9) + dUdz * d(3) |
363 |
virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) |
364 |
endif |
365 |
|
366 |
if (do_pot) then |
367 |
#ifdef IS_MPI |
368 |
pot_row(atom1) = pot_row(atom1) + 2.0d0*eps*R126 |
369 |
pot_col(atom2) = pot_col(atom2) + 2.0d0*eps*R126 |
370 |
#else |
371 |
pot = pot + 4.0*eps*R126 |
372 |
#endif |
373 |
endif |
374 |
|
375 |
return |
376 |
end subroutine do_gb_pair |
377 |
|
378 |
end module gb_pair |
379 |
|
380 |
|