ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2277
Committed: Fri Aug 26 21:30:41 2005 UTC (19 years ago) by chrisfen
File size: 14049 byte(s)
Log Message:
added some probably nonfunctional get*cut routines

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