ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90
(Generate patch)

Comparing trunk/OOPSE-4/src/UseTheForce/DarkSide/gb.F90 (file contents):
Revision 1608 by gezelter, Wed Oct 20 04:02:48 2004 UTC vs.
Revision 2277 by chrisfen, Fri Aug 26 21:30:41 2005 UTC

# Line 1 | Line 1
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 gb_pair
44    use force_globals
45    use definitions
# Line 21 | Line 63 | module gb_pair
63    public :: check_gb_pair_FF
64    public :: set_gb_pair_params
65    public :: do_gb_pair
66 +  public :: getGayBerneCut
67  
68   contains
69  
# Line 46 | Line 89 | contains
89      return
90    end subroutine set_gb_pair_params
91  
92 +  function getGayBerneCut(atomID) result(cutValue)
93 +    integer, intent(in) :: atomID !! nah... we don't need to use this...
94 +    real(kind=dp) :: cutValue
95  
96 +    cutValue = gb_sigma*2.5_dp
97 +  end function getGayBerneCut
98 +
99    subroutine do_gb_pair(atom1, atom2, d, r, r2, sw, vpair, fpair, &
100 <       pot, u_l, f, t, do_pot)
100 >       pot, A, f, t, do_pot)
101      
102      integer, intent(in) :: atom1, atom2
103      integer :: id1, id2
# Line 56 | Line 105 | contains
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 <    real (kind=dp), dimension(3,nLocal) :: u_l
108 >    real (kind=dp), dimension(9,nLocal) :: A
109      real (kind=dp), dimension(3,nLocal) :: f
110      real (kind=dp), dimension(3,nLocal) :: t
111      logical, intent(in) :: do_pot
# Line 100 | Line 149 | contains
149      r4 = r2*r2
150  
151   #ifdef IS_MPI
152 <    ul1(1) = u_l_Row(1,atom1)
153 <    ul1(2) = u_l_Row(2,atom1)
154 <    ul1(3) = u_l_Row(3,atom1)
155 <
156 <    ul2(1) = u_l_Col(1,atom2)
157 <    ul2(2) = u_l_Col(2,atom2)
158 <    ul2(3) = u_l_Col(3,atom2)
152 >    ul1(1) = A_Row(3,atom1)
153 >    ul1(2) = A_Row(6,atom1)
154 >    ul1(3) = A_Row(9,atom1)
155 >
156 >    ul2(1) = A_Col(3,atom2)
157 >    ul2(2) = A_Col(6,atom2)
158 >    ul2(3) = A_Col(9,atom2)
159   #else
160 <    ul1(1) = u_l(1,atom1)
161 <    ul1(2) = u_l(2,atom1)
162 <    ul1(3) = u_l(3,atom1)
160 >    ul1(1) = A(3,atom1)
161 >    ul1(2) = A(6,atom1)
162 >    ul1(3) = A(9,atom1)
163  
164 <    ul2(1) = u_l(1,atom2)
165 <    ul2(2) = u_l(2,atom2)
166 <    ul2(3) = u_l(3,atom2)
164 >    ul2(1) = A(3,atom2)
165 >    ul2(2) = A(6,atom2)
166 >    ul2(3) = A(9,atom2)
167   #endif
168      
169      dru1dx = ul1(1)
# Line 178 | Line 227 | contains
227      dgdy = term1y + line3y
228      dgdz = term1z + line3z
229  
230 <    term1u1x = 2.0d0*(line1a+line2a)*d(1)
231 <    term1u1y = 2.0d0*(line1a+line2a)*d(2)
232 <    term1u1z = 2.0d0*(line1a+line2a)*d(3)
233 <    term1u2x = 2.0d0*(line1a-line2a)*d(1)
234 <    term1u2y = 2.0d0*(line1a-line2a)*d(2)
235 <    term1u2z = 2.0d0*(line1a-line2a)*d(3)
230 >    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      term2a = -line3a/opXdot
238      term2b =  line3b/omXdot
# Line 219 | Line 268 | contains
268      dBigRdx = drdx/gb_sigma + dgdx*gfact
269      dBigRdy = drdy/gb_sigma + dgdy*gfact
270      dBigRdz = drdz/gb_sigma + dgdz*gfact
271 +
272      dBigRdu1x = dgdu1x*gfact
273      dBigRdu1y = dgdu1y*gfact
274      dBigRdu1z = dgdu1z*gfact
275      dBigRdu2x = dgdu2x*gfact
276      dBigRdu2y = dgdu2y*gfact
277      dBigRdu2z = dgdu2z*gfact
278 <  
278 >
279      ! Now, we must do it again for g(ChiPrime) and dgpdx
280  
281      line1a = dotsum/opXpdot
# Line 244 | Line 294 | contains
294      dgpdy = term1y + line3y
295      dgpdz = term1z + line3z
296      
297 <    term1u1x = 2.0d0*(line1a+line2a)*d(1)
298 <    term1u1y = 2.0d0*(line1a+line2a)*d(2)
299 <    term1u1z = 2.0d0*(line1a+line2a)*d(3)
300 <    term1u2x = 2.0d0*(line1a-line2a)*d(1)
301 <    term1u2y = 2.0d0*(line1a-line2a)*d(2)
302 <    term1u2z = 2.0d0*(line1a-line2a)*d(3)
303 <    
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 >
304      term2a = -line3a/opXpdot
305      term2b =  line3b/omXpdot
306      
# Line 274 | Line 324 | contains
324      gmu = gp**gb_mu
325      gpi = 1.0d0 / gp
326      gmum = gmu*gpi
327 <  
278 <    ! write(*,*) atom1, atom2, Chi, u1dotu2
327 >
328      curlyE = 1.0d0/dsqrt(1.0d0 - Chi*Chi*u1dotu2*u1dotu2)
329  
330      dcE = (curlyE**3)*Chi*Chi*u1dotu2
# Line 383 | Line 432 | end module gb_pair
432    end subroutine do_gb_pair
433  
434   end module gb_pair
386
387
388  subroutine set_gb_pair_params(sigma, l2b_ratio, eps, eps_ratio, mu, nu)
389    use definitions, ONLY : dp
390    use gb_pair, ONLY : module_set_gb_pair_params => set_gb_pair_params
391    real( kind = dp ), intent(inout) :: sigma, l2b_ratio, eps, eps_ratio
392    real( kind = dp ), intent(inout) :: mu, nu
393    call module_set_gb_pair_params(sigma, l2b_ratio, eps, eps_ratio, mu, nu)
394 end subroutine set_gb_pair_params

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines