ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2361
Committed: Wed Oct 12 21:00:50 2005 UTC (18 years, 8 months ago) by gezelter
File size: 14314 byte(s)
Log Message:
simplifying long range potential array

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