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

Comparing trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90 (file contents):
Revision 1628 by gezelter, Thu Oct 21 20:15:31 2004 UTC vs.
Revision 2086 by gezelter, Tue Mar 8 21:06:12 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   !! Calculates Long Range forces Lennard-Jones interactions.
44   !! @author Charles F. Vardeman II
45   !! @author Matthew Meineke
46 < !! @version $Id: LJ.F90,v 1.3 2004-10-21 20:15:25 gezelter Exp $, $Date: 2004-10-21 20:15:25 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
46 > !! @version $Id: LJ.F90,v 1.9 2005-03-08 21:06:12 gezelter Exp $, $Date: 2005-03-08 21:06:12 $, $Name: not supported by cvs2svn $, $Revision: 1.9 $
47  
48 +
49   module lj
50    use atype_module
51    use switcheroo
# Line 20 | Line 63 | module lj
63    integer, parameter :: DP = selected_real_kind(15)
64    
65    type, private :: LjType
66 <     integer :: ident
66 >     integer :: c_ident
67 >     integer :: atid
68       real(kind=dp) :: sigma
69       real(kind=dp) :: epsilon
70 +     logical :: soft_pot
71    end type LjType
72    
73    type(LjType), dimension(:), allocatable :: ParameterMap
# Line 36 | Line 81 | module lj
81       real(kind=dp)  :: tp6
82       real(kind=dp)  :: tp12
83       real(kind=dp)  :: delta
84 +     logical :: soft_pot    
85    end type MixParameters
86    
87    type(MixParameters), dimension(:,:), allocatable :: MixingMap
# Line 51 | Line 97 | module lj
97    public :: useGeometricMixing
98    public :: do_lj_pair
99    public :: newLJtype  
100 +  public :: getSigma
101 +  public :: getEpsilon
102    
103   contains
104  
105 <  subroutine newLJtype(ident, sigma, epsilon, status)
106 <    integer,intent(in) :: ident
105 >  subroutine newLJtype(c_ident, sigma, epsilon, soft_pot, status)
106 >    integer,intent(in) :: c_ident
107      real(kind=dp),intent(in) :: sigma
108      real(kind=dp),intent(in) :: epsilon
109 +    integer, intent(in) :: soft_pot
110      integer,intent(out) :: status
111 <    integer :: nAtypes
111 >    integer :: nATypes, myATID
112 >    integer, pointer :: MatchList(:) => null()
113  
114      status = 0
115      
116 <    !! Be simple-minded and assume that we need a ParameterMap that
67 <    !! is the same size as the total number of atom types
116 >    myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
117  
118      if (.not.allocated(ParameterMap)) then
119        
120 +       !call getMatchingElementList(atypes, "is_LennardJones", .true., &
121 +       !     nLJTypes, MatchList)
122         nAtypes = getSize(atypes)
72    
123         if (nAtypes == 0) then
124            status = -1
125            return
126         end if
127 <      
127 >        
128         if (.not. allocated(ParameterMap)) then
129            allocate(ParameterMap(nAtypes))
130         endif
131 <      
131 >      
132      end if
133  
134 <    if (ident .gt. size(ParameterMap)) then
134 >    if (myATID .gt. size(ParameterMap)) then
135         status = -1
136         return
137      endif
138      
139      ! set the values for ParameterMap for this atom type:
140  
141 <    ParameterMap(ident)%ident = ident
142 <    ParameterMap(ident)%epsilon = epsilon
143 <    ParameterMap(ident)%sigma = sigma
144 <    
141 >    ParameterMap(myATID)%c_ident = c_ident
142 >    ParameterMap(myATID)%atid = myATID
143 >    ParameterMap(myATID)%epsilon = epsilon
144 >    ParameterMap(myATID)%sigma = sigma
145 >    if (soft_pot == 1) then
146 >      ParameterMap(myATID)%soft_pot = .true.
147 >    else
148 >      ParameterMap(myATID)%soft_pot = .false.
149 >    endif
150    end subroutine newLJtype
151 +
152 +  function getSigma(atid) result (s)
153 +    integer, intent(in) :: atid
154 +    integer :: localError
155 +    real(kind=dp) :: s
156      
157 +    if (.not.allocated(ParameterMap)) then
158 +       call handleError("LJ", "no ParameterMap was present before first call of getSigma!")
159 +       return
160 +    end if
161 +    
162 +    s = ParameterMap(atid)%sigma
163 +  end function getSigma
164 +
165 +  function getEpsilon(atid) result (e)
166 +    integer, intent(in) :: atid
167 +    integer :: localError
168 +    real(kind=dp) :: e
169 +    
170 +    if (.not.allocated(ParameterMap)) then
171 +       call handleError("LJ", "no ParameterMap was present before first call of getEpsilon!")
172 +       return
173 +    end if
174 +    
175 +    e = ParameterMap(atid)%epsilon
176 +  end function getEpsilon
177 +
178 +
179    subroutine setCutoffLJ(rcut, do_shift, status)
180      logical, intent(in):: do_shift
181      integer :: status, myStatus
# Line 119 | Line 201 | contains
201    end subroutine useGeometricMixing
202    
203    subroutine createMixingMap(status)
204 <    integer :: nAtypes
204 >    integer :: nATIDs
205      integer :: status
206      integer :: i
207      integer :: j
208      real ( kind = dp ) :: Sigma_i, Sigma_j
209      real ( kind = dp ) :: Epsilon_i, Epsilon_j
210      real ( kind = dp ) :: rcut6
211 +    logical :: i_is_LJ, j_is_LJ
212  
213      status = 0
214 +
215 +    if (.not. allocated(ParameterMap)) then
216 +       call handleError("LJ", "no ParameterMap was present before call of createMixingMap!")
217 +       status = -1
218 +       return
219 +    endif
220      
221 <    nAtypes = size(ParameterMap)
221 >    nATIDs = size(ParameterMap)
222      
223 <    if (nAtypes == 0) then
223 >    if (nATIDs == 0) then
224         status = -1
225         return
226      end if
227  
228 +    if (.not. allocated(MixingMap)) then
229 +       allocate(MixingMap(nATIDs, nATIDs))
230 +    endif
231 +
232      if (.not.have_rcut) then
233         status = -1
234         return
235      endif
236 <    
144 <    if (.not. allocated(MixingMap)) then
145 <       allocate(MixingMap(nAtypes, nAtypes))
146 <    endif
147 <    
236 >        
237      rcut6 = LJ_rcut**6
238      
239      ! This loops through all atypes, even those that don't support LJ forces.
240 <    do i = 1, nAtypes
241 <      
242 <       Epsilon_i = ParameterMap(i)%epsilon
243 <       Sigma_i = ParameterMap(i)%sigma
244 <      
156 <       ! do self mixing rule
157 <       MixingMap(i,i)%sigma   = Sigma_i          
158 <       MixingMap(i,i)%sigma6  = Sigma_i ** 6          
159 <       MixingMap(i,i)%tp6     = (MixingMap(i,i)%sigma6)/rcut6          
160 <       MixingMap(i,i)%tp12    = (MixingMap(i,i)%tp6) ** 2
161 <       MixingMap(i,i)%epsilon = Epsilon_i          
162 <       MixingMap(i,i)%delta   = -4.0_DP * MixingMap(i,i)%epsilon * &
163 <            (MixingMap(i,i)%tp12 - MixingMap(i,i)%tp6)
164 <      
165 <       do j = i + 1, nAtypes
240 >    do i = 1, nATIDs
241 >       call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
242 >       if (i_is_LJ) then
243 >          Epsilon_i = ParameterMap(i)%epsilon
244 >          Sigma_i = ParameterMap(i)%sigma
245            
246 <          Epsilon_j = ParameterMap(j)%epsilon
247 <          Sigma_j = ParameterMap(j)%sigma
246 >          ! do self mixing rule
247 >          MixingMap(i,i)%sigma   = Sigma_i          
248 >          MixingMap(i,i)%sigma6  = Sigma_i ** 6          
249 >          MixingMap(i,i)%tp6     = (MixingMap(i,i)%sigma6)/rcut6          
250 >          MixingMap(i,i)%tp12    = (MixingMap(i,i)%tp6) ** 2
251 >          MixingMap(i,i)%epsilon = Epsilon_i          
252 >          MixingMap(i,i)%delta   = -4.0_DP * MixingMap(i,i)%epsilon * &
253 >               (MixingMap(i,i)%tp12 - MixingMap(i,i)%tp6)
254 >          MixingMap(i,i)%soft_pot = ParameterMap(i)%soft_pot
255            
256 <          ! only the distance parameter uses different mixing policies
257 <          if (useGeometricDistanceMixing) then
172 <             ! only for OPLS as far as we can tell
173 <             MixingMap(i,j)%sigma = dsqrt(Sigma_i * Sigma_j)
174 <          else
175 <             ! everyone else
176 <             MixingMap(i,j)%sigma = 0.5_dp * (Sigma_i + Sigma_j)
177 <          endif
256 >          do j = i + 1, nATIDs
257 >             call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ)
258            
259 <          ! energy parameter is always geometric mean:
260 <          MixingMap(i,j)%epsilon = dsqrt(Epsilon_i * Epsilon_j)
261 <                    
262 <          MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6
263 <          MixingMap(i,j)%tp6    = MixingMap(i,j)%sigma6/rcut6
264 <          MixingMap(i,j)%tp12    = (MixingMap(i,j)%tp6) ** 2
265 <          
266 <          MixingMap(i,j)%delta = -4.0_DP * MixingMap(i,j)%epsilon * &
267 <               (MixingMap(i,j)%tp12 - MixingMap(i,j)%tp6)
268 <          
269 <          MixingMap(j,i)%sigma   = MixingMap(i,j)%sigma
270 <          MixingMap(j,i)%sigma6  = MixingMap(i,j)%sigma6
271 <          MixingMap(j,i)%tp6     = MixingMap(i,j)%tp6
272 <          MixingMap(j,i)%tp12    = MixingMap(i,j)%tp12
273 <          MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
274 <          MixingMap(j,i)%delta   = MixingMap(i,j)%delta
275 <          
276 <       end do
259 >             if (j_is_LJ) then
260 >                Epsilon_j = ParameterMap(j)%epsilon
261 >                Sigma_j = ParameterMap(j)%sigma
262 >                
263 >                ! only the distance parameter uses different mixing policies
264 >                if (useGeometricDistanceMixing) then
265 >                   ! only for OPLS as far as we can tell
266 >                   MixingMap(i,j)%sigma = dsqrt(Sigma_i * Sigma_j)
267 >                else
268 >                   ! everyone else
269 >                   MixingMap(i,j)%sigma = 0.5_dp * (Sigma_i + Sigma_j)
270 >                endif
271 >                
272 >                ! energy parameter is always geometric mean:
273 >                MixingMap(i,j)%epsilon = dsqrt(Epsilon_i * Epsilon_j)
274 >                
275 >                MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6
276 >                MixingMap(i,j)%tp6    = MixingMap(i,j)%sigma6/rcut6
277 >                MixingMap(i,j)%tp12    = (MixingMap(i,j)%tp6) ** 2
278 >                
279 >                MixingMap(i,j)%delta = -4.0_DP * MixingMap(i,j)%epsilon * &
280 >                     (MixingMap(i,j)%tp12 - MixingMap(i,j)%tp6)
281 >                
282 >                MixingMap(i,j)%soft_pot = ParameterMap(i)%soft_pot .or. ParameterMap(j)%soft_pot
283 >                
284 >                
285 >                MixingMap(j,i)%sigma   = MixingMap(i,j)%sigma
286 >                MixingMap(j,i)%sigma6  = MixingMap(i,j)%sigma6
287 >                MixingMap(j,i)%tp6     = MixingMap(i,j)%tp6
288 >                MixingMap(j,i)%tp12    = MixingMap(i,j)%tp12
289 >                MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
290 >                MixingMap(j,i)%delta   = MixingMap(i,j)%delta
291 >                MixingMap(j,i)%soft_pot   = MixingMap(i,j)%soft_pot
292 >             endif
293 >          end do
294 >       endif
295      end do
296      
297 +    haveMixingMap = .true.
298 +
299    end subroutine createMixingMap
300          
301    subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
# Line 219 | Line 319 | contains
319      real( kind = dp ) :: t6
320      real( kind = dp ) :: t12
321      real( kind = dp ) :: delta
322 +    logical :: soft_pot
323      integer :: id1, id2, localError
324  
325      if (.not.haveMixingMap) then
# Line 235 | Line 336 | contains
336      sigma6   = MixingMap(atid_Row(atom1),atid_Col(atom2))%sigma6
337      epsilon  = MixingMap(atid_Row(atom1),atid_Col(atom2))%epsilon
338      delta    = MixingMap(atid_Row(atom1),atid_Col(atom2))%delta
339 +    soft_pot =  MixingMap(atid_Row(atom1),atid_Col(atom2))%soft_pot
340   #else
341      sigma6   = MixingMap(atid(atom1),atid(atom2))%sigma6
342      epsilon  = MixingMap(atid(atom1),atid(atom2))%epsilon
343      delta    = MixingMap(atid(atom1),atid(atom2))%delta
344 +    soft_pot    = MixingMap(atid(atom1),atid(atom2))%soft_pot
345   #endif
346  
347      r6 = r2 * r2 * r2
348      
349      t6  = sigma6/ r6
350      t12 = t6 * t6    
248  
249    pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
250    if (LJ_do_shift) then
251       pot_temp = pot_temp + delta
252    endif
351  
352 <    vpair = vpair + pot_temp
352 >    if (soft_pot) then
353 >      pot_temp = 4.0E0_DP * epsilon * t12
354 >      if (LJ_do_shift) then
355 >         pot_temp = pot_temp + delta
356 >      endif
357 >  
358 >      vpair = vpair + pot_temp
359 >        
360 >      dudr = sw * 24.0E0_DP * epsilon * (-2.0E0_DP)*t12 / rij
361 >
362 >    else
363 >      pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
364 >      if (LJ_do_shift) then
365 >         pot_temp = pot_temp + delta
366 >      endif
367 >  
368 >      vpair = vpair + pot_temp
369 >        
370 >      dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
371 >    endif
372        
256    dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
257      
373      drdx = d(1) / rij
374      drdy = d(2) / rij
375      drdz = d(3) / rij
# Line 314 | Line 429 | end module lj
429    !! Calculates the mixing for sigma or epslon
430      
431   end module lj
317
318 subroutine newLJtype(ident, sigma, epsilon, status)
319  use lj, ONLY : module_newLJtype => newLJtype
320  integer, parameter :: DP = selected_real_kind(15)
321  integer,intent(inout) :: ident
322  real(kind=dp),intent(inout) :: sigma
323  real(kind=dp),intent(inout) :: epsilon
324  integer,intent(inout) :: status
325  
326  call module_newLJtype(ident, sigma, epsilon, status)
327  
328 end subroutine newLJtype
329
330 subroutine useGeometricMixing()
331  use lj, ONLY: module_useGeometricMixing => useGeometricMixing
332  
333  call module_useGeometricMixing()
334  return
335 end subroutine useGeometricMixing

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines