ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 834
Committed: Tue Oct 28 20:09:45 2003 UTC (20 years, 8 months ago) by gezelter
File size: 8928 byte(s)
Log Message:
Compatibility 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.13 2003-10-28 20:09:37 gezelter Exp $, $Date: 2003-10-28 20:09:37 $, $Name: not supported by cvs2svn $, $Revision: 1.13 $
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 haveCut = .true.
97
98 if (havePolicy) then
99 status = 0
100 call createMixingList(myStatus)
101 if (myStatus /= 0) then
102 status = -1
103 return
104 end if
105
106 LJ_FF_initialized = .true.
107 end if
108
109 return
110 end subroutine setCutoffLJ
111
112 subroutine createMixingList(status)
113 integer :: nAtypes
114 integer :: status
115 integer :: i
116 integer :: j
117 real ( kind = dp ) :: mySigma_i,mySigma_j
118 real ( kind = dp ) :: myEpsilon_i,myEpsilon_j
119 real ( kind = dp ) :: rcut6
120 status = 0
121
122 nAtypes = getSize(atypes)
123 if (nAtypes == 0) then
124 status = -1
125 return
126 end if
127
128 if (.not. associated(ljMixed)) then
129 allocate(ljMixed(nAtypes, nAtypes))
130 endif
131
132 rcut6 = LJ_rcut**6
133
134 do i = 1, nAtypes
135
136 call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i)
137 call getElementProperty(atypes, i, "lj_sigma", mySigma_i)
138 ! do self mixing rule
139 ljMixed(i,i)%sigma = mySigma_i
140
141 ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
142
143 ljMixed(i,i)%tp6 = (ljMixed(i,i)%sigma6)/rcut6
144
145 ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2
146
147
148 ljMixed(i,i)%epsilon = myEpsilon_i
149
150 ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
151 (ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6)
152
153 do j = i + 1, nAtypes
154 call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j)
155 call getElementProperty(atypes,j,"lj_sigma", mySigma_j)
156
157 ljMixed(i,j)%sigma = &
158 calcLJMix("sigma",mySigma_i, &
159 mySigma_j)
160
161 ljMixed(i,j)%sigma6 = &
162 (ljMixed(i,j)%sigma)**6
163
164
165 ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6
166
167 ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2
168
169
170 ljMixed(i,j)%epsilon = &
171 calcLJMix("epsilon",myEpsilon_i, &
172 myEpsilon_j)
173
174 ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * &
175 (ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6)
176
177
178 ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
179 ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
180 ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6
181 ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12
182 ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
183 ljMixed(j,i)%delta = ljMixed(i,j)%delta
184
185 end do
186 end do
187
188 end subroutine createMixingList
189
190 subroutine do_lj_pair(atom1, atom2, d, rij, r2, pot, f, do_pot, do_stress)
191
192 integer, intent(in) :: atom1, atom2
193 real( kind = dp ), intent(in) :: rij, r2
194 real( kind = dp ) :: pot
195 real( kind = dp ), dimension(3,getNlocal()) :: f
196 real( kind = dp ), intent(in), dimension(3) :: d
197 logical, intent(in) :: do_pot, do_stress
198
199 ! local Variables
200 real( kind = dp ) :: drdx, drdy, drdz
201 real( kind = dp ) :: fx, fy, fz
202 real( kind = dp ) :: pot_temp, dudr
203 real( kind = dp ) :: sigma6
204 real( kind = dp ) :: epsilon
205 real( kind = dp ) :: r6
206 real( kind = dp ) :: t6
207 real( kind = dp ) :: t12
208 real( kind = dp ) :: delta
209 integer :: id1, id2
210
211
212 if (rij.lt.LJ_rcut) then
213
214 ! Look up the correct parameters in the mixing matrix
215 #ifdef IS_MPI
216 sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
217 epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
218 delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
219 #else
220 sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6
221 epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon
222 delta = ljMixed(atid(atom1),atid(atom2))%delta
223 #endif
224
225 r6 = r2 * r2 * r2
226
227 t6 = sigma6/ r6
228 t12 = t6 * t6
229
230 pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta
231
232 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
233
234 drdx = d(1) / rij
235 drdy = d(2) / rij
236 drdz = d(3) / rij
237
238 fx = dudr * drdx
239 fy = dudr * drdy
240 fz = dudr * drdz
241
242 #ifdef IS_MPI
243 if (do_pot) then
244 pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
245 pot_Col(atom2) = pot_Col(atom2) + pot_temp*0.5
246 endif
247
248 f_Row(1,atom1) = f_Row(1,atom1) + fx
249 f_Row(2,atom1) = f_Row(2,atom1) + fy
250 f_Row(3,atom1) = f_Row(3,atom1) + fz
251
252 f_Col(1,atom2) = f_Col(1,atom2) - fx
253 f_Col(2,atom2) = f_Col(2,atom2) - fy
254 f_Col(3,atom2) = f_Col(3,atom2) - fz
255
256 #else
257 if (do_pot) pot = pot + pot_temp
258
259 f(1,atom1) = f(1,atom1) + fx
260 f(2,atom1) = f(2,atom1) + fy
261 f(3,atom1) = f(3,atom1) + fz
262
263 f(1,atom2) = f(1,atom2) - fx
264 f(2,atom2) = f(2,atom2) - fy
265 f(3,atom2) = f(3,atom2) - fz
266 #endif
267
268 if (do_stress) then
269
270 #ifdef IS_MPI
271 id1 = tagRow(atom1)
272 id2 = tagColumn(atom2)
273 #else
274 id1 = atom1
275 id2 = atom2
276 #endif
277
278 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
279
280 ! because the d vector is the rj - ri vector, and
281 ! because fx, fy, fz are the force on atom i, we need a
282 ! negative sign here:
283
284 tau_Temp(1) = tau_Temp(1) - d(1) * fx
285 tau_Temp(2) = tau_Temp(2) - d(1) * fy
286 tau_Temp(3) = tau_Temp(3) - d(1) * fz
287 tau_Temp(4) = tau_Temp(4) - d(2) * fx
288 tau_Temp(5) = tau_Temp(5) - d(2) * fy
289 tau_Temp(6) = tau_Temp(6) - d(2) * fz
290 tau_Temp(7) = tau_Temp(7) - d(3) * fx
291 tau_Temp(8) = tau_Temp(8) - d(3) * fy
292 tau_Temp(9) = tau_Temp(9) - d(3) * fz
293
294 virial_Temp = virial_Temp + &
295 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
296
297 endif
298 endif
299
300 endif
301 return
302
303 end subroutine do_lj_pair
304
305
306 !! Calculates the mixing for sigma or epslon
307
308 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
309 character(len=*) :: thisParam
310 real(kind = dp) :: param1
311 real(kind = dp) :: param2
312 real(kind = dp ) :: myMixParam
313
314 integer, optional :: status
315
316 myMixParam = 0.0_dp
317
318 if (present(status)) status = 0
319 select case (LJ_Mixing_Policy)
320 case (1)
321 select case (thisParam)
322 case ("sigma")
323 myMixParam = 0.5_dp * (param1 + param2)
324 case ("epsilon")
325 myMixParam = sqrt(param1 * param2)
326 case default
327 status = -1
328 end select
329 case default
330 status = -1
331 end select
332 end function calcLJMix
333
334 end module lj