ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_LJ_FF.F90
Revision: 317
Committed: Tue Mar 11 23:13:06 2003 UTC (21 years, 3 months ago) by gezelter
File size: 7639 byte(s)
Log Message:
Bug fixes

File Contents

# Content
1 !! Calculates Long Range forces Lennard-Jones interactions.
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.4 2003-03-11 23:13:06 gezelter Exp $, $Date: 2003-03-11 23:13:06 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
6
7
8
9 module lj
10 use simulation
11 use definitions
12 use forceGlobals
13 use atype_typedefs
14 use vector_class
15 #ifdef IS_MPI
16 use mpiSimulation
17 #endif
18 implicit none
19 PRIVATE
20
21 !! Logical has lj force field module been initialized?
22 logical, save :: isljFFinit = .false.
23
24
25 !! Public methods and data
26 public :: init_ljFF
27 public :: do_lj_pair
28
29 type :: lj_mixed_params
30 !! Lennard-Jones epsilon
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
36 !! Lennard-Jones Sigma to sixth
37 real ( kind = dp ) :: sigma6 = 0.0_dp
38 end type lj_mixed_params
39
40 type (lj_mixed_params), dimension(:,:), pointer :: ljMixed
41
42 contains
43
44 subroutine init_ljFF(ListHead,status)
45 type (atype),pointer :: ListHead
46 integer :: nTypes
47 integer, intent(out) :: status
48
49 integer :: myStatus
50
51 status = 0
52 call createMixingList(ListHead,myStatus)
53 if (myStatus /= 0) then
54 status = -1
55 return
56 end if
57
58 end subroutine init_ljFF
59
60 subroutine createMixingList(ListHead,status)
61 type (atype), pointer :: ListHead
62 integer :: listSize
63 integer :: status
64 integer :: i
65 integer :: j
66 type (atype), pointer :: tmpPtr_i => null()
67 type (atype), pointer :: tmpPtr_j => null()
68 status = 0
69
70 listSize = get_this_ListLen(ListHead)
71 if (listSize == 0) then
72 status = -1
73 return
74 end if
75
76
77 if (.not. associated(ljMixed)) then
78 allocate(ljMixed(listSize,listSize))
79 else
80 status = -1
81 return
82 end if
83
84
85
86 tmpPtr_i => ListHead
87 tmpPtr_j => tmpPtr_i%next
88 do while (associated(tmpPtr_i))
89 i = i + 1
90 ! do self mixing rule
91 ljMixed(i,i)%sigma = tmpPtr_i%lj_sigma
92
93 ljMixed(i,i)%sigma2 = (ljMixed(i,i)%sigma) ** 2
94
95 ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
96 ljMixed(i,i)%epsilon = tmpPtr_i%lj_epsilon
97
98 j = i + 1
99 do while (associated(tmpPtr_j))
100
101 ljMixed(i,j)%sigma = &
102 calcLJMix("sigma",tmpPtr_i%lj_sigma, &
103 tmpPtr_j%lj_sigma)
104
105 ljMixed(i,j)%sigma2 = &
106 (ljMixed(i,j)%sigma)**2
107
108 ljMixed(i,j)%sigma6 = &
109 (ljMixed(i,j)%sigma)**6
110
111 ljMixed(i,j)%epsilon = &
112 calcLJMix("epsilon",tmpPtr_i%lj_epsilon, &
113 tmpPtr_j%lj_epsilon)
114 ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
115 ljMixed(j,i)%sigma2 = ljMixed(i,j)%sigma2
116 ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
117 ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
118
119
120 tmpPtr_j => tmpPtr_j%next
121 j = j + 1
122 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
129 end do
130
131 end subroutine createMixingList
132
133
134 subroutine do_lj_pair(atom1, atom2, d, rij, pot, f)
135
136 integer, intent(in) :: atom1, atom2
137 real( kind = dp ), intent(in) :: rij
138 real( kind = dp ) :: pot
139 real( kind = dp ), dimension(3, getNlocal()) :: f
140 real( kind = dp ), intent(in), dimension(3) :: d
141
142 ! local Variables
143 real( kind = dp ) :: drdx, drdy, drdz
144 real( kind = dp ) :: fx, fy, fz
145 real( kind = dp ) :: pot_temp, dudr
146 real( kind = dp ) :: sigma
147 real( kind = dp ) :: sigma2
148 real( kind = dp ) :: sigma6
149 real( kind = dp ) :: epsilon
150
151 real( kind = dp ) :: rcut
152 real( kind = dp ) :: rcut2
153 real( kind = dp ) :: rcut6
154 real( kind = dp ) :: r2
155 real( kind = dp ) :: r6
156
157 real( kind = dp ) :: t6
158 real( kind = dp ) :: t12
159 real( kind = dp ) :: tp6
160 real( kind = dp ) :: tp12
161 real( kind = dp ) :: delta
162
163 ! Look up the correct parameters in the mixing matrix
164 sigma = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma
165 sigma2 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2
166 sigma6 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6
167 epsilon = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon
168
169 call getRcut(rcut,rcut2=rcut2,rcut6=rcut6)
170
171 r2 = rij * rij
172 r6 = r2 * r2 * r2
173
174 t6 = sigma6/ r6
175 t12 = t6 * t6
176
177 tp6 = sigma6 / rcut6
178 tp12 = tp6*tp6
179
180 delta = -4.0_DP * epsilon * (tp12 - tp6)
181
182 if (rij.le.rcut) then
183
184 pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta
185
186 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
187
188 drdx = -d(1) / rij
189 drdy = -d(2) / rij
190 drdz = -d(3) / rij
191
192 fx = dudr * drdx
193 fy = dudr * drdy
194 fz = dudr * drdz
195
196 #ifdef IS_MPI
197 pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
198 pot_Col(atom2) = pot_Col(atom2) + pot_temp*0.5
199
200 f_Row(1,atom1) = f_Row(1,atom1) + fx
201 f_Row(2,atom1) = f_Row(2,atom1) + fy
202 f_Row(3,atom1) = f_Row(3,atom1) + fz
203
204 f_Col(1,atom2) = f_Col(1,atom2) - fx
205 f_Col(2,atom2) = f_Col(2,atom2) - fy
206 f_Col(3,atom2) = f_Col(3,atom2) - fz
207
208 #else
209 pot = pot + pot_temp
210
211 f(1,atom1) = f(1,atom1) + fx
212 f(2,atom1) = f(2,atom1) + fy
213 f(3,atom1) = f(3,atom1) + fz
214
215 f(1,atom2) = f(1,atom2) - fx
216 f(2,atom2) = f(2,atom2) - fy
217 f(3,atom2) = f(3,atom2) - fz
218 #endif
219
220 if (doStress()) then
221 tau_Temp(1) = tau_Temp(1) + fx * d(1)
222 tau_Temp(2) = tau_Temp(2) + fx * d(2)
223 tau_Temp(3) = tau_Temp(3) + fx * d(3)
224 tau_Temp(4) = tau_Temp(4) + fy * d(1)
225 tau_Temp(5) = tau_Temp(5) + fy * d(2)
226 tau_Temp(6) = tau_Temp(6) + fy * d(3)
227 tau_Temp(7) = tau_Temp(7) + fz * d(1)
228 tau_Temp(8) = tau_Temp(8) + fz * d(2)
229 tau_Temp(9) = tau_Temp(9) + fz * d(3)
230 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
231 endif
232
233 endif
234 return
235 end subroutine do_lj_pair
236
237
238 !! Calculates the mixing for sigma or epslon based on character
239 !! string for initialzition of mixing array.
240 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
241 character(len=*) :: thisParam
242 real(kind = dp) :: param1
243 real(kind = dp) :: param2
244 real(kind = dp ) :: myMixParam
245 character(len = getStringLen()) :: thisMixingRule
246 integer, optional :: status
247
248 !! get the mixing rules from the simulation
249 thisMixingRule = returnMixingRules()
250 myMixParam = 0.0_dp
251
252 if (present(status)) status = 0
253 select case (thisMixingRule)
254 case ("standard")
255 select case (thisParam)
256 case ("sigma")
257 myMixParam = 0.5_dp * (param1 + param2)
258 case ("epsilon")
259 myMixParam = sqrt(param1 * param2)
260 case default
261 status = -1
262 end select
263 case("LJglass")
264 case default
265 status = -1
266 end select
267 end function calcLJMix
268
269 end module lj