ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/gb.F90
Revision: 1006
Committed: Mon Jul 17 18:13:26 2006 UTC (19 years ago) by xsun
File size: 16631 byte(s)
Log Message:
fixed d / sigma0 confusion

File Contents

# User Rev Content
1 gezelter 246 !!
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 676 module gayberne
44 gezelter 115 use force_globals
45     use definitions
46     use simulation
47 kdaily 668 use atype_module
48     use vector_class
49 gezelter 981 use linearalgebra
50 kdaily 668 use status
51     use lj
52 gezelter 982 use fForceOptions
53 gezelter 115 #ifdef IS_MPI
54     use mpiSimulation
55     #endif
56 kdaily 529
57 gezelter 115 implicit none
58    
59 kdaily 668 private
60 gezelter 115
61 gezelter 676 #define __FORTRAN90
62     #include "UseTheForce/DarkSide/fInteractionMap.h"
63 gezelter 115
64 gezelter 981 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 gezelter 676 public :: newGBtype
71 gezelter 981 public :: complete_GB_FF
72 gezelter 115 public :: do_gb_pair
73 chrisfen 578 public :: getGayBerneCut
74 gezelter 676 public :: destroyGBtypes
75 gezelter 115
76 gezelter 676 type :: GBtype
77     integer :: atid
78 gezelter 981 real(kind = dp ) :: d
79     real(kind = dp ) :: l
80 gezelter 676 real(kind = dp ) :: eps
81     real(kind = dp ) :: eps_ratio
82 gezelter 981 real(kind = dp ) :: dw
83     logical :: isLJ
84 gezelter 676 end type GBtype
85 gezelter 981
86 gezelter 676 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 gezelter 981
93 gezelter 676 type(GBList), save :: GBMap
94 gezelter 981
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 gezelter 115 contains
110 gezelter 981
111     subroutine newGBtype(c_ident, d, l, eps, eps_ratio, dw, status)
112 gezelter 676
113     integer, intent(in) :: c_ident
114 gezelter 981 real( kind = dp ), intent(in) :: d, l, eps, eps_ratio, dw
115 gezelter 676 integer, intent(out) :: status
116 gezelter 981
117     integer :: nGBTypes, nLJTypes, ntypes, myATID
118 gezelter 676 integer, pointer :: MatchList(:) => null()
119     integer :: current, i
120 kdaily 668 status = 0
121 gezelter 981
122 gezelter 676 if (.not.associated(GBMap%GBtypes)) then
123 gezelter 981
124 gezelter 676 call getMatchingElementList(atypes, "is_GayBerne", .true., &
125     nGBtypes, MatchList)
126    
127 gezelter 981 call getMatchingElementList(atypes, "is_LennardJones", .true., &
128     nLJTypes, MatchList)
129    
130     GBMap%nGBtypes = nGBtypes + nLJTypes
131    
132     allocate(GBMap%GBtypes(nGBtypes + nLJTypes))
133    
134 gezelter 676 ntypes = getSize(atypes)
135    
136 gezelter 981 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 gezelter 676
180 gezelter 981 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 xsun 1006 GBMap%GBtypes(current)%d = getSigma(myATID) / sqrt(2.0_dp)
191 gezelter 981 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 gezelter 982
202 gezelter 981
203     end subroutine complete_GB_FF
204 kdaily 668
205 gezelter 981 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 kdaily 668
211 gezelter 981 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 gezelter 676 endif
221 kdaily 668
222 gezelter 981 do i = 1, nGBtypes
223 kdaily 668
224 gezelter 981 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 gezelter 676
230 gezelter 981 do j = i, nGBtypes
231 gezelter 676
232 gezelter 981 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 xsun 1006
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 gezelter 981 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    
261     if (i.ne.j) then
262     GBMixingMap(j,i)%sigma0 = GBMixingMap(i,j)%sigma0
263     GBMixingMap(j,i)%dw = GBMixingMap(i,j)%dw
264     GBMixingMap(j,i)%eps0 = GBMixingMap(i,j)%eps0
265     GBMixingMap(j,i)%x2 = GBMixingMap(i,j)%x2
266     GBMixingMap(j,i)%xa2 = GBMixingMap(i,j)%xa2
267     GBMixingMap(j,i)%xai2 = GBMixingMap(i,j)%xai2
268     GBMixingMap(j,i)%xp2 = GBMixingMap(i,j)%xp2
269     GBMixingMap(j,i)%xpap2 = GBMixingMap(i,j)%xpap2
270     GBMixingMap(j,i)%xpapi2 = GBMixingMap(i,j)%xpapi2
271     endif
272     enddo
273     enddo
274     haveMixingMap = .true.
275 gezelter 983 mu = getGayBerneMu()
276     nu = getGayBerneNu()
277 gezelter 981 end subroutine createGBMixingMap
278 gezelter 676
279 gezelter 981
280 chrisfen 580 !! gay berne cutoff should be a parameter in globals, this is a temporary
281     !! work around - this should be fixed when gay berne is up and running
282 gezelter 676
283 chrisfen 578 function getGayBerneCut(atomID) result(cutValue)
284 gezelter 676 integer, intent(in) :: atomID
285     integer :: gbt1
286 gezelter 981 real(kind=dp) :: cutValue, l, d
287 gezelter 115
288 gezelter 676 if (GBMap%currentGBtype == 0) then
289     call handleError("GB", "No members in GBMap")
290     return
291     end if
292    
293     gbt1 = GBMap%atidToGBtype(atomID)
294 gezelter 981 l = GBMap%GBtypes(gbt1)%l
295     d = GBMap%GBtypes(gbt1)%d
296     cutValue = 2.5_dp*max(l,d)
297 gezelter 676
298 chrisfen 578 end function getGayBerneCut
299    
300 gezelter 115 subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
301 gezelter 981 pot, Amat, f, t, do_pot)
302 kdaily 529
303 gezelter 115 integer, intent(in) :: atom1, atom2
304 gezelter 676 integer :: atid1, atid2, gbt1, gbt2, id1, id2
305 gezelter 115 real (kind=dp), intent(inout) :: r, r2
306     real (kind=dp), dimension(3), intent(in) :: d
307     real (kind=dp), dimension(3), intent(inout) :: fpair
308     real (kind=dp) :: pot, sw, vpair
309 gezelter 981 real (kind=dp), dimension(9,nLocal) :: Amat
310 gezelter 115 real (kind=dp), dimension(3,nLocal) :: f
311     real (kind=dp), dimension(3,nLocal) :: t
312     logical, intent(in) :: do_pot
313 gezelter 981 real (kind = dp), dimension(3) :: ul1, ul2, rxu1, rxu2, uxu, rhat
314 gezelter 115
315 gezelter 981 real (kind = dp) :: sigma0, dw, eps0, x2, xa2, xai2, xp2, xpap2, xpapi2
316 kdaily 997 real (kind = dp) :: e1, e2, eps, sigma, s3, s03, au2, bu2, au, bu, a, b, g, g2
317 gezelter 981 real (kind = dp) :: U, BigR, R3, R6, R7, R12, R13, H, Hp, fx, fy, fz
318     real (kind = dp) :: dUdr, dUda, dUdb, dUdg, pref1, pref2
319     logical :: i_is_lj, j_is_lj
320    
321     if (.not.haveMixingMap) then
322     call createGBMixingMap()
323     endif
324    
325 gezelter 676 #ifdef IS_MPI
326     atid1 = atid_Row(atom1)
327     atid2 = atid_Col(atom2)
328     #else
329     atid1 = atid(atom1)
330     atid2 = atid(atom2)
331     #endif
332 gezelter 115
333 gezelter 676 gbt1 = GBMap%atidToGBtype(atid1)
334 gezelter 981 gbt2 = GBMap%atidToGBtype(atid2)
335 gezelter 676
336 gezelter 981 i_is_LJ = GBMap%GBTypes(gbt1)%isLJ
337     j_is_LJ = GBMap%GBTypes(gbt2)%isLJ
338 gezelter 676
339 gezelter 981 sigma0 = GBMixingMap(gbt1, gbt2)%sigma0
340     dw = GBMixingMap(gbt1, gbt2)%dw
341     eps0 = GBMixingMap(gbt1, gbt2)%eps0
342     x2 = GBMixingMap(gbt1, gbt2)%x2
343     xa2 = GBMixingMap(gbt1, gbt2)%xa2
344     xai2 = GBMixingMap(gbt1, gbt2)%xai2
345     xp2 = GBMixingMap(gbt1, gbt2)%xp2
346     xpap2 = GBMixingMap(gbt1, gbt2)%xpap2
347     xpapi2 = GBMixingMap(gbt1, gbt2)%xpapi2
348    
349 gezelter 115 #ifdef IS_MPI
350 gezelter 686 ul1(1) = A_Row(7,atom1)
351     ul1(2) = A_Row(8,atom1)
352 gezelter 246 ul1(3) = A_Row(9,atom1)
353 gezelter 115
354 gezelter 686 ul2(1) = A_Col(7,atom2)
355     ul2(2) = A_Col(8,atom2)
356 gezelter 246 ul2(3) = A_Col(9,atom2)
357 gezelter 115 #else
358 gezelter 981 ul1(1) = Amat(7,atom1)
359     ul1(2) = Amat(8,atom1)
360     ul1(3) = Amat(9,atom1)
361 gezelter 115
362 gezelter 981 ul2(1) = Amat(7,atom2)
363     ul2(2) = Amat(8,atom2)
364     ul2(3) = Amat(9,atom2)
365 gezelter 115 #endif
366 kdaily 529
367 gezelter 981 if (i_is_LJ) then
368     a = 0.0_dp
369     ul1 = 0.0_dp
370     else
371     a = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
372     endif
373 gezelter 115
374 gezelter 981 if (j_is_LJ) then
375     b = 0.0_dp
376     ul2 = 0.0_dp
377     else
378     b = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
379     endif
380 gezelter 115
381 gezelter 981 if (i_is_LJ.or.j_is_LJ) then
382     g = 0.0_dp
383     else
384     g = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
385     endif
386 gezelter 1002
387 gezelter 981 au = a / r
388     bu = b / r
389 kdaily 997
390 gezelter 1002 au2 = au * au
391     bu2 = bu * bu
392     g2 = g * g
393 gezelter 115
394 kdaily 997 H = (xa2 * au2 + xai2 * bu2 - 2.0_dp*x2*au*bu*g) / (1.0_dp - x2*g2)
395     Hp = (xpap2*au2 + xpapi2*bu2 - 2.0_dp*xp2*au*bu*g) / (1.0_dp - xp2*g2)
396 gezelter 1002
397 gezelter 981 sigma = sigma0 / sqrt(1.0_dp - H)
398     e1 = 1.0_dp / sqrt(1.0_dp - x2*g2)
399     e2 = 1.0_dp - Hp
400     eps = eps0 * (e1**nu) * (e2**mu)
401     BigR = dw*sigma0 / (r - sigma + dw*sigma0)
402    
403     R3 = BigR*BigR*BigR
404     R6 = R3*R3
405     R7 = R6 * BigR
406     R12 = R6*R6
407     R13 = R6*R7
408 gezelter 115
409 gezelter 981 U = 4.0_dp * eps * (R12 - R6)
410 gezelter 115
411 gezelter 981 s3 = sigma*sigma*sigma
412     s03 = sigma0*sigma0*sigma0
413 kdaily 529
414 gezelter 981 pref1 = - 8.0_dp * eps * mu * (R12 - R6) / (e2 * r)
415 gezelter 983
416 gezelter 981 pref2 = 8.0_dp * eps * s3 * (6.0_dp*R13 - 3.0_dp*R7) / (dw*r*s03)
417 gezelter 507
418 kdaily 997 dUdr = - (pref1 * Hp + pref2 * (sigma0*sigma0*r/s3 + H))
419 kdaily 529
420 gezelter 981 dUda = pref1 * (xpap2*au - xp2*bu*g) / (1.0_dp - xp2 * g2) &
421     + pref2 * (xa2 * au - x2 *bu*g) / (1.0_dp - x2 * g2)
422 kdaily 997
423 gezelter 981 dUdb = pref1 * (xpapi2*bu - xp2*au*g) / (1.0_dp - xp2 * g2) &
424     + pref2 * (xai2 * bu - x2 *au*g) / (1.0_dp - x2 * g2)
425 gezelter 507
426 gezelter 981 dUdg = 4.0_dp * eps * nu * (R12 - R6) * x2 * g / (1.0_dp - x2*g2) &
427     + 8.0_dp * eps * mu * (R12 - R6) * (xp2*au*bu - Hp*xp2*g) / &
428     (1.0_dp - xp2 * g2) / e2 &
429 gezelter 983 + 8.0_dp * eps * s3 * (3.0_dp * R7 - 6.0_dp * R13) * &
430 gezelter 981 (x2 * au * bu - H * x2 * g) / (1.0_dp - x2 * g2) / (dw * s03)
431 kdaily 997
432 gezelter 981
433     rhat = d / r
434 gezelter 507
435 gezelter 983 fx = dUdr * rhat(1) + dUda * ul1(1) + dUdb * ul2(1)
436     fy = dUdr * rhat(2) + dUda * ul1(2) + dUdb * ul2(2)
437     fz = dUdr * rhat(3) + dUda * ul1(3) + dUdb * ul2(3)
438 gezelter 115
439 gezelter 981 rxu1 = cross_product(d, ul1)
440     rxu2 = cross_product(d, ul2)
441     uxu = cross_product(ul1, ul2)
442 kdaily 997
443     !!$ write(*,*) 'pref = ' , pref1, pref2
444 gezelter 983 !!$ write(*,*) 'rxu1 = ' , rxu1(1), rxu1(2), rxu1(3)
445     !!$ write(*,*) 'rxu2 = ' , rxu2(1), rxu2(2), rxu2(3)
446     !!$ write(*,*) 'uxu = ' , uxu(1), uxu(2), uxu(3)
447 kdaily 997 !!$ write(*,*) 'dUda = ', dUda, dudb, dudg, dudr
448     !!$ write(*,*) 'H = ', H,hp,sigma, e1, e2, BigR
449     !!$ write(*,*) 'chi = ', xa2, xai2, x2
450     !!$ write(*,*) 'chip = ', xpap2, xpapi2, xp2
451     !!$ write(*,*) 'eps = ', eps0, e1, e2, eps
452     !!$ write(*,*) 'U =', U, pref1, pref2
453     !!$ write(*,*) 'f =', fx, fy, fz
454     !!$ write(*,*) 'au =', au, bu, g
455     !!$
456    
457    
458 gezelter 115 #ifdef IS_MPI
459 gezelter 981 f_Row(1,atom1) = f_Row(1,atom1) + fx
460     f_Row(2,atom1) = f_Row(2,atom1) + fy
461     f_Row(3,atom1) = f_Row(3,atom1) + fz
462 kdaily 529
463 gezelter 981 f_Col(1,atom2) = f_Col(1,atom2) - fx
464     f_Col(2,atom2) = f_Col(2,atom2) - fy
465     f_Col(3,atom2) = f_Col(3,atom2) - fz
466 kdaily 529
467 gezelter 981 t_Row(1,atom1) = t_Row(1,atom1) + dUda*rxu1(1) - dUdg*uxu(1)
468     t_Row(2,atom1) = t_Row(2,atom1) + dUda*rxu1(2) - dUdg*uxu(2)
469     t_Row(3,atom1) = t_Row(3,atom1) + dUda*rxu1(3) - dUdg*uxu(3)
470    
471     t_Col(1,atom2) = t_Col(1,atom2) + dUdb*rxu2(1) + dUdg*uxu(1)
472     t_Col(2,atom2) = t_Col(2,atom2) + dUdb*rxu2(2) + dUdg*uxu(2)
473     t_Col(3,atom2) = t_Col(3,atom2) + dUdb*rxu2(3) + dUdg*uxu(3)
474 gezelter 115 #else
475 gezelter 981 f(1,atom1) = f(1,atom1) + fx
476     f(2,atom1) = f(2,atom1) + fy
477     f(3,atom1) = f(3,atom1) + fz
478 kdaily 529
479 gezelter 981 f(1,atom2) = f(1,atom2) - fx
480     f(2,atom2) = f(2,atom2) - fy
481     f(3,atom2) = f(3,atom2) - fz
482 kdaily 529
483 gezelter 981 t(1,atom1) = t(1,atom1) + dUda*rxu1(1) - dUdg*uxu(1)
484     t(2,atom1) = t(2,atom1) + dUda*rxu1(2) - dUdg*uxu(2)
485     t(3,atom1) = t(3,atom1) + dUda*rxu1(3) - dUdg*uxu(3)
486    
487     t(1,atom2) = t(1,atom2) + dUdb*rxu2(1) + dUdg*uxu(1)
488     t(2,atom2) = t(2,atom2) + dUdb*rxu2(2) + dUdg*uxu(2)
489     t(3,atom2) = t(3,atom2) + dUdb*rxu2(3) + dUdg*uxu(3)
490 gezelter 115 #endif
491 kdaily 668
492 gezelter 115 if (do_pot) then
493     #ifdef IS_MPI
494 gezelter 981 pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5d0*U*sw
495     pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5d0*U*sw
496 gezelter 115 #else
497 gezelter 981 pot = pot + U*sw
498 gezelter 115 #endif
499     endif
500 gezelter 676
501 gezelter 981 vpair = vpair + U*sw
502 gezelter 115 #ifdef IS_MPI
503     id1 = AtomRowToGlobal(atom1)
504     id2 = AtomColToGlobal(atom2)
505     #else
506     id1 = atom1
507     id2 = atom2
508     #endif
509 kdaily 529
510 gezelter 115 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
511 kdaily 529
512 gezelter 981 fpair(1) = fpair(1) + fx
513     fpair(2) = fpair(2) + fy
514     fpair(3) = fpair(3) + fz
515 kdaily 529
516 gezelter 115 endif
517 kdaily 529
518 gezelter 115 return
519     end subroutine do_gb_pair
520 gezelter 981
521 gezelter 676 subroutine destroyGBTypes()
522    
523     GBMap%nGBtypes = 0
524     GBMap%currentGBtype = 0
525 kdaily 668
526 gezelter 676 if (associated(GBMap%GBtypes)) then
527     deallocate(GBMap%GBtypes)
528     GBMap%GBtypes => null()
529 kdaily 668 end if
530    
531 gezelter 676 if (associated(GBMap%atidToGBtype)) then
532     deallocate(GBMap%atidToGBtype)
533     GBMap%atidToGBtype => null()
534     end if
535 kdaily 668
536 gezelter 981 haveMixingMap = .false.
537    
538 gezelter 676 end subroutine destroyGBTypes
539 kdaily 668
540 gezelter 676 end module gayberne
541 kdaily 668