ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2379
Committed: Mon Oct 17 22:42:07 2005 UTC (18 years, 10 months ago) by kdaily
File size: 25167 byte(s)
Log Message:
using notation from Cleaver paper

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 gayberne
44 use force_globals
45 use definitions
46 use simulation
47 use atype_module
48 use vector_class
49 use status
50 use lj
51 #ifdef IS_MPI
52 use mpiSimulation
53 #endif
54
55 implicit none
56
57 private
58
59 #define __FORTRAN90
60 #include "UseTheForce/DarkSide/fInteractionMap.h"
61
62 public :: newGBtype
63 public :: do_gb_pair
64 public :: do_gb_lj_pair
65 public :: getGayBerneCut
66 public :: destroyGBtypes
67
68 type :: GBtype
69 integer :: atid
70 real(kind = dp ) :: sigma
71 real(kind = dp ) :: l2b_ratio
72 real(kind = dp ) :: eps
73 real(kind = dp ) :: eps_ratio
74 real(kind = dp ) :: mu
75 real(kind = dp ) :: nu
76 real(kind = dp ) :: sigma_l
77 real(kind = dp ) :: eps_l
78 end type GBtype
79
80 type, private :: GBList
81 integer :: nGBtypes = 0
82 integer :: currentGBtype = 0
83 type(GBtype), pointer :: GBtypes(:) => null()
84 integer, pointer :: atidToGBtype(:) => null()
85 end type GBList
86
87 type(GBList), save :: GBMap
88
89 contains
90
91 subroutine newGBtype(c_ident, sigma, l2b_ratio, eps, eps_ratio, mu, nu, &
92 status)
93
94 integer, intent(in) :: c_ident
95 real( kind = dp ), intent(in) :: sigma, l2b_ratio, eps, eps_ratio
96 real( kind = dp ), intent(in) :: mu, nu
97 integer, intent(out) :: status
98
99 integer :: nGBTypes, ntypes, myATID
100 integer, pointer :: MatchList(:) => null()
101 integer :: current, i
102 status = 0
103
104 if (.not.associated(GBMap%GBtypes)) then
105
106 call getMatchingElementList(atypes, "is_GayBerne", .true., &
107 nGBtypes, MatchList)
108
109 GBMap%nGBtypes = nGBtypes
110
111 allocate(GBMap%GBtypes(nGBtypes))
112
113 ntypes = getSize(atypes)
114
115 allocate(GBMap%atidToGBtype(ntypes))
116
117 !! initialize atidToGBtype to -1 so that we can use this
118 !! array to figure out which atom comes first in the GBLJ
119 !! routine
120
121 do i = 1, ntypes
122 GBMap%atidToGBtype(i) = -1
123 enddo
124
125 endif
126
127 GBMap%currentGBtype = GBMap%currentGBtype + 1
128 current = GBMap%currentGBtype
129
130 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
131 GBMap%atidToGBtype(myATID) = current
132 GBMap%GBtypes(current)%atid = myATID
133 GBMap%GBtypes(current)%sigma = sigma
134 GBMap%GBtypes(current)%l2b_ratio = l2b_ratio
135 GBMap%GBtypes(current)%eps = eps
136 GBMap%GBtypes(current)%eps_ratio = eps_ratio
137 GBMap%GBtypes(current)%mu = mu
138 GBMap%GBtypes(current)%nu = nu
139 GBMap%GBtypes(current)%sigma_l = sigma*l2b_ratio
140 GBMap%GBtypes(current)%eps_l = eps*eps_ratio
141
142 return
143 end subroutine newGBtype
144
145
146 !! gay berne cutoff should be a parameter in globals, this is a temporary
147 !! work around - this should be fixed when gay berne is up and running
148
149 function getGayBerneCut(atomID) result(cutValue)
150 integer, intent(in) :: atomID
151 integer :: gbt1
152 real(kind=dp) :: cutValue, sigma, l2b_ratio
153
154 if (GBMap%currentGBtype == 0) then
155 call handleError("GB", "No members in GBMap")
156 return
157 end if
158
159 gbt1 = GBMap%atidToGBtype(atomID)
160 sigma = GBMap%GBtypes(gbt1)%sigma
161 l2b_ratio = GBMap%GBtypes(gbt1)%l2b_ratio
162
163 cutValue = l2b_ratio*sigma*2.5_dp
164 end function getGayBerneCut
165
166 subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
167 pot, A, f, t, do_pot)
168
169 integer, intent(in) :: atom1, atom2
170 integer :: atid1, atid2, gbt1, gbt2, id1, id2
171 real (kind=dp), intent(inout) :: r, r2
172 real (kind=dp), dimension(3), intent(in) :: d
173 real (kind=dp), dimension(3), intent(inout) :: fpair
174 real (kind=dp) :: pot, sw, vpair
175 real (kind=dp), dimension(9,nLocal) :: A
176 real (kind=dp), dimension(3,nLocal) :: f
177 real (kind=dp), dimension(3,nLocal) :: t
178 logical, intent(in) :: do_pot
179 real (kind = dp), dimension(3) :: ul1
180 real (kind = dp), dimension(3) :: ul2
181
182 real(kind=dp) :: sigma, l2b_ratio, epsilon, eps_ratio, mu, nu, sigma_l, eps_l
183 real(kind=dp) :: chi, chiprime, emu, s2
184 real(kind=dp) :: r4, rdotu1, rdotu2, u1dotu2, g, gp, gpi, gmu, gmum
185 real(kind=dp) :: curlyE, enu, enum, eps, dotsum, dotdiff, ds2, dd2
186 real(kind=dp) :: opXdot, omXdot, opXpdot, omXpdot, pref, gfact
187 real(kind=dp) :: BigR, Ri, Ri2, Ri6, Ri7, Ri12, Ri13, R126, R137
188 real(kind=dp) :: dru1dx, dru1dy, dru1dz
189 real(kind=dp) :: dru2dx, dru2dy, dru2dz
190 real(kind=dp) :: dBigRdx, dBigRdy, dBigRdz
191 real(kind=dp) :: dBigRdu1x, dBigRdu1y, dBigRdu1z
192 real(kind=dp) :: dBigRdu2x, dBigRdu2y, dBigRdu2z
193 real(kind=dp) :: dUdx, dUdy, dUdz
194 real(kind=dp) :: dUdu1x, dUdu1y, dUdu1z, dUdu2x, dUdu2y, dUdu2z
195 real(kind=dp) :: dcE, dcEdu1x, dcEdu1y, dcEdu1z, dcEdu2x, dcEdu2y, dcEdu2z
196 real(kind=dp) :: depsdu1x, depsdu1y, depsdu1z, depsdu2x, depsdu2y, depsdu2z
197 real(kind=dp) :: drdx, drdy, drdz
198 real(kind=dp) :: dgdx, dgdy, dgdz
199 real(kind=dp) :: dgdu1x, dgdu1y, dgdu1z, dgdu2x, dgdu2y, dgdu2z
200 real(kind=dp) :: dgpdx, dgpdy, dgpdz
201 real(kind=dp) :: dgpdu1x, dgpdu1y, dgpdu1z, dgpdu2x, dgpdu2y, dgpdu2z
202 real(kind=dp) :: line1a, line1bx, line1by, line1bz
203 real(kind=dp) :: line2a, line2bx, line2by, line2bz
204 real(kind=dp) :: line3a, line3b, line3, line3x, line3y, line3z
205 real(kind=dp) :: term1x, term1y, term1z, term1u1x, term1u1y, term1u1z
206 real(kind=dp) :: term1u2x, term1u2y, term1u2z
207 real(kind=dp) :: term2a, term2b, term2u1x, term2u1y, term2u1z
208 real(kind=dp) :: term2u2x, term2u2y, term2u2z
209 real(kind=dp) :: yick1, yick2, mess1, mess2
210
211 #ifdef IS_MPI
212 atid1 = atid_Row(atom1)
213 atid2 = atid_Col(atom2)
214 #else
215 atid1 = atid(atom1)
216 atid2 = atid(atom2)
217 #endif
218
219 gbt1 = GBMap%atidToGBtype(atid1)
220 gbt2 = GBMap%atidToGBtype(atid2)
221
222 if (gbt1 .eq. gbt2) then
223 sigma = GBMap%GBtypes(gbt1)%sigma
224 l2b_ratio = GBMap%GBtypes(gbt1)%l2b_ratio
225 epsilon = GBMap%GBtypes(gbt1)%eps
226 eps_ratio = GBMap%GBtypes(gbt1)%eps_ratio
227 mu = GBMap%GBtypes(gbt1)%mu
228 nu = GBMap%GBtypes(gbt1)%nu
229 sigma_l = GBMap%GBtypes(gbt1)%sigma_l
230 eps_l = GBMap%GBtypes(gbt1)%eps_l
231 else
232 call handleError("GB", "GB-pair was called with two different GB types! OOPSE can only handle interactions for one GB type at a time.")
233 endif
234
235 s2 = (l2b_ratio)**2
236 emu = (eps_ratio)**(1.0d0/mu)
237
238 chi = (s2 - 1.0d0)/(s2 + 1.0d0)
239 chiprime = (1.0d0 - emu)/(1.0d0 + emu)
240
241 r4 = r2*r2
242
243 #ifdef IS_MPI
244 ul1(1) = A_Row(3,atom1)
245 ul1(2) = A_Row(6,atom1)
246 ul1(3) = A_Row(9,atom1)
247
248 ul2(1) = A_Col(3,atom2)
249 ul2(2) = A_Col(6,atom2)
250 ul2(3) = A_Col(9,atom2)
251 #else
252 ul1(1) = A(3,atom1)
253 ul1(2) = A(6,atom1)
254 ul1(3) = A(9,atom1)
255
256 ul2(1) = A(3,atom2)
257 ul2(2) = A(6,atom2)
258 ul2(3) = A(9,atom2)
259 #endif
260
261 dru1dx = ul1(1)
262 dru2dx = ul2(1)
263 dru1dy = ul1(2)
264 dru2dy = ul2(2)
265 dru1dz = ul1(3)
266 dru2dz = ul2(3)
267
268 drdx = d(1) / r
269 drdy = d(2) / r
270 drdz = d(3) / r
271
272 ! do some dot products:
273 ! NB the r in these dot products is the actual intermolecular vector,
274 ! and is not the unit vector in that direction.
275
276 rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
277 rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
278 u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
279
280 ! This stuff is all for the calculation of g(Chi) and dgdx
281 ! Line numbers roughly follow the lines in equation A25 of Luckhurst
282 ! et al. Liquid Crystals 8, 451-464 (1990).
283 ! We note however, that there are some major typos in that Appendix
284 ! of the Luckhurst paper, particularly in equations A23, A29 and A31
285 ! We have attempted to correct them below.
286
287 dotsum = rdotu1+rdotu2
288 dotdiff = rdotu1-rdotu2
289 ds2 = dotsum*dotsum
290 dd2 = dotdiff*dotdiff
291
292 opXdot = 1.0d0 + Chi*u1dotu2
293 omXdot = 1.0d0 - Chi*u1dotu2
294 opXpdot = 1.0d0 + ChiPrime*u1dotu2
295 omXpdot = 1.0d0 - ChiPrime*u1dotu2
296
297 line1a = dotsum/opXdot
298 line1bx = dru1dx + dru2dx
299 line1by = dru1dy + dru2dy
300 line1bz = dru1dz + dru2dz
301
302 line2a = dotdiff/omXdot
303 line2bx = dru1dx - dru2dx
304 line2by = dru1dy - dru2dy
305 line2bz = dru1dz - dru2dz
306
307 term1x = -Chi*(line1a*line1bx + line2a*line2bx)/r2
308 term1y = -Chi*(line1a*line1by + line2a*line2by)/r2
309 term1z = -Chi*(line1a*line1bz + line2a*line2bz)/r2
310
311 line3a = ds2/opXdot
312 line3b = dd2/omXdot
313 line3 = Chi*(line3a + line3b)/r4
314 line3x = d(1)*line3
315 line3y = d(2)*line3
316 line3z = d(3)*line3
317
318 dgdx = term1x + line3x
319 dgdy = term1y + line3y
320 dgdz = term1z + line3z
321
322 term1u1x = 2.0d0*(line1a+line2a)*d(1)
323 term1u1y = 2.0d0*(line1a+line2a)*d(2)
324 term1u1z = 2.0d0*(line1a+line2a)*d(3)
325 term1u2x = 2.0d0*(line1a-line2a)*d(1)
326 term1u2y = 2.0d0*(line1a-line2a)*d(2)
327 term1u2z = 2.0d0*(line1a-line2a)*d(3)
328
329 term2a = -line3a/opXdot
330 term2b = line3b/omXdot
331
332 term2u1x = Chi*ul2(1)*(term2a + term2b)
333 term2u1y = Chi*ul2(2)*(term2a + term2b)
334 term2u1z = Chi*ul2(3)*(term2a + term2b)
335 term2u2x = Chi*ul1(1)*(term2a + term2b)
336 term2u2y = Chi*ul1(2)*(term2a + term2b)
337 term2u2z = Chi*ul1(3)*(term2a + term2b)
338
339 pref = -Chi*0.5d0/r2
340
341 dgdu1x = pref*(term1u1x+term2u1x)
342 dgdu1y = pref*(term1u1y+term2u1y)
343 dgdu1z = pref*(term1u1z+term2u1z)
344 dgdu2x = pref*(term1u2x+term2u2x)
345 dgdu2y = pref*(term1u2y+term2u2y)
346 dgdu2z = pref*(term1u2z+term2u2z)
347
348 g = 1.0d0 - Chi*(line3a + line3b)/(2.0d0*r2)
349
350 BigR = (r - sigma*(g**(-0.5d0)) + sigma)/sigma
351 Ri = 1.0d0/BigR
352 Ri2 = Ri*Ri
353 Ri6 = Ri2*Ri2*Ri2
354 Ri7 = Ri6*Ri
355 Ri12 = Ri6*Ri6
356 Ri13 = Ri6*Ri7
357
358 gfact = (g**(-1.5d0))*0.5d0
359
360 dBigRdx = drdx/sigma + dgdx*gfact
361 dBigRdy = drdy/sigma + dgdy*gfact
362 dBigRdz = drdz/sigma + dgdz*gfact
363
364 dBigRdu1x = dgdu1x*gfact
365 dBigRdu1y = dgdu1y*gfact
366 dBigRdu1z = dgdu1z*gfact
367 dBigRdu2x = dgdu2x*gfact
368 dBigRdu2y = dgdu2y*gfact
369 dBigRdu2z = dgdu2z*gfact
370
371 ! Now, we must do it again for g(ChiPrime) and dgpdx
372
373 line1a = dotsum/opXpdot
374 line2a = dotdiff/omXpdot
375 term1x = -ChiPrime*(line1a*line1bx + line2a*line2bx)/r2
376 term1y = -ChiPrime*(line1a*line1by + line2a*line2by)/r2
377 term1z = -ChiPrime*(line1a*line1bz + line2a*line2bz)/r2
378 line3a = ds2/opXpdot
379 line3b = dd2/omXpdot
380 line3 = ChiPrime*(line3a + line3b)/r4
381 line3x = d(1)*line3
382 line3y = d(2)*line3
383 line3z = d(3)*line3
384
385 dgpdx = term1x + line3x
386 dgpdy = term1y + line3y
387 dgpdz = term1z + line3z
388
389 term1u1x = 2.00d0*(line1a+line2a)*d(1)
390 term1u1y = 2.00d0*(line1a+line2a)*d(2)
391 term1u1z = 2.00d0*(line1a+line2a)*d(3)
392 term1u2x = 2.0d0*(line1a-line2a)*d(1)
393 term1u2y = 2.0d0*(line1a-line2a)*d(2)
394 term1u2z = 2.0d0*(line1a-line2a)*d(3)
395
396 term2a = -line3a/opXpdot
397 term2b = line3b/omXpdot
398
399 term2u1x = ChiPrime*ul2(1)*(term2a + term2b)
400 term2u1y = ChiPrime*ul2(2)*(term2a + term2b)
401 term2u1z = ChiPrime*ul2(3)*(term2a + term2b)
402 term2u2x = ChiPrime*ul1(1)*(term2a + term2b)
403 term2u2y = ChiPrime*ul1(2)*(term2a + term2b)
404 term2u2z = ChiPrime*ul1(3)*(term2a + term2b)
405
406 pref = -ChiPrime*0.5d0/r2
407
408 dgpdu1x = pref*(term1u1x+term2u1x)
409 dgpdu1y = pref*(term1u1y+term2u1y)
410 dgpdu1z = pref*(term1u1z+term2u1z)
411 dgpdu2x = pref*(term1u2x+term2u2x)
412 dgpdu2y = pref*(term1u2y+term2u2y)
413 dgpdu2z = pref*(term1u2z+term2u2z)
414
415 gp = 1.0d0 - ChiPrime*(line3a + line3b)/(2.0d0*r2)
416 gmu = gp**mu
417 gpi = 1.0d0 / gp
418 gmum = gmu*gpi
419
420 curlyE = 1.0d0/dsqrt(1.0d0 - Chi*Chi*u1dotu2*u1dotu2)
421 !!$
422 !!$ dcE = -(curlyE**3)*Chi*Chi*u1dotu2
423 dcE = (curlyE**3)*Chi*Chi*u1dotu2
424
425 dcEdu1x = dcE*ul2(1)
426 dcEdu1y = dcE*ul2(2)
427 dcEdu1z = dcE*ul2(3)
428 dcEdu2x = dcE*ul1(1)
429 dcEdu2y = dcE*ul1(2)
430 dcEdu2z = dcE*ul1(3)
431
432 enu = curlyE**nu
433 enum = enu/curlyE
434
435 eps = epsilon*enu*gmu
436
437 yick1 = epsilon*enu*mu*gmum
438 yick2 = epsilon*gmu*nu*enum
439
440 depsdu1x = yick1*dgpdu1x + yick2*dcEdu1x
441 depsdu1y = yick1*dgpdu1y + yick2*dcEdu1y
442 depsdu1z = yick1*dgpdu1z + yick2*dcEdu1z
443 depsdu2x = yick1*dgpdu2x + yick2*dcEdu2x
444 depsdu2y = yick1*dgpdu2y + yick2*dcEdu2y
445 depsdu2z = yick1*dgpdu2z + yick2*dcEdu2z
446
447 R126 = Ri12 - Ri6
448 R137 = 6.0d0*Ri7 - 12.0d0*Ri13
449
450 mess1 = gmu*R137
451 mess2 = R126*mu*gmum
452
453 dUdx = 4.0d0*epsilon*enu*(mess1*dBigRdx + mess2*dgpdx)*sw
454 dUdy = 4.0d0*epsilon*enu*(mess1*dBigRdy + mess2*dgpdy)*sw
455 dUdz = 4.0d0*epsilon*enu*(mess1*dBigRdz + mess2*dgpdz)*sw
456
457 dUdu1x = 4.0d0*(R126*depsdu1x + eps*R137*dBigRdu1x)*sw
458 dUdu1y = 4.0d0*(R126*depsdu1y + eps*R137*dBigRdu1y)*sw
459 dUdu1z = 4.0d0*(R126*depsdu1z + eps*R137*dBigRdu1z)*sw
460 dUdu2x = 4.0d0*(R126*depsdu2x + eps*R137*dBigRdu2x)*sw
461 dUdu2y = 4.0d0*(R126*depsdu2y + eps*R137*dBigRdu2y)*sw
462 dUdu2z = 4.0d0*(R126*depsdu2z + eps*R137*dBigRdu2z)*sw
463
464 #ifdef IS_MPI
465 f_Row(1,atom1) = f_Row(1,atom1) + dUdx
466 f_Row(2,atom1) = f_Row(2,atom1) + dUdy
467 f_Row(3,atom1) = f_Row(3,atom1) + dUdz
468
469 f_Col(1,atom2) = f_Col(1,atom2) - dUdx
470 f_Col(2,atom2) = f_Col(2,atom2) - dUdy
471 f_Col(3,atom2) = f_Col(3,atom2) - dUdz
472
473 t_Row(1,atom1) = t_Row(1,atom1)- ul1(3)*dUdu1y + ul1(2)*dUdu1z
474 t_Row(2,atom1) = t_Row(2,atom1)- ul1(1)*dUdu1z + ul1(3)*dUdu1x
475 t_Row(3,atom1) = t_Row(3,atom1)- ul1(2)*dUdu1x + ul1(1)*dUdu1y
476
477 t_Col(1,atom2) = t_Col(1,atom2) - ul2(3)*dUdu2y + ul2(2)*dUdu2z
478 t_Col(2,atom2) = t_Col(2,atom2) - ul2(1)*dUdu2z + ul2(3)*dUdu2x
479 t_Col(3,atom2) = t_Col(3,atom2) - ul2(2)*dUdu2x + ul2(1)*dUdu2y
480 #else
481 f(1,atom1) = f(1,atom1) + dUdx
482 f(2,atom1) = f(2,atom1) + dUdy
483 f(3,atom1) = f(3,atom1) + dUdz
484
485 f(1,atom2) = f(1,atom2) - dUdx
486 f(2,atom2) = f(2,atom2) - dUdy
487 f(3,atom2) = f(3,atom2) - dUdz
488
489 t(1,atom1) = t(1,atom1)- ul1(3)*dUdu1y + ul1(2)*dUdu1z
490 t(2,atom1) = t(2,atom1)- ul1(1)*dUdu1z + ul1(3)*dUdu1x
491 t(3,atom1) = t(3,atom1)- ul1(2)*dUdu1x + ul1(1)*dUdu1y
492
493 t(1,atom2) = t(1,atom2)- ul2(3)*dUdu2y + ul2(2)*dUdu2z
494 t(2,atom2) = t(2,atom2)- ul2(1)*dUdu2z + ul2(3)*dUdu2x
495 t(3,atom2) = t(3,atom2)- ul2(2)*dUdu2x + ul2(1)*dUdu2y
496 #endif
497
498 if (do_pot) then
499 #ifdef IS_MPI
500 pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 2.0d0*eps*R126*sw
501 pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 2.0d0*eps*R126*sw
502 #else
503 pot = pot + 4.0*eps*R126*sw
504 #endif
505 endif
506
507 vpair = vpair + 4.0*eps*R126
508 #ifdef IS_MPI
509 id1 = AtomRowToGlobal(atom1)
510 id2 = AtomColToGlobal(atom2)
511 #else
512 id1 = atom1
513 id2 = atom2
514 #endif
515
516 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
517
518 fpair(1) = fpair(1) + dUdx
519 fpair(2) = fpair(2) + dUdy
520 fpair(3) = fpair(3) + dUdz
521
522 endif
523
524 return
525 end subroutine do_gb_pair
526
527 subroutine do_gb_lj_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
528 pot, A, f, t, do_pot)
529
530 integer, intent(in) :: atom1, atom2
531 integer :: id1, id2
532 real (kind=dp), intent(inout) :: r, r2
533 real (kind=dp), dimension(3), intent(in) :: d
534 real (kind=dp), dimension(3), intent(inout) :: fpair
535 real (kind=dp) :: pot, sw, vpair
536 real (kind=dp), dimension(9,nLocal) :: A
537 real (kind=dp), dimension(3,nLocal) :: f
538 real (kind=dp), dimension(3,nLocal) :: t
539 logical, intent(in) :: do_pot
540 real (kind = dp), dimension(3) :: ul
541
542 real(kind=dp) :: gb_sigma, gb_eps, gb_eps_ratio, gb_mu, gb_l2b_ratio
543 real(kind=dp) :: s0, l2, d2, lj2
544 real(kind=dp) :: eE, eS, eab, eabf, moom, mum1
545 real(kind=dp) :: dx, dy, dz, drdx, drdy, drdz, rdotu
546 real(kind=dp) :: mess, sab, dsabdct, depmudct
547 real(kind=dp) :: epmu, depmudx, depmudy, depmudz
548 real(kind=dp) :: depmudux, depmuduy, depmuduz
549 real(kind=dp) :: BigR, dBigRdx, dBigRdy, dBigRdz
550 real(kind=dp) :: dBigRdux, dBigRduy, dBigRduz
551 real(kind=dp) :: dUdx, dUdy, dUdz, dUdux, dUduy, dUduz, e0
552 real(kind=dp) :: Ri, Ri3, Ri6, Ri7, Ri12, Ri13, R126, R137, prefactor
553 real(kind=dp) :: chipoalphap2, chioalpha2, ec, epsnot
554 real(kind=dp) :: drdotudx, drdotudy, drdotudz
555 real(kind=dp) :: drdotudux, drdotuduy, drdotuduz
556 real(kind=dp) :: ljeps, ljsigma
557 integer :: ljt1, ljt2, atid1, atid2, gbt1, gbt2
558 logical :: gb_first
559
560 #ifdef IS_MPI
561 atid1 = atid_Row(atom1)
562 atid2 = atid_Col(atom2)
563 #else
564 atid1 = atid(atom1)
565 atid2 = atid(atom2)
566 #endif
567
568 gbt1 = GBMap%atidToGBtype(atid1)
569 gbt2 = GBMap%atidToGBtype(atid2)
570
571 if (gbt1 .eq. -1) then
572 gb_first = .false.
573 if (gbt2 .eq. -1) then
574 call handleError("GB", "GBLJ was called without a GB type.")
575 endif
576 else
577 gb_first = .true.
578 if (gbt2 .ne. -1) then
579 call handleError("GB", "GBLJ was called with two GB types (instead of one).")
580 endif
581 endif
582
583 ri =1/r
584
585 dx = d(1)
586 dy = d(2)
587 dz = d(3)
588
589 drdx = dx *ri
590 drdy = dy *ri
591 drdz = dz *ri
592
593 if(gb_first)then
594 #ifdef IS_MPI
595 ul(1) = A_Row(3,atom1)
596 ul(2) = A_Row(6,atom1)
597 ul(3) = A_Row(9,atom1)
598 #else
599 ul(1) = A(3,atom1)
600 ul(2) = A(6,atom1)
601 ul(3) = A(9,atom1)
602 #endif
603 gb_sigma = GBMap%GBtypes(gbt1)%sigma
604 gb_l2b_ratio = GBMap%GBtypes(gbt1)%l2b_ratio
605 gb_eps = GBMap%GBtypes(gbt1)%eps
606 gb_eps_ratio = GBMap%GBtypes(gbt1)%eps_ratio
607 gb_mu = GBMap%GBtypes(gbt1)%mu
608
609 ljsigma = getSigma(atid2)
610 ljeps = getEpsilon(atid2)
611 else
612 #ifdef IS_MPI
613 ul(1) = A_Col(3,atom2)
614 ul(2) = A_Col(6,atom2)
615 ul(3) = A_Col(9,atom2)
616 #else
617 ul(1) = A(3,atom2)
618 ul(2) = A(6,atom2)
619 ul(3) = A(9,atom2)
620 #endif
621 gb_sigma = GBMap%GBtypes(gbt2)%sigma
622 gb_l2b_ratio = GBMap%GBtypes(gbt2)%l2b_ratio
623 gb_eps = GBMap%GBtypes(gbt2)%eps
624 gb_eps_ratio = GBMap%GBtypes(gbt2)%eps_ratio
625 gb_mu = GBMap%GBtypes(gbt2)%mu
626
627 ljsigma = getSigma(atid1)
628 ljeps = getEpsilon(atid1)
629 endif
630
631 write(*,*) 'd u', dx, dy, dz, ul(1), ul(2), ul(3)
632
633 rdotu = (dx*ul(1)+dy*ul(2)+dz*ul(3))*ri
634
635 drdotudx = ul(1)*ri-rdotu*dx*ri*ri
636 drdotudy = ul(2)*ri-rdotu*dy*ri*ri
637 drdotudz = ul(3)*ri-rdotu*dz*ri*ri
638 drdotudux = drdx
639 drdotuduy = drdy
640 drdotuduz = drdz
641
642 l2 = (gb_sigma*gb_l2b_ratio)**2
643 d2 = gb_sigma**2
644 lj2 = ljsigma**2
645 s0 = dsqrt(d2 + lj2)
646
647 chioalpha2 = (l2 - d2)/(l2 + lj2)
648
649 eE = dsqrt(gb_eps*ljeps)
650 eS = dsqrt(gb_eps*gb_eps_ratio*ljeps)
651 moom = 1.0d0 / gb_mu
652 mum1 = gb_mu-1
653 chipoalphap2 = 1 - (eE/eS)**moom
654
655 !! mess matches cleaver (eq 20)
656
657 mess = 1-rdotu*rdotu*chioalpha2
658 sab = 1.0d0/dsqrt(mess)
659
660 write(*,*) 's', s0, sab, rdotu, chioalpha2
661 dsabdct = s0*sab*sab*sab*rdotu*chioalpha2
662
663 eab = 1-chipoalphap2*rdotu*rdotu
664 eabf = eS*(eab**gb_mu)
665
666 write(*,*) 'e', eS, chipoalphap2, gb_mu, rdotu, eab, mum1
667
668 depmudct = -2*eS*chipoalphap2*gb_mu*rdotu*(eab**mum1)
669
670 BigR = (r - sab*s0 + s0)/s0
671 dBigRdx = (drdx -dsabdct*drdotudx)/s0
672 dBigRdy = (drdy -dsabdct*drdotudy)/s0
673 dBigRdz = (drdz -dsabdct*drdotudz)/s0
674 dBigRdux = (-dsabdct*drdotudux)/s0
675 dBigRduy = (-dsabdct*drdotuduy)/s0
676 dBigRduz = (-dsabdct*drdotuduz)/s0
677
678 write(*,*) 'ds dep', dsabdct, depmudct
679 write(*,*) 'drdotudu', drdotudux, drdotuduy, drdotuduz
680 depmudx = depmudct*drdotudx
681 depmudy = depmudct*drdotudy
682 depmudz = depmudct*drdotudz
683 depmudux = depmudct*drdotudux
684 depmuduy = depmudct*drdotuduy
685 depmuduz = depmudct*drdotuduz
686
687 Ri = 1.0d0/BigR
688 Ri3 = Ri*Ri*Ri
689 Ri6 = Ri3*Ri3
690 Ri7 = Ri6*Ri
691 Ri12 = Ri6*Ri6
692 Ri13 = Ri6*Ri7
693 R126 = Ri12 - Ri6
694 R137 = 6.0d0*Ri7 - 12.0d0*Ri13
695
696 prefactor = 4.0d0
697
698 dUdx = prefactor*(eabf*R137*dBigRdx + R126*depmudx)
699 dUdy = prefactor*(eabf*R137*dBigRdy + R126*depmudy)
700 dUdz = prefactor*(eabf*R137*dBigRdz + R126*depmudz)
701 write(*,*) 'dRdu', dbigrdux, dbigrduy, dbigrduz
702 write(*,*) 'dEdu', depmudux, depmuduy, depmuduz
703 dUdux = prefactor*(eabf*R137*dBigRdux + R126*depmudux)
704 dUduy = prefactor*(eabf*R137*dBigRduy + R126*depmuduy)
705 dUduz = prefactor*(eabf*R137*dBigRduz + R126*depmuduz)
706
707 #ifdef IS_MPI
708 f_Row(1,atom1) = f_Row(1,atom1) - dUdx
709 f_Row(2,atom1) = f_Row(2,atom1) - dUdy
710 f_Row(3,atom1) = f_Row(3,atom1) - dUdz
711
712 f_Col(1,atom2) = f_Col(1,atom2) + dUdx
713 f_Col(2,atom2) = f_Col(2,atom2) + dUdy
714 f_Col(3,atom2) = f_Col(3,atom2) + dUdz
715
716 if (gb_first) then
717 t_Row(1,atom1) = t_Row(1,atom1) + ul(2)*dUduz - ul(3)*dUduy
718 t_Row(2,atom1) = t_Row(2,atom1) + ul(3)*dUdux - ul(1)*dUduz
719 t_Row(3,atom1) = t_Row(3,atom1) + ul(1)*dUduy - ul(2)*dUdux
720 else
721 t_Col(1,atom2) = t_Col(1,atom2) + ul(2)*dUduz - ul(3)*dUduy
722 t_Col(2,atom2) = t_Col(2,atom2) + ul(3)*dUdux - ul(1)*dUduz
723 t_Col(3,atom2) = t_Col(3,atom2) + ul(1)*dUduy - ul(2)*dUdux
724 endif
725 #else
726 f(1,atom1) = f(1,atom1) + dUdx
727 f(2,atom1) = f(2,atom1) + dUdy
728 f(3,atom1) = f(3,atom1) + dUdz
729
730 f(1,atom2) = f(1,atom2) - dUdx
731 f(2,atom2) = f(2,atom2) - dUdy
732 f(3,atom2) = f(3,atom2) - dUdz
733
734 ! torques are cross products:
735
736 write(*,*) 'dU', dUdux, duduy, duduz
737
738 if (gb_first) then
739 t(1,atom1) = t(1,atom1) + ul(2)*dUduz - ul(3)*dUduy
740 t(2,atom1) = t(2,atom1) + ul(3)*dUdux - ul(1)*dUduz
741 t(3,atom1) = t(3,atom1) + ul(1)*dUduy - ul(2)*dUdux
742 write(*,*) 'T1', t(1,atom1), t(2,atom1), t(3,atom1)
743 else
744 t(1,atom2) = t(1,atom2) + ul(2)*dUduz - ul(3)*dUduy
745 t(2,atom2) = t(2,atom2) + ul(3)*dUdux - ul(1)*dUduz
746 t(3,atom2) = t(3,atom2) + ul(1)*dUduy - ul(2)*dUdux
747
748 write(*,*) 'T2', t(1,atom2), t(2,atom2), t(3,atom2)
749 endif
750
751 #endif
752
753 if (do_pot) then
754 #ifdef IS_MPI
755 pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 2.0d0*eps*R126*sw
756 pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 2.0d0*eps*R126*sw
757 #else
758 pot = pot + prefactor*eabf*R126*sw
759 #endif
760 endif
761
762 vpair = vpair + 4.0*epmu*R126
763 #ifdef IS_MPI
764 id1 = AtomRowToGlobal(atom1)
765 id2 = AtomColToGlobal(atom2)
766 #else
767 id1 = atom1
768 id2 = atom2
769 #endif
770
771 If (Molmembershiplist(Id1) .Ne. Molmembershiplist(Id2)) Then
772
773 Fpair(1) = Fpair(1) + Dudx
774 Fpair(2) = Fpair(2) + Dudy
775 Fpair(3) = Fpair(3) + Dudz
776
777 Endif
778
779 return
780
781 end subroutine do_gb_lj_pair
782
783 subroutine destroyGBTypes()
784
785 GBMap%nGBtypes = 0
786 GBMap%currentGBtype = 0
787
788 if (associated(GBMap%GBtypes)) then
789 deallocate(GBMap%GBtypes)
790 GBMap%GBtypes => null()
791 end if
792
793 if (associated(GBMap%atidToGBtype)) then
794 deallocate(GBMap%atidToGBtype)
795 GBMap%atidToGBtype => null()
796 end if
797
798 end subroutine destroyGBTypes
799
800 end module gayberne
801