ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
Revision: 2958
Committed: Mon Jul 24 14:51:09 2006 UTC (18 years, 1 month ago) by xsun
File size: 16087 byte(s)
Log Message:
fixed a getSigma bug

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 linearalgebra
50 use status
51 use lj
52 use fForceOptions
53 #ifdef IS_MPI
54 use mpiSimulation
55 #endif
56
57 implicit none
58
59 private
60
61 #define __FORTRAN90
62 #include "UseTheForce/DarkSide/fInteractionMap.h"
63
64 logical, save :: haveGBMap = .false.
65 logical, save :: haveMixingMap = .false.
66 real(kind=dp), save :: mu = 2.0_dp
67 real(kind=dp), save :: nu = 1.0_dp
68
69
70 public :: newGBtype
71 public :: complete_GB_FF
72 public :: do_gb_pair
73 public :: getGayBerneCut
74 public :: destroyGBtypes
75
76 type :: GBtype
77 integer :: atid
78 real(kind = dp ) :: d
79 real(kind = dp ) :: l
80 real(kind = dp ) :: eps
81 real(kind = dp ) :: eps_ratio
82 real(kind = dp ) :: dw
83 logical :: isLJ
84 end type GBtype
85
86 type, private :: GBList
87 integer :: nGBtypes = 0
88 integer :: currentGBtype = 0
89 type(GBtype), pointer :: GBtypes(:) => null()
90 integer, pointer :: atidToGBtype(:) => null()
91 end type GBList
92
93 type(GBList), save :: GBMap
94
95 type :: GBMixParameters
96 real(kind=DP) :: sigma0
97 real(kind=DP) :: eps0
98 real(kind=DP) :: dw
99 real(kind=DP) :: x2
100 real(kind=DP) :: xa2
101 real(kind=DP) :: xai2
102 real(kind=DP) :: xp2
103 real(kind=DP) :: xpap2
104 real(kind=DP) :: xpapi2
105 end type GBMixParameters
106
107 type(GBMixParameters), dimension(:,:), allocatable :: GBMixingMap
108
109 contains
110
111 subroutine newGBtype(c_ident, d, l, eps, eps_ratio, dw, status)
112
113 integer, intent(in) :: c_ident
114 real( kind = dp ), intent(in) :: d, l, eps, eps_ratio, dw
115 integer, intent(out) :: status
116
117 integer :: nGBTypes, nLJTypes, ntypes, myATID
118 integer, pointer :: MatchList(:) => null()
119 integer :: current, i
120 status = 0
121
122 if (.not.associated(GBMap%GBtypes)) then
123
124 call getMatchingElementList(atypes, "is_GayBerne", .true., &
125 nGBtypes, MatchList)
126
127 call getMatchingElementList(atypes, "is_LennardJones", .true., &
128 nLJTypes, MatchList)
129
130 GBMap%nGBtypes = nGBtypes + nLJTypes
131
132 allocate(GBMap%GBtypes(nGBtypes + nLJTypes))
133
134 ntypes = getSize(atypes)
135
136 allocate(GBMap%atidToGBtype(ntypes))
137 endif
138
139 GBMap%currentGBtype = GBMap%currentGBtype + 1
140 current = GBMap%currentGBtype
141
142 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
143
144 GBMap%atidToGBtype(myATID) = current
145 GBMap%GBtypes(current)%atid = myATID
146 GBMap%GBtypes(current)%d = d
147 GBMap%GBtypes(current)%l = l
148 GBMap%GBtypes(current)%eps = eps
149 GBMap%GBtypes(current)%eps_ratio = eps_ratio
150 GBMap%GBtypes(current)%dw = dw
151 GBMap%GBtypes(current)%isLJ = .false.
152
153 return
154 end subroutine newGBtype
155
156 subroutine complete_GB_FF(status)
157 integer :: status
158 integer :: i, j, l, m, lm, function_type
159 real(kind=dp) :: thisDP, sigma
160 integer :: alloc_stat, iTheta, iPhi, nSteps, nAtypes, myATID, current
161 logical :: thisProperty
162
163 status = 0
164 if (GBMap%currentGBtype == 0) then
165 call handleError("complete_GB_FF", "No members in GBMap")
166 status = -1
167 return
168 end if
169
170 nAtypes = getSize(atypes)
171
172 if (nAtypes == 0) then
173 status = -1
174 return
175 end if
176
177 ! atypes comes from c side
178 do i = 1, nAtypes
179
180 myATID = getFirstMatchingElement(atypes, 'c_ident', i)
181 call getElementProperty(atypes, myATID, "is_LennardJones", thisProperty)
182
183 if (thisProperty) then
184 GBMap%currentGBtype = GBMap%currentGBtype + 1
185 current = GBMap%currentGBtype
186
187 GBMap%atidToGBtype(myATID) = current
188 GBMap%GBtypes(current)%atid = myATID
189 GBMap%GBtypes(current)%isLJ = .true.
190 GBMap%GBtypes(current)%d = getSigma(myATID) / sqrt(2.0_dp)
191 GBMap%GBtypes(current)%l = GBMap%GBtypes(current)%d
192 GBMap%GBtypes(current)%eps = getEpsilon(myATID)
193 GBMap%GBtypes(current)%eps_ratio = 1.0_dp
194 GBMap%GBtypes(current)%dw = 1.0_dp
195
196 endif
197
198 end do
199
200 haveGBMap = .true.
201
202
203 end subroutine complete_GB_FF
204
205 subroutine createGBMixingMap()
206 integer :: nGBtypes, i, j
207 real (kind = dp) :: d1, l1, e1, er1, dw1
208 real (kind = dp) :: d2, l2, e2, er2, dw2
209 real (kind = dp) :: er, ermu, xp, ap2
210
211 if (GBMap%currentGBtype == 0) then
212 call handleError("GB", "No members in GBMap")
213 return
214 end if
215
216 nGBtypes = GBMap%nGBtypes
217
218 if (.not. allocated(GBMixingMap)) then
219 allocate(GBMixingMap(nGBtypes, nGBtypes))
220 endif
221
222 do i = 1, nGBtypes
223
224 d1 = GBMap%GBtypes(i)%d
225 l1 = GBMap%GBtypes(i)%l
226 e1 = GBMap%GBtypes(i)%eps
227 er1 = GBMap%GBtypes(i)%eps_ratio
228 dw1 = GBMap%GBtypes(i)%dw
229
230 do j = 1, nGBtypes
231
232 d2 = GBMap%GBtypes(j)%d
233 l2 = GBMap%GBtypes(j)%l
234 e2 = GBMap%GBtypes(j)%eps
235 er2 = GBMap%GBtypes(j)%eps_ratio
236 dw2 = GBMap%GBtypes(j)%dw
237
238 ! Cleaver paper uses sqrt of squares to get sigma0 for
239 ! mixed interactions.
240
241 GBMixingMap(i,j)%sigma0 = sqrt(d1*d1 + d2*d2)
242 GBMixingMap(i,j)%xa2 = (l1*l1 - d1*d1)/(l1*l1 + d2*d2)
243 GBMixingMap(i,j)%xai2 = (l2*l2 - d2*d2)/(l2*l2 + d1*d1)
244 GBMixingMap(i,j)%x2 = (l1*l1 - d1*d1) * (l2*l2 - d2*d2) / &
245 ((l2*l2 + d1*d1) * (l1*l1 + d2*d2))
246
247 ! assumed LB mixing rules for now:
248
249 GBMixingMap(i,j)%dw = 0.5_dp * (dw1 + dw2)
250 GBMixingMap(i,j)%eps0 = sqrt(e1 * e2)
251
252 er = sqrt(er1 * er2)
253 ermu = er**(1.0_dp / mu)
254 xp = (1.0_dp - ermu) / (1.0_dp + ermu)
255 ap2 = 1.0_dp / (1.0_dp + ermu)
256
257 GBMixingMap(i,j)%xp2 = xp*xp
258 GBMixingMap(i,j)%xpap2 = xp*ap2
259 GBMixingMap(i,j)%xpapi2 = xp/ap2
260 enddo
261 enddo
262 haveMixingMap = .true.
263 mu = getGayBerneMu()
264 nu = getGayBerneNu()
265 end subroutine createGBMixingMap
266
267
268 !! gay berne cutoff should be a parameter in globals, this is a temporary
269 !! work around - this should be fixed when gay berne is up and running
270
271 function getGayBerneCut(atomID) result(cutValue)
272 integer, intent(in) :: atomID
273 integer :: gbt1
274 real(kind=dp) :: cutValue, l, d
275
276 if (GBMap%currentGBtype == 0) then
277 call handleError("GB", "No members in GBMap")
278 return
279 end if
280
281 gbt1 = GBMap%atidToGBtype(atomID)
282 l = GBMap%GBtypes(gbt1)%l
283 d = GBMap%GBtypes(gbt1)%d
284
285 ! sigma is actually sqrt(2)*l for prolate ellipsoids
286
287 cutValue = 2.5_dp*sqrt(2.0_dp)*max(l,d)
288
289 end function getGayBerneCut
290
291 subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
292 pot, Amat, f, t, do_pot)
293
294 integer, intent(in) :: atom1, atom2
295 integer :: atid1, atid2, gbt1, gbt2, id1, id2
296 real (kind=dp), intent(inout) :: r, r2
297 real (kind=dp), dimension(3), intent(in) :: d
298 real (kind=dp), dimension(3), intent(inout) :: fpair
299 real (kind=dp) :: pot, sw, vpair
300 real (kind=dp), dimension(9,nLocal) :: Amat
301 real (kind=dp), dimension(3,nLocal) :: f
302 real (kind=dp), dimension(3,nLocal) :: t
303 logical, intent(in) :: do_pot
304 real (kind = dp), dimension(3) :: ul1, ul2, rxu1, rxu2, uxu, rhat
305
306 real (kind = dp) :: sigma0, dw, eps0, x2, xa2, xai2, xp2, xpap2, xpapi2
307 real (kind = dp) :: e1, e2, eps, sigma, s3, s03, au2, bu2, au, bu, a, b, g, g2
308 real (kind = dp) :: U, BigR, R3, R6, R7, R12, R13, H, Hp, fx, fy, fz
309 real (kind = dp) :: dUdr, dUda, dUdb, dUdg, pref1, pref2
310 logical :: i_is_lj, j_is_lj
311
312 if (.not.haveMixingMap) then
313 call createGBMixingMap()
314 endif
315
316 #ifdef IS_MPI
317 atid1 = atid_Row(atom1)
318 atid2 = atid_Col(atom2)
319 #else
320 atid1 = atid(atom1)
321 atid2 = atid(atom2)
322 #endif
323
324 gbt1 = GBMap%atidToGBtype(atid1)
325 gbt2 = GBMap%atidToGBtype(atid2)
326
327 i_is_LJ = GBMap%GBTypes(gbt1)%isLJ
328 j_is_LJ = GBMap%GBTypes(gbt2)%isLJ
329
330 sigma0 = GBMixingMap(gbt1, gbt2)%sigma0
331 dw = GBMixingMap(gbt1, gbt2)%dw
332 eps0 = GBMixingMap(gbt1, gbt2)%eps0
333 x2 = GBMixingMap(gbt1, gbt2)%x2
334 xa2 = GBMixingMap(gbt1, gbt2)%xa2
335 xai2 = GBMixingMap(gbt1, gbt2)%xai2
336 xp2 = GBMixingMap(gbt1, gbt2)%xp2
337 xpap2 = GBMixingMap(gbt1, gbt2)%xpap2
338 xpapi2 = GBMixingMap(gbt1, gbt2)%xpapi2
339
340 #ifdef IS_MPI
341 ul1(1) = A_Row(7,atom1)
342 ul1(2) = A_Row(8,atom1)
343 ul1(3) = A_Row(9,atom1)
344
345 ul2(1) = A_Col(7,atom2)
346 ul2(2) = A_Col(8,atom2)
347 ul2(3) = A_Col(9,atom2)
348 #else
349 ul1(1) = Amat(7,atom1)
350 ul1(2) = Amat(8,atom1)
351 ul1(3) = Amat(9,atom1)
352
353 ul2(1) = Amat(7,atom2)
354 ul2(2) = Amat(8,atom2)
355 ul2(3) = Amat(9,atom2)
356 #endif
357
358 if (i_is_LJ) then
359 a = 0.0_dp
360 ul1 = 0.0_dp
361 else
362 a = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
363 endif
364
365 if (j_is_LJ) then
366 b = 0.0_dp
367 ul2 = 0.0_dp
368 else
369 b = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
370 endif
371
372 if (i_is_LJ.or.j_is_LJ) then
373 g = 0.0_dp
374 else
375 g = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
376 endif
377
378 au = a / r
379 bu = b / r
380
381 au2 = au * au
382 bu2 = bu * bu
383 g2 = g * g
384
385 H = (xa2 * au2 + xai2 * bu2 - 2.0_dp*x2*au*bu*g) / (1.0_dp - x2*g2)
386 Hp = (xpap2*au2 + xpapi2*bu2 - 2.0_dp*xp2*au*bu*g) / (1.0_dp - xp2*g2)
387
388 sigma = sigma0 / sqrt(1.0_dp - H)
389 e1 = 1.0_dp / sqrt(1.0_dp - x2*g2)
390 e2 = 1.0_dp - Hp
391 eps = eps0 * (e1**nu) * (e2**mu)
392 BigR = dw*sigma0 / (r - sigma + dw*sigma0)
393
394 R3 = BigR*BigR*BigR
395 R6 = R3*R3
396 R7 = R6 * BigR
397 R12 = R6*R6
398 R13 = R6*R7
399
400 U = 4.0_dp * eps * (R12 - R6)
401
402 s3 = sigma*sigma*sigma
403 s03 = sigma0*sigma0*sigma0
404
405 pref1 = - 8.0_dp * eps * mu * (R12 - R6) / (e2 * r)
406
407 pref2 = 8.0_dp * eps * s3 * (6.0_dp*R13 - 3.0_dp*R7) / (dw*r*s03)
408
409 dUdr = - (pref1 * Hp + pref2 * (sigma0*sigma0*r/s3 + H))
410
411 dUda = pref1 * (xpap2*au - xp2*bu*g) / (1.0_dp - xp2 * g2) &
412 + pref2 * (xa2 * au - x2 *bu*g) / (1.0_dp - x2 * g2)
413
414 dUdb = pref1 * (xpapi2*bu - xp2*au*g) / (1.0_dp - xp2 * g2) &
415 + pref2 * (xai2 * bu - x2 *au*g) / (1.0_dp - x2 * g2)
416
417 dUdg = 4.0_dp * eps * nu * (R12 - R6) * x2 * g / (1.0_dp - x2*g2) &
418 + 8.0_dp * eps * mu * (R12 - R6) * (xp2*au*bu - Hp*xp2*g) / &
419 (1.0_dp - xp2 * g2) / e2 &
420 + 8.0_dp * eps * s3 * (3.0_dp * R7 - 6.0_dp * R13) * &
421 (x2 * au * bu - H * x2 * g) / (1.0_dp - x2 * g2) / (dw * s03)
422
423
424 rhat = d / r
425
426 fx = dUdr * rhat(1) + dUda * ul1(1) + dUdb * ul2(1)
427 fy = dUdr * rhat(2) + dUda * ul1(2) + dUdb * ul2(2)
428 fz = dUdr * rhat(3) + dUda * ul1(3) + dUdb * ul2(3)
429
430 rxu1 = cross_product(d, ul1)
431 rxu2 = cross_product(d, ul2)
432 uxu = cross_product(ul1, ul2)
433
434 !!$ write(*,*) 'pref = ' , pref1, pref2
435 !!$ write(*,*) 'rxu1 = ' , rxu1(1), rxu1(2), rxu1(3)
436 !!$ write(*,*) 'rxu2 = ' , rxu2(1), rxu2(2), rxu2(3)
437 !!$ write(*,*) 'uxu = ' , uxu(1), uxu(2), uxu(3)
438 !!$ write(*,*) 'dUda = ', dUda, dudb, dudg, dudr
439 !!$ write(*,*) 'H = ', H,hp,sigma, e1, e2, BigR
440 !!$ write(*,*) 'chi = ', xa2, xai2, x2
441 !!$ write(*,*) 'chip = ', xpap2, xpapi2, xp2
442 !!$ write(*,*) 'eps = ', eps0, e1, e2, eps
443 !!$ write(*,*) 'U =', U, pref1, pref2
444 !!$ write(*,*) 'f =', fx, fy, fz
445 !!$ write(*,*) 'au =', au, bu, g
446 !!$
447
448
449 #ifdef IS_MPI
450 f_Row(1,atom1) = f_Row(1,atom1) + fx
451 f_Row(2,atom1) = f_Row(2,atom1) + fy
452 f_Row(3,atom1) = f_Row(3,atom1) + fz
453
454 f_Col(1,atom2) = f_Col(1,atom2) - fx
455 f_Col(2,atom2) = f_Col(2,atom2) - fy
456 f_Col(3,atom2) = f_Col(3,atom2) - fz
457
458 t_Row(1,atom1) = t_Row(1,atom1) + dUda*rxu1(1) - dUdg*uxu(1)
459 t_Row(2,atom1) = t_Row(2,atom1) + dUda*rxu1(2) - dUdg*uxu(2)
460 t_Row(3,atom1) = t_Row(3,atom1) + dUda*rxu1(3) - dUdg*uxu(3)
461
462 t_Col(1,atom2) = t_Col(1,atom2) + dUdb*rxu2(1) + dUdg*uxu(1)
463 t_Col(2,atom2) = t_Col(2,atom2) + dUdb*rxu2(2) + dUdg*uxu(2)
464 t_Col(3,atom2) = t_Col(3,atom2) + dUdb*rxu2(3) + dUdg*uxu(3)
465 #else
466 f(1,atom1) = f(1,atom1) + fx
467 f(2,atom1) = f(2,atom1) + fy
468 f(3,atom1) = f(3,atom1) + fz
469
470 f(1,atom2) = f(1,atom2) - fx
471 f(2,atom2) = f(2,atom2) - fy
472 f(3,atom2) = f(3,atom2) - fz
473
474 t(1,atom1) = t(1,atom1) + dUda*rxu1(1) - dUdg*uxu(1)
475 t(2,atom1) = t(2,atom1) + dUda*rxu1(2) - dUdg*uxu(2)
476 t(3,atom1) = t(3,atom1) + dUda*rxu1(3) - dUdg*uxu(3)
477
478 t(1,atom2) = t(1,atom2) + dUdb*rxu2(1) + dUdg*uxu(1)
479 t(2,atom2) = t(2,atom2) + dUdb*rxu2(2) + dUdg*uxu(2)
480 t(3,atom2) = t(3,atom2) + dUdb*rxu2(3) + dUdg*uxu(3)
481 #endif
482
483 if (do_pot) then
484 #ifdef IS_MPI
485 pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5d0*U*sw
486 pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5d0*U*sw
487 #else
488 pot = pot + U*sw
489 #endif
490 endif
491
492 vpair = vpair + U*sw
493 #ifdef IS_MPI
494 id1 = AtomRowToGlobal(atom1)
495 id2 = AtomColToGlobal(atom2)
496 #else
497 id1 = atom1
498 id2 = atom2
499 #endif
500
501 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
502
503 fpair(1) = fpair(1) + fx
504 fpair(2) = fpair(2) + fy
505 fpair(3) = fpair(3) + fz
506
507 endif
508
509 return
510 end subroutine do_gb_pair
511
512 subroutine destroyGBTypes()
513
514 GBMap%nGBtypes = 0
515 GBMap%currentGBtype = 0
516
517 if (associated(GBMap%GBtypes)) then
518 deallocate(GBMap%GBtypes)
519 GBMap%GBtypes => null()
520 end if
521
522 if (associated(GBMap%atidToGBtype)) then
523 deallocate(GBMap%atidToGBtype)
524 GBMap%atidToGBtype => null()
525 end if
526
527 haveMixingMap = .false.
528
529 end subroutine destroyGBTypes
530
531 end module gayberne
532