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: LJ.F90,v 1.2 2004-10-21 15:25:30 chuckv Exp $, $Date: 2004-10-21 15:25:30 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $ |
6 |
|
7 |
module lj |
8 |
use atype_module |
9 |
use switcheroo |
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 |
integer, parameter :: DP = selected_real_kind(15) |
21 |
|
22 |
#define __FORTRAN90 |
23 |
#include "UseTheForce/fForceField.h" |
24 |
|
25 |
integer, save :: LJ_Mixing_Policy |
26 |
real(kind=DP), save :: LJ_rcut |
27 |
logical, save :: havePolicy = .false. |
28 |
logical, save :: haveCut = .false. |
29 |
logical, save :: LJ_do_shift = .false. |
30 |
|
31 |
!! Logical has lj force field module been initialized? |
32 |
|
33 |
logical, save :: LJ_FF_initialized = .false. |
34 |
|
35 |
!! Public methods and data |
36 |
public :: init_LJ_FF |
37 |
public :: setCutoffLJ |
38 |
public :: do_lj_pair |
39 |
public :: newLJtype |
40 |
|
41 |
!! structure for lj type parameters |
42 |
type, private :: ljType |
43 |
integer :: lj_ident |
44 |
real(kind=dp) :: lj_sigma |
45 |
real(kind=dp) :: lj_epsilon |
46 |
end type ljType |
47 |
|
48 |
!! List of lj type parameters |
49 |
type, private :: ljTypeList |
50 |
integer :: n_lj_types = 0 |
51 |
integer :: currentAddition = 0 |
52 |
type(ljType), pointer :: ljParams(:) => null() |
53 |
end type ljTypeList |
54 |
|
55 |
!! The list of lj Parameters |
56 |
type (ljTypeList), save :: ljParameterList |
57 |
|
58 |
|
59 |
type :: lj_mixed_params |
60 |
!! Lennard-Jones epsilon |
61 |
real ( kind = dp ) :: epsilon = 0.0_dp |
62 |
!! Lennard-Jones Sigma |
63 |
real ( kind = dp ) :: sigma = 0.0_dp |
64 |
!! Lennard-Jones Sigma to sixth |
65 |
real ( kind = dp ) :: sigma6 = 0.0_dp |
66 |
!! |
67 |
real ( kind = dp ) :: tp6 |
68 |
real ( kind = dp ) :: tp12 |
69 |
real ( kind = dp ) :: delta = 0.0_dp |
70 |
end type lj_mixed_params |
71 |
|
72 |
type (lj_mixed_params), dimension(:,:), pointer :: ljMixed |
73 |
|
74 |
|
75 |
|
76 |
contains |
77 |
|
78 |
subroutine newLJtype(ident,lj_sigma,lj_epsilon,status) |
79 |
integer,intent(in) :: ident |
80 |
real(kind=dp),intent(in) :: lj_sigma |
81 |
real(kind=dp),intent(in) :: lj_epsilon |
82 |
integer,intent(out) :: status |
83 |
|
84 |
integer,pointer :: Matchlist(:) => null() |
85 |
integer :: current |
86 |
integer :: nAtypes |
87 |
status = 0 |
88 |
|
89 |
!! Assume that atypes has already been set and get the total number of types in atypes |
90 |
|
91 |
|
92 |
|
93 |
! check to see if this is the first time into |
94 |
if (.not.associated(ljParameterList%ljParams)) then |
95 |
call getMatchingElementList(atypes, "is_lj", .true., nAtypes, MatchList) |
96 |
ljParameterList%n_lj_types = nAtypes |
97 |
if (nAtypes == 0) then |
98 |
status = -1 |
99 |
return |
100 |
end if |
101 |
allocate(ljParameterList%ljParams(nAtypes)) |
102 |
end if |
103 |
|
104 |
ljParameterList%currentAddition = ljParameterList%currentAddition + 1 |
105 |
current = ljParameterList%currentAddition |
106 |
|
107 |
! set the values for ljParameterList |
108 |
ljParameterList%ljParams(current)%lj_ident = ident |
109 |
ljParameterList%ljParams(current)%lj_epsilon = lj_epsilon |
110 |
ljParameterList%ljParams(current)%lj_sigma = lj_sigma |
111 |
|
112 |
end subroutine newLJtype |
113 |
|
114 |
subroutine init_LJ_FF(mix_Policy, status) |
115 |
integer, intent(in) :: mix_Policy |
116 |
integer, intent(out) :: status |
117 |
integer :: myStatus |
118 |
|
119 |
if (mix_Policy == LB_MIXING_RULE) then |
120 |
LJ_Mixing_Policy = LB_MIXING_RULE |
121 |
else |
122 |
if (mix_Policy == EXPLICIT_MIXING_RULE) then |
123 |
LJ_Mixing_Policy = EXPLICIT_MIXING_RULE |
124 |
else |
125 |
write(*,*) 'Unknown Mixing Policy!' |
126 |
status = -1 |
127 |
return |
128 |
endif |
129 |
endif |
130 |
|
131 |
havePolicy = .true. |
132 |
|
133 |
if (haveCut) then |
134 |
status = 0 |
135 |
call createMixingList(myStatus) |
136 |
if (myStatus /= 0) then |
137 |
status = -1 |
138 |
return |
139 |
end if |
140 |
|
141 |
LJ_FF_initialized = .true. |
142 |
end if |
143 |
|
144 |
end subroutine init_LJ_FF |
145 |
|
146 |
subroutine setCutoffLJ(rcut, do_shift, status) |
147 |
logical, intent(in):: do_shift |
148 |
integer :: status, myStatus |
149 |
real(kind=dp) :: rcut |
150 |
|
151 |
#define __FORTRAN90 |
152 |
#include "UseTheForce/fSwitchingFunction.h" |
153 |
|
154 |
status = 0 |
155 |
|
156 |
LJ_rcut = rcut |
157 |
LJ_do_shift = do_shift |
158 |
call set_switch(LJ_SWITCH, rcut, rcut) |
159 |
haveCut = .true. |
160 |
|
161 |
if (havePolicy) then |
162 |
status = 0 |
163 |
call createMixingList(myStatus) |
164 |
if (myStatus /= 0) then |
165 |
status = -1 |
166 |
return |
167 |
end if |
168 |
|
169 |
LJ_FF_initialized = .true. |
170 |
end if |
171 |
|
172 |
return |
173 |
end subroutine setCutoffLJ |
174 |
|
175 |
subroutine createMixingList(status) |
176 |
integer :: nAtypes |
177 |
integer :: status |
178 |
integer :: i |
179 |
integer :: j |
180 |
real ( kind = dp ) :: mySigma_i,mySigma_j |
181 |
real ( kind = dp ) :: myEpsilon_i,myEpsilon_j |
182 |
real ( kind = dp ) :: rcut6 |
183 |
logical :: I_isLJ, J_isLJ |
184 |
status = 0 |
185 |
|
186 |
! we only allocate this array to the number of lj_atypes |
187 |
nAtypes = size(ljParameterList%ljParams) |
188 |
if (nAtypes == 0) then |
189 |
status = -1 |
190 |
return |
191 |
end if |
192 |
|
193 |
if (.not. associated(ljMixed)) then |
194 |
allocate(ljMixed(nAtypes, nAtypes)) |
195 |
endif |
196 |
|
197 |
rcut6 = LJ_rcut**6 |
198 |
|
199 |
! This loops through all atypes, even those that don't support LJ forces. |
200 |
do i = 1, nAtypes |
201 |
|
202 |
myEpsilon_i = ljParameterList%ljParams(i)%lj_epsilon |
203 |
mySigma_i = ljParameterList%ljParams(i)%lj_sigma |
204 |
|
205 |
! do self mixing rule |
206 |
ljMixed(i,i)%sigma = mySigma_i |
207 |
|
208 |
ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6 |
209 |
|
210 |
ljMixed(i,i)%tp6 = (ljMixed(i,i)%sigma6)/rcut6 |
211 |
|
212 |
ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2 |
213 |
|
214 |
|
215 |
ljMixed(i,i)%epsilon = myEpsilon_i |
216 |
|
217 |
ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * & |
218 |
(ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6) |
219 |
|
220 |
do j = i + 1, nAtypes |
221 |
|
222 |
myEpsilon_j = ljParameterList%ljParams(j)%lj_epsilon |
223 |
mySigma_j = ljParameterList%ljParams(j)%lj_sigma |
224 |
|
225 |
|
226 |
ljMixed(i,j)%sigma = & |
227 |
calcLJMix("sigma",mySigma_i, & |
228 |
mySigma_j) |
229 |
|
230 |
ljMixed(i,j)%sigma6 = & |
231 |
(ljMixed(i,j)%sigma)**6 |
232 |
|
233 |
|
234 |
ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6 |
235 |
|
236 |
ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2 |
237 |
|
238 |
|
239 |
ljMixed(i,j)%epsilon = & |
240 |
calcLJMix("epsilon",myEpsilon_i, & |
241 |
myEpsilon_j) |
242 |
|
243 |
ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * & |
244 |
(ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6) |
245 |
|
246 |
|
247 |
ljMixed(j,i)%sigma = ljMixed(i,j)%sigma |
248 |
ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6 |
249 |
ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6 |
250 |
ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12 |
251 |
ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon |
252 |
ljMixed(j,i)%delta = ljMixed(i,j)%delta |
253 |
|
254 |
end do |
255 |
end do |
256 |
|
257 |
end subroutine createMixingList |
258 |
|
259 |
subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
260 |
pot, f, do_pot) |
261 |
|
262 |
integer, intent(in) :: atom1, atom2 |
263 |
real( kind = dp ), intent(in) :: rij, r2 |
264 |
real( kind = dp ) :: pot, sw, vpair |
265 |
real( kind = dp ), dimension(3,nLocal) :: f |
266 |
real( kind = dp ), intent(in), dimension(3) :: d |
267 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
268 |
logical, intent(in) :: do_pot |
269 |
|
270 |
! local Variables |
271 |
real( kind = dp ) :: drdx, drdy, drdz |
272 |
real( kind = dp ) :: fx, fy, fz |
273 |
real( kind = dp ) :: pot_temp, dudr |
274 |
real( kind = dp ) :: sigma6 |
275 |
real( kind = dp ) :: epsilon |
276 |
real( kind = dp ) :: r6 |
277 |
real( kind = dp ) :: t6 |
278 |
real( kind = dp ) :: t12 |
279 |
real( kind = dp ) :: delta |
280 |
integer :: id1, id2 |
281 |
|
282 |
! Look up the correct parameters in the mixing matrix |
283 |
#ifdef IS_MPI |
284 |
sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6 |
285 |
epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon |
286 |
delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta |
287 |
#else |
288 |
sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6 |
289 |
epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon |
290 |
delta = ljMixed(atid(atom1),atid(atom2))%delta |
291 |
#endif |
292 |
|
293 |
r6 = r2 * r2 * r2 |
294 |
|
295 |
t6 = sigma6/ r6 |
296 |
t12 = t6 * t6 |
297 |
|
298 |
pot_temp = 4.0E0_DP * epsilon * (t12 - t6) |
299 |
if (LJ_do_shift) then |
300 |
pot_temp = pot_temp + delta |
301 |
endif |
302 |
|
303 |
vpair = vpair + pot_temp |
304 |
|
305 |
dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij |
306 |
|
307 |
drdx = d(1) / rij |
308 |
drdy = d(2) / rij |
309 |
drdz = d(3) / rij |
310 |
|
311 |
fx = dudr * drdx |
312 |
fy = dudr * drdy |
313 |
fz = dudr * drdz |
314 |
|
315 |
|
316 |
#ifdef IS_MPI |
317 |
if (do_pot) then |
318 |
pot_Row(atom1) = pot_Row(atom1) + sw*pot_temp*0.5 |
319 |
pot_Col(atom2) = pot_Col(atom2) + sw*pot_temp*0.5 |
320 |
endif |
321 |
|
322 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
323 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
324 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
325 |
|
326 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
327 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
328 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
329 |
|
330 |
#else |
331 |
if (do_pot) pot = pot + sw*pot_temp |
332 |
|
333 |
f(1,atom1) = f(1,atom1) + fx |
334 |
f(2,atom1) = f(2,atom1) + fy |
335 |
f(3,atom1) = f(3,atom1) + fz |
336 |
|
337 |
f(1,atom2) = f(1,atom2) - fx |
338 |
f(2,atom2) = f(2,atom2) - fy |
339 |
f(3,atom2) = f(3,atom2) - fz |
340 |
#endif |
341 |
|
342 |
#ifdef IS_MPI |
343 |
id1 = AtomRowToGlobal(atom1) |
344 |
id2 = AtomColToGlobal(atom2) |
345 |
#else |
346 |
id1 = atom1 |
347 |
id2 = atom2 |
348 |
#endif |
349 |
|
350 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
351 |
|
352 |
fpair(1) = fpair(1) + fx |
353 |
fpair(2) = fpair(2) + fy |
354 |
fpair(3) = fpair(3) + fz |
355 |
|
356 |
endif |
357 |
|
358 |
return |
359 |
|
360 |
end subroutine do_lj_pair |
361 |
|
362 |
|
363 |
!! Calculates the mixing for sigma or epslon |
364 |
|
365 |
function calcLJMix(thisParam,param1,param2,status) result(myMixParam) |
366 |
character(len=*) :: thisParam |
367 |
real(kind = dp) :: param1 |
368 |
real(kind = dp) :: param2 |
369 |
real(kind = dp ) :: myMixParam |
370 |
|
371 |
integer, optional :: status |
372 |
|
373 |
myMixParam = 0.0_dp |
374 |
|
375 |
if (present(status)) status = 0 |
376 |
select case (LJ_Mixing_Policy) |
377 |
case (1) |
378 |
select case (thisParam) |
379 |
case ("sigma") |
380 |
myMixParam = 0.5_dp * (param1 + param2) |
381 |
case ("epsilon") |
382 |
myMixParam = sqrt(param1 * param2) |
383 |
case default |
384 |
status = -1 |
385 |
end select |
386 |
case default |
387 |
status = -1 |
388 |
end select |
389 |
end function calcLJMix |
390 |
|
391 |
end module lj |
392 |
|
393 |
subroutine newLJtype(ident,lj_sigma,lj_epsilon,status) |
394 |
use lj, ONLY : module_newLJtype => newLJtype |
395 |
integer, parameter :: DP = selected_real_kind(15) |
396 |
integer,intent(inout) :: ident |
397 |
real(kind=dp),intent(inout) :: lj_sigma |
398 |
real(kind=dp),intent(inout) :: lj_epsilon |
399 |
integer,intent(inout) :: status |
400 |
|
401 |
call module_newLJtype(ident,lj_sigma,lj_epsilon,status) |
402 |
|
403 |
end subroutine newLJtype |
404 |
|