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