ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_LJ_FF.F90
(Generate patch)

Comparing trunk/OOPSE_old/src/mdtools/libmdCode/calc_LJ_FF.F90 (file contents):
Revision 309 by gezelter, Mon Mar 10 23:19:23 2003 UTC vs.
Revision 320 by chuckv, Wed Mar 12 14:55:29 2003 UTC

# Line 2 | Line 2
2   !! Corresponds to the force field defined in lj_FF.cpp
3   !! @author Charles F. Vardeman II
4   !! @author Matthew Meineke
5 < !! @version $Id: calc_LJ_FF.F90,v 1.3 2003-03-10 23:19:23 gezelter Exp $, $Date: 2003-03-10 23:19:23 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
5 > !! @version $Id: calc_LJ_FF.F90,v 1.5 2003-03-12 14:55:29 chuckv Exp $, $Date: 2003-03-12 14:55:29 $, $Name: not supported by cvs2svn $, $Revision: 1.5 $
6  
7  
8  
# Line 11 | Line 11 | module lj
11    use definitions
12    use forceGlobals
13    use atype_typedefs
14 <  use generic_atypes
14 >  use vector_class
15   #ifdef IS_MPI
16    use mpiSimulation
17   #endif
# Line 31 | Line 31 | module lj
31       real ( kind = dp )  :: epsilon = 0.0_dp
32       !! Lennard-Jones Sigma
33       real ( kind = dp )  :: sigma = 0.0_dp
34     !! Lennard-Jones Sigma Squared
35     real ( kind = dp )  :: sigma2 = 0.0_dp
34       !! Lennard-Jones Sigma to sixth
35       real ( kind = dp )  :: sigma6 = 0.0_dp
36 +     !!
37 +     real ( kind = dp )  :: tp6
38 +     real ( kind = dp )  :: tp12
39 +
40 +     real ( kind = dp )  :: delta  = 0.0_dp
41 +
42    end type lj_mixed_params
43    
44    type (lj_mixed_params), dimension(:,:), pointer :: ljMixed
45    
46   contains
47  
48 <  subroutine init_ljFF(ListHead,status)
45 <    type (atype),pointer :: ListHead
46 <    integer :: nTypes
48 >  subroutine init_ljFF(status)
49      integer, intent(out) :: status
48    
50      integer :: myStatus
51      
52      status = 0
53 <    call createMixingList(ListHead,myStatus)
53 >    call createMixingList(myStatus)
54      if (myStatus /= 0) then
55         status = -1
56         return
# Line 57 | Line 58 | contains
58      
59    end subroutine init_ljFF
60    
61 <  subroutine createMixingList(ListHead,status)
62 <    type (atype), pointer :: ListHead
62 <    integer :: listSize
61 >  subroutine createMixingList(status)
62 >    integer :: nAtypes
63      integer :: status
64      integer :: i
65      integer :: j
66 <    type (atype), pointer :: tmpPtr_i => null()
67 <    type (atype), pointer :: tmpPtr_j => null()
66 >    real ( kind = dp ) :: mySigma_i,mySigma_j
67 >    real ( kind = dp ) :: myEpsilon_i,myEpsilon_j
68 >    real ( kind = dp ) :: rcut, rcut6
69      status = 0
70      
71 <    listSize = get_this_ListLen(ListHead)
71 >    nAtypes = getSize(atype)
72      if (listSize == 0) then
73         status = -1
74         return
# Line 80 | Line 81 | contains
81         status = -1
82         return
83      end if
84 <    
85 <    
86 <    
87 <    tmpPtr_i => ListHead
88 <    tmpPtr_j => tmpPtr_i%next
88 <    do while (associated(tmpPtr_i))
89 <       i = i + 1
84 >    call getRcut(rcut, rcut6=rcut6)
85 >    do i = 1, nAtypes
86 >
87 >       call getElementProperty(atypes,i,"lj_epsilon",myEpsilon_i)
88 >       call getElementProperty(atypes,i,"lj_sigma",  mySigma_i)
89         ! do self mixing rule
90 <       ljMixed(i,i)%sigma   = tmpPtr_i%lj_sigma
90 >       ljMixed(i,i)%sigma   = mySigma_i
91        
93       ljMixed(i,i)%sigma2  = (ljMixed(i,i)%sigma) ** 2
94      
92         ljMixed(i,i)%sigma6  = (ljMixed(i,i)%sigma) ** 6                                                                                                  
93 <       ljMixed(i,i)%epsilon = tmpPtr_i%lj_epsilon
94 <      
95 <       j = i + 1
96 <       do while (associated(tmpPtr_j))
93 >       ljMixed(i,i)%tp6     = ljMixed(i,i)%sigma6/rcut6
94 >
95 >       ljMixed(i,i)%tp12    = (ljMixed(i,i)%tp6) ** 2
96 >
97 >       ljMixed(i,i)%epsilon = myEpsilon_i
98 >
99 >       ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
100 >            (ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6)
101 >      
102 >       do j = i + 1, nAtypes
103 >          call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j)
104 >          call getElementProperty(atypes,j,"lj_sigma",  mySigma_j)
105            
106            ljMixed(i,j)%sigma  =  &
107 <               calcLJMix("sigma",tmpPtr_i%lj_sigma, &
108 <               tmpPtr_j%lj_sigma)
107 >               calcLJMix("sigma",mySigma_i, &
108 >               mySigma_j)
109            
105          ljMixed(i,j)%sigma2  = &
106               (ljMixed(i,j)%sigma)**2
107          
110            ljMixed(i,j)%sigma6 = &
111                 (ljMixed(i,j)%sigma)**6
112            
113 +          
114 +          ljMixed(i,j)%tp6     = ljMixed(i,j)%sigma6/rcut6
115 +          
116 +          ljMixed(i,j)%tp12    = (ljMixed(i,j)%tp6) ** 2
117 +          
118 +          
119            ljMixed(i,j)%epsilon = &
120 <               calcLJMix("epsilon",tmpPtr_i%lj_epsilon, &
121 <               tmpPtr_j%lj_epsilon)
122 <          ljMixed(j,i)%sigma   = ljMixed(i,j)%sigma
123 <          ljMixed(j,i)%sigma2  = ljMixed(i,j)%sigma2
120 >               calcLJMix("epsilon",myEpsilon_i, &
121 >               myEpsilon_j)
122 >          
123 >          ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * &
124 >               (ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6)
125 >          
126 >          
127 >          ljMixed(j,i)%sigma   = ljMixed(i,j)%sigma
128            ljMixed(j,i)%sigma6  = ljMixed(i,j)%sigma6
129 +          ljMixed(j,i)%tp6     = ljMixed(i,j)%tp6
130 +          ljMixed(j,i)%tp12    = ljMixed(i,j)%tp12
131            ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
132 +          ljMixed(j,i)%delta   = ljMixed(i,j)%delta
133            
119          
120          tmpPtr_j => tmpPtr_j%next
121          j = j + 1
134         end do
123       ! advance pointers
124       tmpPtr_i => tmpPtr_i%next
125       if (associated(tmpPtr_i)) then
126          tmpPtr_j => tmpPtr_i%next
127       endif
128      
135      end do
136 +
137      
138    end subroutine createMixingList
139    
140  
141 <  subroutine do_lj_pair(atom1, atom2, atype1, atype2, dx, dy, dz, rij, &
135 <       pot, f)
141 >  subroutine do_lj_pair(atom1, atom2, d, rij, pot, f)
142  
143      integer, intent(in) ::  atom1, atom2
144 <    real( kind = dp ), intent(in) :: dx, dy, dz, rij
144 >    real( kind = dp ), intent(in) :: rij
145      real( kind = dp ) :: pot
146      real( kind = dp ), dimension(3, getNlocal()) :: f    
147 <    type (atype), pointer :: atype1
142 <    type (atype), pointer :: atype2
147 >    real( kind = dp ), intent(in), dimension(3) :: d
148  
149      ! local Variables
150      real( kind = dp ) :: drdx, drdy, drdz
# Line 163 | Line 168 | contains
168      real( kind = dp ) :: delta
169  
170      ! Look up the correct parameters in the mixing matrix
171 <    sigma    = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma
172 <    sigma2   = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2
173 <    sigma6   = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6
174 <    epsilon  = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon
171 > #ifdef IS_MPI
172 >    sigma    = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma
173 >    sigma6   = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
174 >    tp6      = ljMixed(atid_Row(atom1),atid_col(atom2))%tp6
175 >    tp12     = ljMixed(atid_Row(atom1),atid_Col(atom2))%tp12
176 >    epsilon  = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
177 >    delta    = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
178 > #else
179 >    sigma    = ljMixed(atid(atom1),atid(atom2))%sigma
180 >    sigma6   = ljMixed(atid(atom1),atid(atom2))%sigma6
181 >    tp6      = ljMixed(atid(atom1),atid(atom2))%tp6
182 >    tp12     = ljMixed(atid(atom1),atid(atom2))%tp12
183 >    epsilon  = ljMixed(atid(atom1),atid(atom2))%epsilon
184 >    delta    = ljMixed(atid(atom1),atid(atom2))%delta
185 > #endif
186 >
187    
188 <    call getRcut(rcut,rcut2=rcut2,rcut6=rcut6)
188 >    call getRcut(rcut)
189      
190      r2 = rij * rij
191      r6 = r2 * r2 * r2
192  
193      t6  = sigma6/ r6
194      t12 = t6 * t6
195 <                                                                              
179 <    tp6 = sigma6 / rcut6
180 <    tp12 = tp6*tp6
181 <                                                                              
182 <    delta = -4.0_DP * epsilon * (tp12 - tp6)
183 <    
195 >      
196      if (rij.le.rcut) then
197        
198         pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta      
199  
200         dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
201  
202 <       drdx = -dx / rij
203 <       drdy = -dy / rij
204 <       drdz = -dz / rij
202 >       drdx = -d(1) / rij
203 >       drdy = -d(2) / rij
204 >       drdz = -d(3) / rij
205        
206         fx = dudr * drdx
207         fy = dudr * drdy
# Line 220 | Line 232 | contains
232   #endif
233        
234         if (doStress()) then
235 <          tau_Temp(1) = tau_Temp(1) + fx * dx
236 <          tau_Temp(2) = tau_Temp(2) + fx * dy
237 <          tau_Temp(3) = tau_Temp(3) + fx * dz
238 <          tau_Temp(4) = tau_Temp(4) + fy * dx
239 <          tau_Temp(5) = tau_Temp(5) + fy * dy
240 <          tau_Temp(6) = tau_Temp(6) + fy * dz
241 <          tau_Temp(7) = tau_Temp(7) + fz * dx
242 <          tau_Temp(8) = tau_Temp(8) + fz * dy
243 <          tau_Temp(9) = tau_Temp(9) + fz * dz
235 >          tau_Temp(1) = tau_Temp(1) + fx * d(1)
236 >          tau_Temp(2) = tau_Temp(2) + fx * d(2)
237 >          tau_Temp(3) = tau_Temp(3) + fx * d(3)
238 >          tau_Temp(4) = tau_Temp(4) + fy * d(1)
239 >          tau_Temp(5) = tau_Temp(5) + fy * d(2)
240 >          tau_Temp(6) = tau_Temp(6) + fy * d(3)
241 >          tau_Temp(7) = tau_Temp(7) + fz * d(1)
242 >          tau_Temp(8) = tau_Temp(8) + fz * d(2)
243 >          tau_Temp(9) = tau_Temp(9) + fz * d(3)
244            virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
245         endif
246  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines