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 1948 by gezelter, Fri Jan 14 20:31:16 2005 UTC vs.
Revision 2062 by tim, Fri Feb 25 21:22:00 2005 UTC

# Line 43 | Line 43
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.7 2005-01-14 20:31:16 gezelter Exp $, $Date: 2005-01-14 20:31:16 $, $Name: not supported by cvs2svn $, $Revision: 1.7 $
46 > !! @version $Id: LJ.F90,v 1.8 2005-02-25 21:22:00 tim Exp $, $Date: 2005-02-25 21:22:00 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $
47  
48  
49   module lj
# Line 67 | Line 67 | module lj
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 80 | 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 100 | Line 102 | contains
102    
103   contains
104  
105 <  subroutine newLJtype(c_ident, sigma, epsilon, status)
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, myATID
112      integer, pointer :: MatchList(:) => null()
# Line 139 | Line 142 | contains
142      ParameterMap(myATID)%atid = myATID
143      ParameterMap(myATID)%epsilon = epsilon
144      ParameterMap(myATID)%sigma = sigma
145 <    
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)
# Line 236 | Line 243 | contains
243         MixingMap(i,i)%epsilon = Epsilon_i          
244         MixingMap(i,i)%delta   = -4.0_DP * MixingMap(i,i)%epsilon * &
245              (MixingMap(i,i)%tp12 - MixingMap(i,i)%tp6)
246 +       MixingMap(i,i)%soft_pot = ParameterMap(i)%soft_pot
247        
248         do j = i + 1, nATIDs
249            
# Line 260 | Line 268 | contains
268            
269            MixingMap(i,j)%delta = -4.0_DP * MixingMap(i,j)%epsilon * &
270                 (MixingMap(i,j)%tp12 - MixingMap(i,j)%tp6)
271 +
272 +          MixingMap(i,j)%soft_pot = ParameterMap(i)%soft_pot .or. ParameterMap(j)%soft_pot
273 +
274            
275            MixingMap(j,i)%sigma   = MixingMap(i,j)%sigma
276            MixingMap(j,i)%sigma6  = MixingMap(i,j)%sigma6
# Line 267 | Line 278 | contains
278            MixingMap(j,i)%tp12    = MixingMap(i,j)%tp12
279            MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
280            MixingMap(j,i)%delta   = MixingMap(i,j)%delta
281 +          MixingMap(j,i)%soft_pot   = MixingMap(i,j)%soft_pot
282            
283         end do
284      end do
# Line 296 | Line 308 | contains
308      real( kind = dp ) :: t6
309      real( kind = dp ) :: t12
310      real( kind = dp ) :: delta
311 +    logical :: soft_pot
312      integer :: id1, id2, localError
313  
314      if (.not.haveMixingMap) then
# Line 312 | Line 325 | contains
325      sigma6   = MixingMap(atid_Row(atom1),atid_Col(atom2))%sigma6
326      epsilon  = MixingMap(atid_Row(atom1),atid_Col(atom2))%epsilon
327      delta    = MixingMap(atid_Row(atom1),atid_Col(atom2))%delta
328 +    soft_pot =  MixingMap(atid_Row(atom1),atid_Col(atom2))%soft_pot
329   #else
330      sigma6   = MixingMap(atid(atom1),atid(atom2))%sigma6
331      epsilon  = MixingMap(atid(atom1),atid(atom2))%epsilon
332      delta    = MixingMap(atid(atom1),atid(atom2))%delta
333 +    soft_pot    = MixingMap(atid(atom1),atid(atom2))%soft_pot
334   #endif
335  
336      r6 = r2 * r2 * r2
337      
338      t6  = sigma6/ r6
339      t12 = t6 * t6    
325  
326    pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
327    if (LJ_do_shift) then
328       pot_temp = pot_temp + delta
329    endif
340  
341 <    vpair = vpair + pot_temp
341 >    if (soft_pot) then
342 >      pot_temp = 4.0E0_DP * epsilon * t12
343 >      if (LJ_do_shift) then
344 >         pot_temp = pot_temp + delta
345 >      endif
346 >  
347 >      vpair = vpair + pot_temp
348 >        
349 >      dudr = sw * 24.0E0_DP * epsilon * (-2.0E0_DP)*t12 / rij
350 >
351 >    else
352 >      pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
353 >      if (LJ_do_shift) then
354 >         pot_temp = pot_temp + delta
355 >      endif
356 >  
357 >      vpair = vpair + pot_temp
358 >        
359 >      dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
360 >    endif
361        
333    dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
334      
362      drdx = d(1) / rij
363      drdy = d(2) / rij
364      drdz = d(3) / rij

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines