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