ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 2 months ago) by gezelter
File size: 13719 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
40 !!
41
42
43 module gb_pair
44 use force_globals
45 use definitions
46 use simulation
47 #ifdef IS_MPI
48 use mpiSimulation
49 #endif
50
51 implicit none
52
53 PRIVATE
54
55 logical, save :: gb_pair_initialized = .false.
56 real(kind=dp), save :: gb_sigma
57 real(kind=dp), save :: gb_l2b_ratio
58 real(kind=dp), save :: gb_eps
59 real(kind=dp), save :: gb_eps_ratio
60 real(kind=dp), save :: gb_mu
61 real(kind=dp), save :: gb_nu
62
63 public :: check_gb_pair_FF
64 public :: set_gb_pair_params
65 public :: do_gb_pair
66
67 contains
68
69 subroutine check_gb_pair_FF(status)
70 integer :: status
71 status = -1
72 if (gb_pair_initialized) status = 0
73 return
74 end subroutine check_gb_pair_FF
75
76 subroutine set_gb_pair_params(sigma, l2b_ratio, eps, eps_ratio, mu, nu)
77 real( kind = dp ), intent(in) :: sigma, l2b_ratio, eps, eps_ratio
78 real( kind = dp ), intent(in) :: mu, nu
79
80 gb_sigma = sigma
81 gb_l2b_ratio = l2b_ratio
82 gb_eps = eps
83 gb_eps_ratio = eps_ratio
84 gb_mu = mu
85 gb_nu = nu
86
87 gb_pair_initialized = .true.
88 return
89 end subroutine set_gb_pair_params
90
91
92 subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
93 pot, A, f, t, do_pot)
94
95 integer, intent(in) :: atom1, atom2
96 integer :: id1, id2
97 real (kind=dp), intent(inout) :: r, r2
98 real (kind=dp), dimension(3), intent(in) :: d
99 real (kind=dp), dimension(3), intent(inout) :: fpair
100 real (kind=dp) :: pot, sw, vpair
101 real (kind=dp), dimension(9,nLocal) :: A
102 real (kind=dp), dimension(3,nLocal) :: f
103 real (kind=dp), dimension(3,nLocal) :: t
104 logical, intent(in) :: do_pot
105 real (kind = dp), dimension(3) :: ul1
106 real (kind = dp), dimension(3) :: ul2
107
108 real(kind=dp) :: chi, chiprime, emu, s2
109 real(kind=dp) :: r4, rdotu1, rdotu2, u1dotu2, g, gp, gpi, gmu, gmum
110 real(kind=dp) :: curlyE, enu, enum, eps, dotsum, dotdiff, ds2, dd2
111 real(kind=dp) :: opXdot, omXdot, opXpdot, omXpdot, pref, gfact
112 real(kind=dp) :: BigR, Ri, Ri2, Ri6, Ri7, Ri12, Ri13, R126, R137
113 real(kind=dp) :: dru1dx, dru1dy, dru1dz
114 real(kind=dp) :: dru2dx, dru2dy, dru2dz
115 real(kind=dp) :: dBigRdx, dBigRdy, dBigRdz
116 real(kind=dp) :: dBigRdu1x, dBigRdu1y, dBigRdu1z
117 real(kind=dp) :: dBigRdu2x, dBigRdu2y, dBigRdu2z
118 real(kind=dp) :: dUdx, dUdy, dUdz
119 real(kind=dp) :: dUdu1x, dUdu1y, dUdu1z, dUdu2x, dUdu2y, dUdu2z
120 real(kind=dp) :: dcE, dcEdu1x, dcEdu1y, dcEdu1z, dcEdu2x, dcEdu2y, dcEdu2z
121 real(kind=dp) :: depsdu1x, depsdu1y, depsdu1z, depsdu2x, depsdu2y, depsdu2z
122 real(kind=dp) :: drdx, drdy, drdz
123 real(kind=dp) :: dgdx, dgdy, dgdz
124 real(kind=dp) :: dgdu1x, dgdu1y, dgdu1z, dgdu2x, dgdu2y, dgdu2z
125 real(kind=dp) :: dgpdx, dgpdy, dgpdz
126 real(kind=dp) :: dgpdu1x, dgpdu1y, dgpdu1z, dgpdu2x, dgpdu2y, dgpdu2z
127 real(kind=dp) :: line1a, line1bx, line1by, line1bz
128 real(kind=dp) :: line2a, line2bx, line2by, line2bz
129 real(kind=dp) :: line3a, line3b, line3, line3x, line3y, line3z
130 real(kind=dp) :: term1x, term1y, term1z, term1u1x, term1u1y, term1u1z
131 real(kind=dp) :: term1u2x, term1u2y, term1u2z
132 real(kind=dp) :: term2a, term2b, term2u1x, term2u1y, term2u1z
133 real(kind=dp) :: term2u2x, term2u2y, term2u2z
134 real(kind=dp) :: yick1, yick2, mess1, mess2
135
136 s2 = (gb_l2b_ratio)**2
137 emu = (gb_eps_ratio)**(1.0d0/gb_mu)
138
139 chi = (s2 - 1.0d0)/(s2 + 1.0d0)
140 chiprime = (1.0d0 - emu)/(1.0d0 + emu)
141
142 r4 = r2*r2
143
144 #ifdef IS_MPI
145 ul1(1) = A_Row(3,atom1)
146 ul1(2) = A_Row(6,atom1)
147 ul1(3) = A_Row(9,atom1)
148
149 ul2(1) = A_Col(3,atom2)
150 ul2(2) = A_Col(6,atom2)
151 ul2(3) = A_Col(9,atom2)
152 #else
153 ul1(1) = A(3,atom1)
154 ul1(2) = A(6,atom1)
155 ul1(3) = A(9,atom1)
156
157 ul2(1) = A(3,atom2)
158 ul2(2) = A(6,atom2)
159 ul2(3) = A(9,atom2)
160 #endif
161
162 dru1dx = ul1(1)
163 dru2dx = ul2(1)
164 dru1dy = ul1(2)
165 dru2dy = ul2(2)
166 dru1dz = ul1(3)
167 dru2dz = ul2(3)
168
169 drdx = d(1) / r
170 drdy = d(2) / r
171 drdz = d(3) / r
172
173 ! do some dot products:
174 ! NB the r in these dot products is the actual intermolecular vector,
175 ! and is not the unit vector in that direction.
176
177 rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
178 rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
179 u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
180
181 ! This stuff is all for the calculation of g(Chi) and dgdx
182 ! Line numbers roughly follow the lines in equation A25 of Luckhurst
183 ! et al. Liquid Crystals 8, 451-464 (1990).
184 ! We note however, that there are some major typos in that Appendix
185 ! of the Luckhurst paper, particularly in equations A23, A29 and A31
186 ! We have attempted to correct them below.
187
188 dotsum = rdotu1+rdotu2
189 dotdiff = rdotu1-rdotu2
190 ds2 = dotsum*dotsum
191 dd2 = dotdiff*dotdiff
192
193 opXdot = 1.0d0 + Chi*u1dotu2
194 omXdot = 1.0d0 - Chi*u1dotu2
195 opXpdot = 1.0d0 + ChiPrime*u1dotu2
196 omXpdot = 1.0d0 - ChiPrime*u1dotu2
197
198 line1a = dotsum/opXdot
199 line1bx = dru1dx + dru2dx
200 line1by = dru1dy + dru2dy
201 line1bz = dru1dz + dru2dz
202
203 line2a = dotdiff/omXdot
204 line2bx = dru1dx - dru2dx
205 line2by = dru1dy - dru2dy
206 line2bz = dru1dz - dru2dz
207
208 term1x = -Chi*(line1a*line1bx + line2a*line2bx)/r2
209 term1y = -Chi*(line1a*line1by + line2a*line2by)/r2
210 term1z = -Chi*(line1a*line1bz + line2a*line2bz)/r2
211
212 line3a = ds2/opXdot
213 line3b = dd2/omXdot
214 line3 = Chi*(line3a + line3b)/r4
215 line3x = d(1)*line3
216 line3y = d(2)*line3
217 line3z = d(3)*line3
218
219 dgdx = term1x + line3x
220 dgdy = term1y + line3y
221 dgdz = term1z + line3z
222
223 term1u1x = 2.0d0*(line1a+line2a)*d(1)
224 term1u1y = 2.0d0*(line1a+line2a)*d(2)
225 term1u1z = 2.0d0*(line1a+line2a)*d(3)
226 term1u2x = 2.0d0*(line1a-line2a)*d(1)
227 term1u2y = 2.0d0*(line1a-line2a)*d(2)
228 term1u2z = 2.0d0*(line1a-line2a)*d(3)
229
230 term2a = -line3a/opXdot
231 term2b = line3b/omXdot
232
233 term2u1x = Chi*ul2(1)*(term2a + term2b)
234 term2u1y = Chi*ul2(2)*(term2a + term2b)
235 term2u1z = Chi*ul2(3)*(term2a + term2b)
236 term2u2x = Chi*ul1(1)*(term2a + term2b)
237 term2u2y = Chi*ul1(2)*(term2a + term2b)
238 term2u2z = Chi*ul1(3)*(term2a + term2b)
239
240 pref = -Chi*0.5d0/r2
241
242 dgdu1x = pref*(term1u1x+term2u1x)
243 dgdu1y = pref*(term1u1y+term2u1y)
244 dgdu1z = pref*(term1u1z+term2u1z)
245 dgdu2x = pref*(term1u2x+term2u2x)
246 dgdu2y = pref*(term1u2y+term2u2y)
247 dgdu2z = pref*(term1u2z+term2u2z)
248
249 g = 1.0d0 - Chi*(line3a + line3b)/(2.0d0*r2)
250
251 BigR = (r - gb_sigma*(g**(-0.5d0)) + gb_sigma)/gb_sigma
252 Ri = 1.0d0/BigR
253 Ri2 = Ri*Ri
254 Ri6 = Ri2*Ri2*Ri2
255 Ri7 = Ri6*Ri
256 Ri12 = Ri6*Ri6
257 Ri13 = Ri6*Ri7
258
259 gfact = (g**(-1.5d0))*0.5d0
260
261 dBigRdx = drdx/gb_sigma + dgdx*gfact
262 dBigRdy = drdy/gb_sigma + dgdy*gfact
263 dBigRdz = drdz/gb_sigma + dgdz*gfact
264 dBigRdu1x = dgdu1x*gfact
265 dBigRdu1y = dgdu1y*gfact
266 dBigRdu1z = dgdu1z*gfact
267 dBigRdu2x = dgdu2x*gfact
268 dBigRdu2y = dgdu2y*gfact
269 dBigRdu2z = dgdu2z*gfact
270
271 ! Now, we must do it again for g(ChiPrime) and dgpdx
272
273 line1a = dotsum/opXpdot
274 line2a = dotdiff/omXpdot
275 term1x = -ChiPrime*(line1a*line1bx + line2a*line2bx)/r2
276 term1y = -ChiPrime*(line1a*line1by + line2a*line2by)/r2
277 term1z = -ChiPrime*(line1a*line1bz + line2a*line2bz)/r2
278 line3a = ds2/opXpdot
279 line3b = dd2/omXpdot
280 line3 = ChiPrime*(line3a + line3b)/r4
281 line3x = d(1)*line3
282 line3y = d(2)*line3
283 line3z = d(3)*line3
284
285 dgpdx = term1x + line3x
286 dgpdy = term1y + line3y
287 dgpdz = term1z + line3z
288
289 term1u1x = 2.0d0*(line1a+line2a)*d(1)
290 term1u1y = 2.0d0*(line1a+line2a)*d(2)
291 term1u1z = 2.0d0*(line1a+line2a)*d(3)
292 term1u2x = 2.0d0*(line1a-line2a)*d(1)
293 term1u2y = 2.0d0*(line1a-line2a)*d(2)
294 term1u2z = 2.0d0*(line1a-line2a)*d(3)
295
296 term2a = -line3a/opXpdot
297 term2b = line3b/omXpdot
298
299 term2u1x = ChiPrime*ul2(1)*(term2a + term2b)
300 term2u1y = ChiPrime*ul2(2)*(term2a + term2b)
301 term2u1z = ChiPrime*ul2(3)*(term2a + term2b)
302 term2u2x = ChiPrime*ul1(1)*(term2a + term2b)
303 term2u2y = ChiPrime*ul1(2)*(term2a + term2b)
304 term2u2z = ChiPrime*ul1(3)*(term2a + term2b)
305
306 pref = -ChiPrime*0.5d0/r2
307
308 dgpdu1x = pref*(term1u1x+term2u1x)
309 dgpdu1y = pref*(term1u1y+term2u1y)
310 dgpdu1z = pref*(term1u1z+term2u1z)
311 dgpdu2x = pref*(term1u2x+term2u2x)
312 dgpdu2y = pref*(term1u2y+term2u2y)
313 dgpdu2z = pref*(term1u2z+term2u2z)
314
315 gp = 1.0d0 - ChiPrime*(line3a + line3b)/(2.0d0*r2)
316 gmu = gp**gb_mu
317 gpi = 1.0d0 / gp
318 gmum = gmu*gpi
319
320 ! write(*,*) atom1, atom2, Chi, u1dotu2
321 curlyE = 1.0d0/dsqrt(1.0d0 - Chi*Chi*u1dotu2*u1dotu2)
322
323 dcE = (curlyE**3)*Chi*Chi*u1dotu2
324
325 dcEdu1x = dcE*ul2(1)
326 dcEdu1y = dcE*ul2(2)
327 dcEdu1z = dcE*ul2(3)
328 dcEdu2x = dcE*ul1(1)
329 dcEdu2y = dcE*ul1(2)
330 dcEdu2z = dcE*ul1(3)
331
332 enu = curlyE**gb_nu
333 enum = enu/curlyE
334
335 eps = gb_eps*enu*gmu
336
337 yick1 = gb_eps*enu*gb_mu*gmum
338 yick2 = gb_eps*gmu*gb_nu*enum
339
340 depsdu1x = yick1*dgpdu1x + yick2*dcEdu1x
341 depsdu1y = yick1*dgpdu1y + yick2*dcEdu1y
342 depsdu1z = yick1*dgpdu1z + yick2*dcEdu1z
343 depsdu2x = yick1*dgpdu2x + yick2*dcEdu2x
344 depsdu2y = yick1*dgpdu2y + yick2*dcEdu2y
345 depsdu2z = yick1*dgpdu2z + yick2*dcEdu2z
346
347 R126 = Ri12 - Ri6
348 R137 = 6.0d0*Ri7 - 12.0d0*Ri13
349
350 mess1 = gmu*R137
351 mess2 = R126*gb_mu*gmum
352
353 dUdx = 4.0d0*gb_eps*enu*(mess1*dBigRdx + mess2*dgpdx)*sw
354 dUdy = 4.0d0*gb_eps*enu*(mess1*dBigRdy + mess2*dgpdy)*sw
355 dUdz = 4.0d0*gb_eps*enu*(mess1*dBigRdz + mess2*dgpdz)*sw
356
357 dUdu1x = 4.0d0*(R126*depsdu1x + eps*R137*dBigRdu1x)*sw
358 dUdu1y = 4.0d0*(R126*depsdu1y + eps*R137*dBigRdu1y)*sw
359 dUdu1z = 4.0d0*(R126*depsdu1z + eps*R137*dBigRdu1z)*sw
360 dUdu2x = 4.0d0*(R126*depsdu2x + eps*R137*dBigRdu2x)*sw
361 dUdu2y = 4.0d0*(R126*depsdu2y + eps*R137*dBigRdu2y)*sw
362 dUdu2z = 4.0d0*(R126*depsdu2z + eps*R137*dBigRdu2z)*sw
363
364 #ifdef IS_MPI
365 f_Row(1,atom1) = f_Row(1,atom1) + dUdx
366 f_Row(2,atom1) = f_Row(2,atom1) + dUdy
367 f_Row(3,atom1) = f_Row(3,atom1) + dUdz
368
369 f_Col(1,atom2) = f_Col(1,atom2) - dUdx
370 f_Col(2,atom2) = f_Col(2,atom2) - dUdy
371 f_Col(3,atom2) = f_Col(3,atom2) - dUdz
372
373 t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y
374 t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z
375 t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x
376
377 t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y
378 t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z
379 t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x
380 #else
381 f(1,atom1) = f(1,atom1) + dUdx
382 f(2,atom1) = f(2,atom1) + dUdy
383 f(3,atom1) = f(3,atom1) + dUdz
384
385 f(1,atom2) = f(1,atom2) - dUdx
386 f(2,atom2) = f(2,atom2) - dUdy
387 f(3,atom2) = f(3,atom2) - dUdz
388
389 t(1,atom1) = t(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y
390 t(2,atom1) = t(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z
391 t(3,atom1) = t(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x
392
393 t(1,atom2) = t(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y
394 t(2,atom2) = t(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z
395 t(3,atom2) = t(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x
396 #endif
397
398 if (do_pot) then
399 #ifdef IS_MPI
400 pot_row(atom1) = pot_row(atom1) + 2.0d0*eps*R126*sw
401 pot_col(atom2) = pot_col(atom2) + 2.0d0*eps*R126*sw
402 #else
403 pot = pot + 4.0*eps*R126*sw
404 #endif
405 endif
406
407 vpair = vpair + 4.0*eps*R126
408 #ifdef IS_MPI
409 id1 = AtomRowToGlobal(atom1)
410 id2 = AtomColToGlobal(atom2)
411 #else
412 id1 = atom1
413 id2 = atom2
414 #endif
415
416 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
417
418 fpair(1) = fpair(1) + dUdx
419 fpair(2) = fpair(2) + dUdy
420 fpair(3) = fpair(3) + dUdz
421
422 endif
423
424 return
425 end subroutine do_gb_pair
426
427 end module gb_pair