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