ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_gb.F90
Revision: 483
Committed: Wed Apr 9 04:06:43 2003 UTC (21 years, 3 months ago) by gezelter
File size: 12054 byte(s)
Log Message:
fixes for NPT and NVT

File Contents

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