ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 1041
Committed: Mon Feb 9 14:48:57 2004 UTC (20 years, 4 months ago) by chrisfen
File size: 9069 byte(s)
Log Message:
Stripped out the hardwired LJ_rcut

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.17 2004-02-09 14:48:57 chrisfen Exp $, $Date: 2004-02-09 14:48:57 $, $Name: not supported by cvs2svn $, $Revision: 1.17 $
6
7 module lj
8 use definitions
9 use atype_module
10 use vector_class
11 use simulation
12 #ifdef IS_MPI
13 use mpiSimulation
14 #endif
15 use force_globals
16
17 implicit none
18 PRIVATE
19
20 #define __FORTRAN90
21 #include "fForceField.h"
22
23 integer, save :: LJ_Mixing_Policy
24 real(kind=DP), save :: LJ_rcut
25 logical, save :: havePolicy = .false., haveCut = .false.
26
27
28 !! Logical has lj force field module been initialized?
29
30 logical, save :: LJ_FF_initialized = .false.
31
32 !! Public methods and data
33 public :: init_LJ_FF
34 public :: setCutoffLJ
35 public :: do_lj_pair
36
37 type :: lj_mixed_params
38 !! Lennard-Jones epsilon
39 real ( kind = dp ) :: epsilon = 0.0_dp
40 !! Lennard-Jones Sigma
41 real ( kind = dp ) :: sigma = 0.0_dp
42 !! Lennard-Jones Sigma to sixth
43 real ( kind = dp ) :: sigma6 = 0.0_dp
44 !!
45 real ( kind = dp ) :: tp6
46 real ( kind = dp ) :: tp12
47
48 real ( kind = dp ) :: delta = 0.0_dp
49
50
51 end type lj_mixed_params
52
53 type (lj_mixed_params), dimension(:,:), pointer :: ljMixed
54
55 contains
56
57 subroutine init_LJ_FF(mix_Policy, status)
58 integer, intent(in) :: mix_Policy
59 integer, intent(out) :: status
60 integer :: myStatus
61
62 if (mix_Policy == LB_MIXING_RULE) then
63 LJ_Mixing_Policy = LB_MIXING_RULE
64 else
65 if (mix_Policy == EXPLICIT_MIXING_RULE) then
66 LJ_Mixing_Policy = EXPLICIT_MIXING_RULE
67 else
68 write(*,*) 'Unknown Mixing Policy!'
69 status = -1
70 return
71 endif
72 endif
73
74 havePolicy = .true.
75
76 if (haveCut) then
77 status = 0
78 call createMixingList(myStatus)
79 if (myStatus /= 0) then
80 status = -1
81 return
82 end if
83
84 LJ_FF_initialized = .true.
85 end if
86
87 end subroutine init_LJ_FF
88
89 subroutine setCutoffLJ(rcut, status)
90 integer :: status, myStatus
91 real(kind=dp) :: rcut
92
93 status = 0
94
95 LJ_rcut = rcut
96 !!$ ! ATTENTION! This is a hardwiring of rcut!
97 !!$ LJ_rcut = 9.0d0
98 haveCut = .true.
99
100 if (havePolicy) then
101 status = 0
102 call createMixingList(myStatus)
103 if (myStatus /= 0) then
104 status = -1
105 return
106 end if
107
108 LJ_FF_initialized = .true.
109 end if
110
111 return
112 end subroutine setCutoffLJ
113
114 subroutine createMixingList(status)
115 integer :: nAtypes
116 integer :: status
117 integer :: i
118 integer :: j
119 real ( kind = dp ) :: mySigma_i,mySigma_j
120 real ( kind = dp ) :: myEpsilon_i,myEpsilon_j
121 real ( kind = dp ) :: rcut6
122 status = 0
123
124 nAtypes = getSize(atypes)
125 if (nAtypes == 0) then
126 status = -1
127 return
128 end if
129
130 if (.not. associated(ljMixed)) then
131 allocate(ljMixed(nAtypes, nAtypes))
132 endif
133
134 rcut6 = LJ_rcut**6
135
136 ! This loops through all atypes, even those that don't support LJ forces.
137 do i = 1, nAtypes
138
139 call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i)
140 call getElementProperty(atypes, i, "lj_sigma", mySigma_i)
141 ! do self mixing rule
142 ljMixed(i,i)%sigma = mySigma_i
143
144 ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
145
146 ljMixed(i,i)%tp6 = (ljMixed(i,i)%sigma6)/rcut6
147
148 ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2
149
150
151 ljMixed(i,i)%epsilon = myEpsilon_i
152
153 ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
154 (ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6)
155
156 do j = i + 1, nAtypes
157 call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j)
158 call getElementProperty(atypes,j,"lj_sigma", mySigma_j)
159
160 ljMixed(i,j)%sigma = &
161 calcLJMix("sigma",mySigma_i, &
162 mySigma_j)
163
164 ljMixed(i,j)%sigma6 = &
165 (ljMixed(i,j)%sigma)**6
166
167
168 ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6
169
170 ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2
171
172
173 ljMixed(i,j)%epsilon = &
174 calcLJMix("epsilon",myEpsilon_i, &
175 myEpsilon_j)
176
177 ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * &
178 (ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6)
179
180
181 ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
182 ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
183 ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6
184 ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12
185 ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
186 ljMixed(j,i)%delta = ljMixed(i,j)%delta
187
188 end do
189 end do
190
191 end subroutine createMixingList
192
193 subroutine do_lj_pair(atom1, atom2, d, rij, r2, pot, f, do_pot, do_stress)
194
195 integer, intent(in) :: atom1, atom2
196 real( kind = dp ), intent(in) :: rij, r2
197 real( kind = dp ) :: pot
198 real( kind = dp ), dimension(3,nLocal) :: f
199 real( kind = dp ), intent(in), dimension(3) :: d
200 logical, intent(in) :: do_pot, do_stress
201
202 ! local Variables
203 real( kind = dp ) :: drdx, drdy, drdz
204 real( kind = dp ) :: fx, fy, fz
205 real( kind = dp ) :: pot_temp, dudr
206 real( kind = dp ) :: sigma6
207 real( kind = dp ) :: epsilon
208 real( kind = dp ) :: r6
209 real( kind = dp ) :: t6
210 real( kind = dp ) :: t12
211 real( kind = dp ) :: delta
212 integer :: id1, id2
213
214 if (rij.lt.LJ_rcut) then
215
216 ! Look up the correct parameters in the mixing matrix
217 #ifdef IS_MPI
218 sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
219 epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
220 delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
221 #else
222 sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6
223 epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon
224 delta = ljMixed(atid(atom1),atid(atom2))%delta
225 #endif
226
227 r6 = r2 * r2 * r2
228
229 t6 = sigma6/ r6
230 t12 = t6 * t6
231
232 pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta
233
234 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
235
236 drdx = d(1) / rij
237 drdy = d(2) / rij
238 drdz = d(3) / rij
239
240 fx = dudr * drdx
241 fy = dudr * drdy
242 fz = dudr * drdz
243
244 #ifdef IS_MPI
245 if (do_pot) then
246 pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
247 pot_Col(atom2) = pot_Col(atom2) + pot_temp*0.5
248 endif
249
250 f_Row(1,atom1) = f_Row(1,atom1) + fx
251 f_Row(2,atom1) = f_Row(2,atom1) + fy
252 f_Row(3,atom1) = f_Row(3,atom1) + fz
253
254 f_Col(1,atom2) = f_Col(1,atom2) - fx
255 f_Col(2,atom2) = f_Col(2,atom2) - fy
256 f_Col(3,atom2) = f_Col(3,atom2) - fz
257
258 #else
259 if (do_pot) pot = pot + pot_temp
260
261 f(1,atom1) = f(1,atom1) + fx
262 f(2,atom1) = f(2,atom1) + fy
263 f(3,atom1) = f(3,atom1) + fz
264
265 f(1,atom2) = f(1,atom2) - fx
266 f(2,atom2) = f(2,atom2) - fy
267 f(3,atom2) = f(3,atom2) - fz
268 #endif
269
270 if (do_stress) then
271
272 #ifdef IS_MPI
273 id1 = tagRow(atom1)
274 id2 = tagColumn(atom2)
275 #else
276 id1 = atom1
277 id2 = atom2
278 #endif
279
280 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
281
282 ! because the d vector is the rj - ri vector, and
283 ! because fx, fy, fz are the force on atom i, we need a
284 ! negative sign here:
285
286 tau_Temp(1) = tau_Temp(1) - d(1) * fx
287 tau_Temp(2) = tau_Temp(2) - d(1) * fy
288 tau_Temp(3) = tau_Temp(3) - d(1) * fz
289 tau_Temp(4) = tau_Temp(4) - d(2) * fx
290 tau_Temp(5) = tau_Temp(5) - d(2) * fy
291 tau_Temp(6) = tau_Temp(6) - d(2) * fz
292 tau_Temp(7) = tau_Temp(7) - d(3) * fx
293 tau_Temp(8) = tau_Temp(8) - d(3) * fy
294 tau_Temp(9) = tau_Temp(9) - d(3) * fz
295
296 virial_Temp = virial_Temp + &
297 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
298
299 endif
300 endif
301
302 endif
303 return
304
305 end subroutine do_lj_pair
306
307
308 !! Calculates the mixing for sigma or epslon
309
310 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
311 character(len=*) :: thisParam
312 real(kind = dp) :: param1
313 real(kind = dp) :: param2
314 real(kind = dp ) :: myMixParam
315
316 integer, optional :: status
317
318 myMixParam = 0.0_dp
319
320 if (present(status)) status = 0
321 select case (LJ_Mixing_Policy)
322 case (1)
323 select case (thisParam)
324 case ("sigma")
325 myMixParam = 0.5_dp * (param1 + param2)
326 case ("epsilon")
327 myMixParam = sqrt(param1 * param2)
328 case default
329 status = -1
330 end select
331 case default
332 status = -1
333 end select
334 end function calcLJMix
335
336 end module lj