ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2355
Committed: Wed Oct 12 18:59:16 2005 UTC (18 years, 8 months ago) by chuckv
File size: 14264 byte(s)
Log Message:
Breaky Breaky: Add Support for seperating potential contributions.

File Contents

# User Rev Content
1 gezelter 1930 !!
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 gezelter 1608 module gb_pair
44     use force_globals
45     use definitions
46     use simulation
47     #ifdef IS_MPI
48     use mpiSimulation
49     #endif
50 kdaily 2226
51 gezelter 1608 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 chrisfen 2277 public :: getGayBerneCut
67 gezelter 1608
68     contains
69    
70     subroutine check_gb_pair_FF(status)
71     integer :: status
72     status = -1
73     if (gb_pair_initialized) status = 0
74     return
75     end subroutine check_gb_pair_FF
76    
77     subroutine set_gb_pair_params(sigma, l2b_ratio, eps, eps_ratio, mu, nu)
78     real( kind = dp ), intent(in) :: sigma, l2b_ratio, eps, eps_ratio
79     real( kind = dp ), intent(in) :: mu, nu
80 kdaily 2226
81 gezelter 1608 gb_sigma = sigma
82     gb_l2b_ratio = l2b_ratio
83     gb_eps = eps
84     gb_eps_ratio = eps_ratio
85     gb_mu = mu
86     gb_nu = nu
87    
88     gb_pair_initialized = .true.
89     return
90     end subroutine set_gb_pair_params
91    
92 chrisfen 2279 !! gay berne cutoff should be a parameter in globals, this is a temporary
93     !! work around - this should be fixed when gay berne is up and running
94 chrisfen 2277 function getGayBerneCut(atomID) result(cutValue)
95     integer, intent(in) :: atomID !! nah... we don't need to use this...
96     real(kind=dp) :: cutValue
97 gezelter 1608
98 chrisfen 2279 cutValue = gb_l2b_ratio*gb_sigma*2.5_dp
99 chrisfen 2277 end function getGayBerneCut
100    
101 gezelter 1608 subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
102 gezelter 1930 pot, A, f, t, do_pot)
103 kdaily 2226
104 gezelter 1608 integer, intent(in) :: atom1, atom2
105     integer :: id1, id2
106     real (kind=dp), intent(inout) :: r, r2
107     real (kind=dp), dimension(3), intent(in) :: d
108     real (kind=dp), dimension(3), intent(inout) :: fpair
109     real (kind=dp) :: pot, sw, vpair
110 gezelter 1930 real (kind=dp), dimension(9,nLocal) :: A
111 gezelter 1608 real (kind=dp), dimension(3,nLocal) :: f
112     real (kind=dp), dimension(3,nLocal) :: t
113     logical, intent(in) :: do_pot
114     real (kind = dp), dimension(3) :: ul1
115     real (kind = dp), dimension(3) :: ul2
116    
117     real(kind=dp) :: chi, chiprime, emu, s2
118     real(kind=dp) :: r4, rdotu1, rdotu2, u1dotu2, g, gp, gpi, gmu, gmum
119     real(kind=dp) :: curlyE, enu, enum, eps, dotsum, dotdiff, ds2, dd2
120     real(kind=dp) :: opXdot, omXdot, opXpdot, omXpdot, pref, gfact
121     real(kind=dp) :: BigR, Ri, Ri2, Ri6, Ri7, Ri12, Ri13, R126, R137
122     real(kind=dp) :: dru1dx, dru1dy, dru1dz
123     real(kind=dp) :: dru2dx, dru2dy, dru2dz
124     real(kind=dp) :: dBigRdx, dBigRdy, dBigRdz
125     real(kind=dp) :: dBigRdu1x, dBigRdu1y, dBigRdu1z
126     real(kind=dp) :: dBigRdu2x, dBigRdu2y, dBigRdu2z
127     real(kind=dp) :: dUdx, dUdy, dUdz
128     real(kind=dp) :: dUdu1x, dUdu1y, dUdu1z, dUdu2x, dUdu2y, dUdu2z
129     real(kind=dp) :: dcE, dcEdu1x, dcEdu1y, dcEdu1z, dcEdu2x, dcEdu2y, dcEdu2z
130     real(kind=dp) :: depsdu1x, depsdu1y, depsdu1z, depsdu2x, depsdu2y, depsdu2z
131     real(kind=dp) :: drdx, drdy, drdz
132     real(kind=dp) :: dgdx, dgdy, dgdz
133     real(kind=dp) :: dgdu1x, dgdu1y, dgdu1z, dgdu2x, dgdu2y, dgdu2z
134     real(kind=dp) :: dgpdx, dgpdy, dgpdz
135     real(kind=dp) :: dgpdu1x, dgpdu1y, dgpdu1z, dgpdu2x, dgpdu2y, dgpdu2z
136     real(kind=dp) :: line1a, line1bx, line1by, line1bz
137     real(kind=dp) :: line2a, line2bx, line2by, line2bz
138     real(kind=dp) :: line3a, line3b, line3, line3x, line3y, line3z
139     real(kind=dp) :: term1x, term1y, term1z, term1u1x, term1u1y, term1u1z
140     real(kind=dp) :: term1u2x, term1u2y, term1u2z
141     real(kind=dp) :: term2a, term2b, term2u1x, term2u1y, term2u1z
142     real(kind=dp) :: term2u2x, term2u2y, term2u2z
143     real(kind=dp) :: yick1, yick2, mess1, mess2
144 kdaily 2226
145 gezelter 1608 s2 = (gb_l2b_ratio)**2
146     emu = (gb_eps_ratio)**(1.0d0/gb_mu)
147    
148     chi = (s2 - 1.0d0)/(s2 + 1.0d0)
149     chiprime = (1.0d0 - emu)/(1.0d0 + emu)
150    
151     r4 = r2*r2
152    
153     #ifdef IS_MPI
154 gezelter 1930 ul1(1) = A_Row(3,atom1)
155     ul1(2) = A_Row(6,atom1)
156     ul1(3) = A_Row(9,atom1)
157 gezelter 1608
158 gezelter 1930 ul2(1) = A_Col(3,atom2)
159     ul2(2) = A_Col(6,atom2)
160     ul2(3) = A_Col(9,atom2)
161 gezelter 1608 #else
162 gezelter 1930 ul1(1) = A(3,atom1)
163     ul1(2) = A(6,atom1)
164     ul1(3) = A(9,atom1)
165 gezelter 1608
166 gezelter 1930 ul2(1) = A(3,atom2)
167     ul2(2) = A(6,atom2)
168     ul2(3) = A(9,atom2)
169 gezelter 1608 #endif
170 kdaily 2226
171 gezelter 1608 dru1dx = ul1(1)
172     dru2dx = ul2(1)
173     dru1dy = ul1(2)
174     dru2dy = ul2(2)
175     dru1dz = ul1(3)
176     dru2dz = ul2(3)
177 kdaily 2226
178 gezelter 1608 drdx = d(1) / r
179     drdy = d(2) / r
180     drdz = d(3) / r
181 kdaily 2226
182 gezelter 1608 ! do some dot products:
183     ! NB the r in these dot products is the actual intermolecular vector,
184     ! and is not the unit vector in that direction.
185 kdaily 2226
186 gezelter 1608 rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
187     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
188     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
189    
190     ! This stuff is all for the calculation of g(Chi) and dgdx
191     ! Line numbers roughly follow the lines in equation A25 of Luckhurst
192     ! et al. Liquid Crystals 8, 451-464 (1990).
193     ! We note however, that there are some major typos in that Appendix
194     ! of the Luckhurst paper, particularly in equations A23, A29 and A31
195     ! We have attempted to correct them below.
196 kdaily 2226
197 gezelter 1608 dotsum = rdotu1+rdotu2
198     dotdiff = rdotu1-rdotu2
199     ds2 = dotsum*dotsum
200     dd2 = dotdiff*dotdiff
201 kdaily 2226
202 gezelter 1608 opXdot = 1.0d0 + Chi*u1dotu2
203     omXdot = 1.0d0 - Chi*u1dotu2
204     opXpdot = 1.0d0 + ChiPrime*u1dotu2
205     omXpdot = 1.0d0 - ChiPrime*u1dotu2
206 kdaily 2226
207 gezelter 1608 line1a = dotsum/opXdot
208     line1bx = dru1dx + dru2dx
209     line1by = dru1dy + dru2dy
210     line1bz = dru1dz + dru2dz
211 kdaily 2226
212 gezelter 1608 line2a = dotdiff/omXdot
213     line2bx = dru1dx - dru2dx
214     line2by = dru1dy - dru2dy
215     line2bz = dru1dz - dru2dz
216 kdaily 2226
217 gezelter 1608 term1x = -Chi*(line1a*line1bx + line2a*line2bx)/r2
218     term1y = -Chi*(line1a*line1by + line2a*line2by)/r2
219     term1z = -Chi*(line1a*line1bz + line2a*line2bz)/r2
220 kdaily 2226
221 gezelter 1608 line3a = ds2/opXdot
222     line3b = dd2/omXdot
223     line3 = Chi*(line3a + line3b)/r4
224     line3x = d(1)*line3
225     line3y = d(2)*line3
226     line3z = d(3)*line3
227 kdaily 2226
228 gezelter 1608 dgdx = term1x + line3x
229     dgdy = term1y + line3y
230     dgdz = term1z + line3z
231    
232 kdaily 2226 term1u1x = (line1a+line2a)*dru1dx
233     term1u1y = (line1a+line2a)*dru1dy
234     term1u1z = (line1a+line2a)*dru1dz
235     term1u2x = (line1a-line2a)*dru2dx
236     term1u2y = (line1a-line2a)*dru2dy
237     term1u2z = (line1a-line2a)*dru2dz
238    
239 gezelter 1608 term2a = -line3a/opXdot
240     term2b = line3b/omXdot
241 kdaily 2226
242 gezelter 1608 term2u1x = Chi*ul2(1)*(term2a + term2b)
243     term2u1y = Chi*ul2(2)*(term2a + term2b)
244     term2u1z = Chi*ul2(3)*(term2a + term2b)
245     term2u2x = Chi*ul1(1)*(term2a + term2b)
246     term2u2y = Chi*ul1(2)*(term2a + term2b)
247     term2u2z = Chi*ul1(3)*(term2a + term2b)
248 kdaily 2226
249 gezelter 1608 pref = -Chi*0.5d0/r2
250    
251     dgdu1x = pref*(term1u1x+term2u1x)
252     dgdu1y = pref*(term1u1y+term2u1y)
253     dgdu1z = pref*(term1u1z+term2u1z)
254     dgdu2x = pref*(term1u2x+term2u2x)
255     dgdu2y = pref*(term1u2y+term2u2y)
256     dgdu2z = pref*(term1u2z+term2u2z)
257    
258     g = 1.0d0 - Chi*(line3a + line3b)/(2.0d0*r2)
259 kdaily 2226
260 gezelter 1608 BigR = (r - gb_sigma*(g**(-0.5d0)) + gb_sigma)/gb_sigma
261     Ri = 1.0d0/BigR
262     Ri2 = Ri*Ri
263     Ri6 = Ri2*Ri2*Ri2
264     Ri7 = Ri6*Ri
265     Ri12 = Ri6*Ri6
266     Ri13 = Ri6*Ri7
267    
268     gfact = (g**(-1.5d0))*0.5d0
269    
270     dBigRdx = drdx/gb_sigma + dgdx*gfact
271     dBigRdy = drdy/gb_sigma + dgdy*gfact
272     dBigRdz = drdz/gb_sigma + dgdz*gfact
273 kdaily 2226
274 gezelter 1608 dBigRdu1x = dgdu1x*gfact
275     dBigRdu1y = dgdu1y*gfact
276     dBigRdu1z = dgdu1z*gfact
277     dBigRdu2x = dgdu2x*gfact
278     dBigRdu2y = dgdu2y*gfact
279     dBigRdu2z = dgdu2z*gfact
280 gezelter 2204
281 gezelter 1608 ! Now, we must do it again for g(ChiPrime) and dgpdx
282    
283     line1a = dotsum/opXpdot
284     line2a = dotdiff/omXpdot
285     term1x = -ChiPrime*(line1a*line1bx + line2a*line2bx)/r2
286     term1y = -ChiPrime*(line1a*line1by + line2a*line2by)/r2
287     term1z = -ChiPrime*(line1a*line1bz + line2a*line2bz)/r2
288     line3a = ds2/opXpdot
289     line3b = dd2/omXpdot
290     line3 = ChiPrime*(line3a + line3b)/r4
291     line3x = d(1)*line3
292     line3y = d(2)*line3
293     line3z = d(3)*line3
294 kdaily 2226
295 gezelter 1608 dgpdx = term1x + line3x
296     dgpdy = term1y + line3y
297     dgpdz = term1z + line3z
298 kdaily 2226
299     term1u1x = (line1a+line2a)*dru1dx
300     term1u1y = (line1a+line2a)*dru1dy
301     term1u1z = (line1a+line2a)*dru1dz
302     term1u2x = (line1a-line2a)*dru2dx
303     term1u2y = (line1a-line2a)*dru2dy
304     term1u2z = (line1a-line2a)*dru2dz
305 gezelter 2204
306 gezelter 1608 term2a = -line3a/opXpdot
307     term2b = line3b/omXpdot
308 kdaily 2226
309 gezelter 1608 term2u1x = ChiPrime*ul2(1)*(term2a + term2b)
310     term2u1y = ChiPrime*ul2(2)*(term2a + term2b)
311     term2u1z = ChiPrime*ul2(3)*(term2a + term2b)
312     term2u2x = ChiPrime*ul1(1)*(term2a + term2b)
313     term2u2y = ChiPrime*ul1(2)*(term2a + term2b)
314     term2u2z = ChiPrime*ul1(3)*(term2a + term2b)
315 kdaily 2226
316 gezelter 1608 pref = -ChiPrime*0.5d0/r2
317 kdaily 2226
318 gezelter 1608 dgpdu1x = pref*(term1u1x+term2u1x)
319     dgpdu1y = pref*(term1u1y+term2u1y)
320     dgpdu1z = pref*(term1u1z+term2u1z)
321     dgpdu2x = pref*(term1u2x+term2u2x)
322     dgpdu2y = pref*(term1u2y+term2u2y)
323     dgpdu2z = pref*(term1u2z+term2u2z)
324 kdaily 2226
325 gezelter 1608 gp = 1.0d0 - ChiPrime*(line3a + line3b)/(2.0d0*r2)
326     gmu = gp**gb_mu
327     gpi = 1.0d0 / gp
328     gmum = gmu*gpi
329 gezelter 2204
330 gezelter 1608 curlyE = 1.0d0/dsqrt(1.0d0 - Chi*Chi*u1dotu2*u1dotu2)
331    
332     dcE = (curlyE**3)*Chi*Chi*u1dotu2
333 kdaily 2226
334 gezelter 1608 dcEdu1x = dcE*ul2(1)
335     dcEdu1y = dcE*ul2(2)
336     dcEdu1z = dcE*ul2(3)
337     dcEdu2x = dcE*ul1(1)
338     dcEdu2y = dcE*ul1(2)
339     dcEdu2z = dcE*ul1(3)
340 kdaily 2226
341 gezelter 1608 enu = curlyE**gb_nu
342     enum = enu/curlyE
343 kdaily 2226
344 gezelter 1608 eps = gb_eps*enu*gmu
345    
346     yick1 = gb_eps*enu*gb_mu*gmum
347     yick2 = gb_eps*gmu*gb_nu*enum
348    
349     depsdu1x = yick1*dgpdu1x + yick2*dcEdu1x
350     depsdu1y = yick1*dgpdu1y + yick2*dcEdu1y
351     depsdu1z = yick1*dgpdu1z + yick2*dcEdu1z
352     depsdu2x = yick1*dgpdu2x + yick2*dcEdu2x
353     depsdu2y = yick1*dgpdu2y + yick2*dcEdu2y
354     depsdu2z = yick1*dgpdu2z + yick2*dcEdu2z
355 kdaily 2226
356 gezelter 1608 R126 = Ri12 - Ri6
357     R137 = 6.0d0*Ri7 - 12.0d0*Ri13
358 kdaily 2226
359 gezelter 1608 mess1 = gmu*R137
360     mess2 = R126*gb_mu*gmum
361 kdaily 2226
362 gezelter 1608 dUdx = 4.0d0*gb_eps*enu*(mess1*dBigRdx + mess2*dgpdx)*sw
363     dUdy = 4.0d0*gb_eps*enu*(mess1*dBigRdy + mess2*dgpdy)*sw
364     dUdz = 4.0d0*gb_eps*enu*(mess1*dBigRdz + mess2*dgpdz)*sw
365 kdaily 2226
366 gezelter 1608 dUdu1x = 4.0d0*(R126*depsdu1x + eps*R137*dBigRdu1x)*sw
367     dUdu1y = 4.0d0*(R126*depsdu1y + eps*R137*dBigRdu1y)*sw
368     dUdu1z = 4.0d0*(R126*depsdu1z + eps*R137*dBigRdu1z)*sw
369     dUdu2x = 4.0d0*(R126*depsdu2x + eps*R137*dBigRdu2x)*sw
370     dUdu2y = 4.0d0*(R126*depsdu2y + eps*R137*dBigRdu2y)*sw
371     dUdu2z = 4.0d0*(R126*depsdu2z + eps*R137*dBigRdu2z)*sw
372 kdaily 2226
373 gezelter 1608 #ifdef IS_MPI
374     f_Row(1,atom1) = f_Row(1,atom1) + dUdx
375     f_Row(2,atom1) = f_Row(2,atom1) + dUdy
376     f_Row(3,atom1) = f_Row(3,atom1) + dUdz
377 kdaily 2226
378 gezelter 1608 f_Col(1,atom2) = f_Col(1,atom2) - dUdx
379     f_Col(2,atom2) = f_Col(2,atom2) - dUdy
380     f_Col(3,atom2) = f_Col(3,atom2) - dUdz
381 kdaily 2226
382 gezelter 1608 t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y
383     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z
384     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x
385 kdaily 2226
386 gezelter 1608 t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y
387     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z
388     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x
389     #else
390     f(1,atom1) = f(1,atom1) + dUdx
391     f(2,atom1) = f(2,atom1) + dUdy
392     f(3,atom1) = f(3,atom1) + dUdz
393 kdaily 2226
394 gezelter 1608 f(1,atom2) = f(1,atom2) - dUdx
395     f(2,atom2) = f(2,atom2) - dUdy
396     f(3,atom2) = f(3,atom2) - dUdz
397 kdaily 2226
398 gezelter 1608 t(1,atom1) = t(1,atom1) - ul1(2)*dUdu1z + ul1(3)*dUdu1y
399     t(2,atom1) = t(2,atom1) - ul1(3)*dUdu1x + ul1(1)*dUdu1z
400     t(3,atom1) = t(3,atom1) - ul1(1)*dUdu1y + ul1(2)*dUdu1x
401 kdaily 2226
402 gezelter 1608 t(1,atom2) = t(1,atom2) - ul2(2)*dUdu2z + ul2(3)*dUdu2y
403     t(2,atom2) = t(2,atom2) - ul2(3)*dUdu2x + ul2(1)*dUdu2z
404     t(3,atom2) = t(3,atom2) - ul2(1)*dUdu2y + ul2(2)*dUdu2x
405     #endif
406 kdaily 2226
407 gezelter 1608 if (do_pot) then
408     #ifdef IS_MPI
409 chuckv 2355 pot_row(GAYBERNE_POT,atom1) = pot_row(GAYBERNE_POT,atom1) + 2.0d0*eps*R126*sw
410     pot_col(GAYBERNE_POT,atom2) = pot_col(GAYBERNE_POT,atom2) + 2.0d0*eps*R126*sw
411 gezelter 1608 #else
412     pot = pot + 4.0*eps*R126*sw
413     #endif
414     endif
415    
416     vpair = vpair + 4.0*eps*R126
417     #ifdef IS_MPI
418     id1 = AtomRowToGlobal(atom1)
419     id2 = AtomColToGlobal(atom2)
420     #else
421     id1 = atom1
422     id2 = atom2
423     #endif
424 kdaily 2226
425 gezelter 1608 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
426 kdaily 2226
427 gezelter 1608 fpair(1) = fpair(1) + dUdx
428     fpair(2) = fpair(2) + dUdy
429     fpair(3) = fpair(3) + dUdz
430 kdaily 2226
431 gezelter 1608 endif
432 kdaily 2226
433 gezelter 1608 return
434     end subroutine do_gb_pair
435    
436     end module gb_pair