1 |
|
!! Calculates Long Range forces Lennard-Jones interactions. |
2 |
– |
!! Corresponds to the force field defined in lj_FF.cpp |
2 |
|
!! @author Charles F. Vardeman II |
3 |
|
!! @author Matthew Meineke |
4 |
< |
!! @version $Id: LJ.F90,v 1.1 2004-10-20 04:02:48 gezelter Exp $, $Date: 2004-10-20 04:02:48 $, $Name: not supported by cvs2svn $, $Revision: 1.1 $ |
4 |
> |
!! @version $Id: LJ.F90,v 1.3 2004-10-21 20:15:25 gezelter Exp $, $Date: 2004-10-21 20:15:25 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $ |
5 |
|
|
6 |
|
module lj |
8 |
– |
use definitions |
7 |
|
use atype_module |
8 |
|
use switcheroo |
9 |
|
use vector_class |
10 |
|
use simulation |
11 |
+ |
use status |
12 |
|
#ifdef IS_MPI |
13 |
|
use mpiSimulation |
14 |
|
#endif |
16 |
|
|
17 |
|
implicit none |
18 |
|
PRIVATE |
20 |
– |
|
21 |
– |
#define __FORTRAN90 |
22 |
– |
#include "UseTheForce/fForceField.h" |
23 |
– |
|
24 |
– |
integer, save :: LJ_Mixing_Policy |
25 |
– |
real(kind=DP), save :: LJ_rcut |
26 |
– |
logical, save :: havePolicy = .false. |
27 |
– |
logical, save :: haveCut = .false. |
28 |
– |
logical, save :: LJ_do_shift = .false. |
19 |
|
|
20 |
< |
!! Logical has lj force field module been initialized? |
20 |
> |
integer, parameter :: DP = selected_real_kind(15) |
21 |
|
|
22 |
< |
logical, save :: LJ_FF_initialized = .false. |
22 |
> |
type, private :: LjType |
23 |
> |
integer :: ident |
24 |
> |
real(kind=dp) :: sigma |
25 |
> |
real(kind=dp) :: epsilon |
26 |
> |
end type LjType |
27 |
|
|
28 |
+ |
type(LjType), dimension(:), allocatable :: ParameterMap |
29 |
+ |
|
30 |
+ |
logical, save :: haveMixingMap = .false. |
31 |
+ |
|
32 |
+ |
type :: MixParameters |
33 |
+ |
real(kind=DP) :: sigma |
34 |
+ |
real(kind=DP) :: epsilon |
35 |
+ |
real(kind=dp) :: sigma6 |
36 |
+ |
real(kind=dp) :: tp6 |
37 |
+ |
real(kind=dp) :: tp12 |
38 |
+ |
real(kind=dp) :: delta |
39 |
+ |
end type MixParameters |
40 |
+ |
|
41 |
+ |
type(MixParameters), dimension(:,:), allocatable :: MixingMap |
42 |
+ |
|
43 |
+ |
real(kind=DP), save :: LJ_rcut |
44 |
+ |
logical, save :: have_rcut = .false. |
45 |
+ |
logical, save :: LJ_do_shift = .false. |
46 |
+ |
logical, save :: useGeometricDistanceMixing = .false. |
47 |
+ |
|
48 |
|
!! Public methods and data |
49 |
< |
public :: init_LJ_FF |
49 |
> |
|
50 |
|
public :: setCutoffLJ |
51 |
+ |
public :: useGeometricMixing |
52 |
|
public :: do_lj_pair |
53 |
+ |
public :: newLJtype |
54 |
|
|
39 |
– |
type :: lj_mixed_params |
40 |
– |
!! Lennard-Jones epsilon |
41 |
– |
real ( kind = dp ) :: epsilon = 0.0_dp |
42 |
– |
!! Lennard-Jones Sigma |
43 |
– |
real ( kind = dp ) :: sigma = 0.0_dp |
44 |
– |
!! Lennard-Jones Sigma to sixth |
45 |
– |
real ( kind = dp ) :: sigma6 = 0.0_dp |
46 |
– |
!! |
47 |
– |
real ( kind = dp ) :: tp6 |
48 |
– |
real ( kind = dp ) :: tp12 |
49 |
– |
real ( kind = dp ) :: delta = 0.0_dp |
50 |
– |
end type lj_mixed_params |
51 |
– |
|
52 |
– |
type (lj_mixed_params), dimension(:,:), pointer :: ljMixed |
53 |
– |
|
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 |
57 |
> |
subroutine newLJtype(ident, sigma, epsilon, status) |
58 |
> |
integer,intent(in) :: ident |
59 |
> |
real(kind=dp),intent(in) :: sigma |
60 |
> |
real(kind=dp),intent(in) :: epsilon |
61 |
> |
integer,intent(out) :: status |
62 |
> |
integer :: nAtypes |
63 |
> |
|
64 |
> |
status = 0 |
65 |
|
|
66 |
< |
if (mix_Policy == LB_MIXING_RULE) then |
67 |
< |
LJ_Mixing_Policy = LB_MIXING_RULE |
63 |
< |
else |
64 |
< |
if (mix_Policy == EXPLICIT_MIXING_RULE) then |
65 |
< |
LJ_Mixing_Policy = EXPLICIT_MIXING_RULE |
66 |
< |
else |
67 |
< |
write(*,*) 'Unknown Mixing Policy!' |
68 |
< |
status = -1 |
69 |
< |
return |
70 |
< |
endif |
71 |
< |
endif |
66 |
> |
!! Be simple-minded and assume that we need a ParameterMap that |
67 |
> |
!! is the same size as the total number of atom types |
68 |
|
|
69 |
< |
havePolicy = .true. |
70 |
< |
|
71 |
< |
if (haveCut) then |
72 |
< |
status = 0 |
73 |
< |
call createMixingList(myStatus) |
78 |
< |
if (myStatus /= 0) then |
69 |
> |
if (.not.allocated(ParameterMap)) then |
70 |
> |
|
71 |
> |
nAtypes = getSize(atypes) |
72 |
> |
|
73 |
> |
if (nAtypes == 0) then |
74 |
|
status = -1 |
75 |
|
return |
76 |
|
end if |
77 |
|
|
78 |
< |
LJ_FF_initialized = .true. |
78 |
> |
if (.not. allocated(ParameterMap)) then |
79 |
> |
allocate(ParameterMap(nAtypes)) |
80 |
> |
endif |
81 |
> |
|
82 |
|
end if |
83 |
< |
|
84 |
< |
end subroutine init_LJ_FF |
85 |
< |
|
83 |
> |
|
84 |
> |
if (ident .gt. size(ParameterMap)) then |
85 |
> |
status = -1 |
86 |
> |
return |
87 |
> |
endif |
88 |
> |
|
89 |
> |
! set the values for ParameterMap for this atom type: |
90 |
> |
|
91 |
> |
ParameterMap(ident)%ident = ident |
92 |
> |
ParameterMap(ident)%epsilon = epsilon |
93 |
> |
ParameterMap(ident)%sigma = sigma |
94 |
> |
|
95 |
> |
end subroutine newLJtype |
96 |
> |
|
97 |
|
subroutine setCutoffLJ(rcut, do_shift, status) |
98 |
|
logical, intent(in):: do_shift |
99 |
|
integer :: status, myStatus |
107 |
|
LJ_rcut = rcut |
108 |
|
LJ_do_shift = do_shift |
109 |
|
call set_switch(LJ_SWITCH, rcut, rcut) |
110 |
< |
haveCut = .true. |
102 |
< |
|
103 |
< |
if (havePolicy) then |
104 |
< |
status = 0 |
105 |
< |
call createMixingList(myStatus) |
106 |
< |
if (myStatus /= 0) then |
107 |
< |
status = -1 |
108 |
< |
return |
109 |
< |
end if |
110 |
< |
|
111 |
< |
LJ_FF_initialized = .true. |
112 |
< |
end if |
110 |
> |
have_rcut = .true. |
111 |
|
|
112 |
|
return |
113 |
|
end subroutine setCutoffLJ |
114 |
+ |
|
115 |
+ |
subroutine useGeometricMixing() |
116 |
+ |
useGeometricDistanceMixing = .true. |
117 |
+ |
haveMixingMap = .false. |
118 |
+ |
return |
119 |
+ |
end subroutine useGeometricMixing |
120 |
|
|
121 |
< |
subroutine createMixingList(status) |
121 |
> |
subroutine createMixingMap(status) |
122 |
|
integer :: nAtypes |
123 |
|
integer :: status |
124 |
|
integer :: i |
125 |
|
integer :: j |
126 |
< |
real ( kind = dp ) :: mySigma_i,mySigma_j |
127 |
< |
real ( kind = dp ) :: myEpsilon_i,myEpsilon_j |
126 |
> |
real ( kind = dp ) :: Sigma_i, Sigma_j |
127 |
> |
real ( kind = dp ) :: Epsilon_i, Epsilon_j |
128 |
|
real ( kind = dp ) :: rcut6 |
129 |
< |
logical :: I_isLJ, J_isLJ |
129 |
> |
|
130 |
|
status = 0 |
131 |
|
|
132 |
< |
nAtypes = getSize(atypes) |
132 |
> |
nAtypes = size(ParameterMap) |
133 |
> |
|
134 |
|
if (nAtypes == 0) then |
135 |
|
status = -1 |
136 |
|
return |
137 |
|
end if |
133 |
– |
|
134 |
– |
if (.not. associated(ljMixed)) then |
135 |
– |
allocate(ljMixed(nAtypes, nAtypes)) |
136 |
– |
endif |
138 |
|
|
139 |
< |
rcut6 = LJ_rcut**6 |
140 |
< |
|
141 |
< |
! This loops through all atypes, even those that don't support LJ forces. |
142 |
< |
do i = 1, nAtypes |
143 |
< |
|
144 |
< |
call getElementProperty(atypes, i, "is_LJ", I_isLJ) |
145 |
< |
|
146 |
< |
if (I_isLJ) then |
147 |
< |
|
148 |
< |
call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i) |
149 |
< |
call getElementProperty(atypes, i, "lj_sigma", mySigma_i) |
150 |
< |
! do self mixing rule |
151 |
< |
ljMixed(i,i)%sigma = mySigma_i |
152 |
< |
|
153 |
< |
ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6 |
154 |
< |
|
155 |
< |
ljMixed(i,i)%tp6 = (ljMixed(i,i)%sigma6)/rcut6 |
156 |
< |
|
157 |
< |
ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2 |
139 |
> |
if (.not.have_rcut) then |
140 |
> |
status = -1 |
141 |
> |
return |
142 |
> |
endif |
143 |
> |
|
144 |
> |
if (.not. allocated(MixingMap)) then |
145 |
> |
allocate(MixingMap(nAtypes, nAtypes)) |
146 |
> |
endif |
147 |
> |
|
148 |
> |
rcut6 = LJ_rcut**6 |
149 |
> |
|
150 |
> |
! This loops through all atypes, even those that don't support LJ forces. |
151 |
> |
do i = 1, nAtypes |
152 |
> |
|
153 |
> |
Epsilon_i = ParameterMap(i)%epsilon |
154 |
> |
Sigma_i = ParameterMap(i)%sigma |
155 |
> |
|
156 |
> |
! do self mixing rule |
157 |
> |
MixingMap(i,i)%sigma = Sigma_i |
158 |
> |
MixingMap(i,i)%sigma6 = Sigma_i ** 6 |
159 |
> |
MixingMap(i,i)%tp6 = (MixingMap(i,i)%sigma6)/rcut6 |
160 |
> |
MixingMap(i,i)%tp12 = (MixingMap(i,i)%tp6) ** 2 |
161 |
> |
MixingMap(i,i)%epsilon = Epsilon_i |
162 |
> |
MixingMap(i,i)%delta = -4.0_DP * MixingMap(i,i)%epsilon * & |
163 |
> |
(MixingMap(i,i)%tp12 - MixingMap(i,i)%tp6) |
164 |
> |
|
165 |
> |
do j = i + 1, nAtypes |
166 |
|
|
167 |
+ |
Epsilon_j = ParameterMap(j)%epsilon |
168 |
+ |
Sigma_j = ParameterMap(j)%sigma |
169 |
|
|
170 |
< |
ljMixed(i,i)%epsilon = myEpsilon_i |
170 |
> |
! only the distance parameter uses different mixing policies |
171 |
> |
if (useGeometricDistanceMixing) then |
172 |
> |
! only for OPLS as far as we can tell |
173 |
> |
MixingMap(i,j)%sigma = dsqrt(Sigma_i * Sigma_j) |
174 |
> |
else |
175 |
> |
! everyone else |
176 |
> |
MixingMap(i,j)%sigma = 0.5_dp * (Sigma_i + Sigma_j) |
177 |
> |
endif |
178 |
|
|
179 |
< |
ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * & |
180 |
< |
(ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6) |
179 |
> |
! energy parameter is always geometric mean: |
180 |
> |
MixingMap(i,j)%epsilon = dsqrt(Epsilon_i * Epsilon_j) |
181 |
> |
|
182 |
> |
MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6 |
183 |
> |
MixingMap(i,j)%tp6 = MixingMap(i,j)%sigma6/rcut6 |
184 |
> |
MixingMap(i,j)%tp12 = (MixingMap(i,j)%tp6) ** 2 |
185 |
|
|
186 |
< |
do j = i + 1, nAtypes |
187 |
< |
|
188 |
< |
call getElementProperty(atypes, j, "is_LJ", J_isLJ) |
189 |
< |
|
190 |
< |
if (J_isLJ) then |
191 |
< |
|
192 |
< |
call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j) |
193 |
< |
call getElementProperty(atypes,j,"lj_sigma", mySigma_j) |
194 |
< |
|
195 |
< |
ljMixed(i,j)%sigma = & |
196 |
< |
calcLJMix("sigma",mySigma_i, & |
175 |
< |
mySigma_j) |
176 |
< |
|
177 |
< |
ljMixed(i,j)%sigma6 = & |
178 |
< |
(ljMixed(i,j)%sigma)**6 |
179 |
< |
|
180 |
< |
|
181 |
< |
ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6 |
182 |
< |
|
183 |
< |
ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2 |
184 |
< |
|
185 |
< |
|
186 |
< |
ljMixed(i,j)%epsilon = & |
187 |
< |
calcLJMix("epsilon",myEpsilon_i, & |
188 |
< |
myEpsilon_j) |
189 |
< |
|
190 |
< |
ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * & |
191 |
< |
(ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6) |
192 |
< |
|
193 |
< |
|
194 |
< |
ljMixed(j,i)%sigma = ljMixed(i,j)%sigma |
195 |
< |
ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6 |
196 |
< |
ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6 |
197 |
< |
ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12 |
198 |
< |
ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon |
199 |
< |
ljMixed(j,i)%delta = ljMixed(i,j)%delta |
200 |
< |
endif |
201 |
< |
end do |
202 |
< |
endif |
186 |
> |
MixingMap(i,j)%delta = -4.0_DP * MixingMap(i,j)%epsilon * & |
187 |
> |
(MixingMap(i,j)%tp12 - MixingMap(i,j)%tp6) |
188 |
> |
|
189 |
> |
MixingMap(j,i)%sigma = MixingMap(i,j)%sigma |
190 |
> |
MixingMap(j,i)%sigma6 = MixingMap(i,j)%sigma6 |
191 |
> |
MixingMap(j,i)%tp6 = MixingMap(i,j)%tp6 |
192 |
> |
MixingMap(j,i)%tp12 = MixingMap(i,j)%tp12 |
193 |
> |
MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon |
194 |
> |
MixingMap(j,i)%delta = MixingMap(i,j)%delta |
195 |
> |
|
196 |
> |
end do |
197 |
|
end do |
198 |
|
|
199 |
< |
end subroutine createMixingList |
200 |
< |
|
199 |
> |
end subroutine createMixingMap |
200 |
> |
|
201 |
|
subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
202 |
|
pot, f, do_pot) |
203 |
|
|
219 |
|
real( kind = dp ) :: t6 |
220 |
|
real( kind = dp ) :: t12 |
221 |
|
real( kind = dp ) :: delta |
222 |
< |
integer :: id1, id2 |
222 |
> |
integer :: id1, id2, localError |
223 |
|
|
224 |
+ |
if (.not.haveMixingMap) then |
225 |
+ |
localError = 0 |
226 |
+ |
call createMixingMap(localError) |
227 |
+ |
if ( localError .ne. 0 ) then |
228 |
+ |
call handleError("LJ", "MixingMap creation failed!") |
229 |
+ |
return |
230 |
+ |
end if |
231 |
+ |
endif |
232 |
+ |
|
233 |
|
! Look up the correct parameters in the mixing matrix |
234 |
|
#ifdef IS_MPI |
235 |
< |
sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6 |
236 |
< |
epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon |
237 |
< |
delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta |
235 |
> |
sigma6 = MixingMap(atid_Row(atom1),atid_Col(atom2))%sigma6 |
236 |
> |
epsilon = MixingMap(atid_Row(atom1),atid_Col(atom2))%epsilon |
237 |
> |
delta = MixingMap(atid_Row(atom1),atid_Col(atom2))%delta |
238 |
|
#else |
239 |
< |
sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6 |
240 |
< |
epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon |
241 |
< |
delta = ljMixed(atid(atom1),atid(atom2))%delta |
239 |
> |
sigma6 = MixingMap(atid(atom1),atid(atom2))%sigma6 |
240 |
> |
epsilon = MixingMap(atid(atom1),atid(atom2))%epsilon |
241 |
> |
delta = MixingMap(atid(atom1),atid(atom2))%delta |
242 |
|
#endif |
243 |
|
|
244 |
|
r6 = r2 * r2 * r2 |
312 |
|
|
313 |
|
|
314 |
|
!! Calculates the mixing for sigma or epslon |
315 |
+ |
|
316 |
+ |
end module lj |
317 |
+ |
|
318 |
+ |
subroutine newLJtype(ident, sigma, epsilon, status) |
319 |
+ |
use lj, ONLY : module_newLJtype => newLJtype |
320 |
+ |
integer, parameter :: DP = selected_real_kind(15) |
321 |
+ |
integer,intent(inout) :: ident |
322 |
+ |
real(kind=dp),intent(inout) :: sigma |
323 |
+ |
real(kind=dp),intent(inout) :: epsilon |
324 |
+ |
integer,intent(inout) :: status |
325 |
|
|
326 |
< |
function calcLJMix(thisParam,param1,param2,status) result(myMixParam) |
327 |
< |
character(len=*) :: thisParam |
328 |
< |
real(kind = dp) :: param1 |
316 |
< |
real(kind = dp) :: param2 |
317 |
< |
real(kind = dp ) :: myMixParam |
326 |
> |
call module_newLJtype(ident, sigma, epsilon, status) |
327 |
> |
|
328 |
> |
end subroutine newLJtype |
329 |
|
|
330 |
< |
integer, optional :: status |
331 |
< |
|
321 |
< |
myMixParam = 0.0_dp |
322 |
< |
|
323 |
< |
if (present(status)) status = 0 |
324 |
< |
select case (LJ_Mixing_Policy) |
325 |
< |
case (1) |
326 |
< |
select case (thisParam) |
327 |
< |
case ("sigma") |
328 |
< |
myMixParam = 0.5_dp * (param1 + param2) |
329 |
< |
case ("epsilon") |
330 |
< |
myMixParam = sqrt(param1 * param2) |
331 |
< |
case default |
332 |
< |
status = -1 |
333 |
< |
end select |
334 |
< |
case default |
335 |
< |
status = -1 |
336 |
< |
end select |
337 |
< |
end function calcLJMix |
330 |
> |
subroutine useGeometricMixing() |
331 |
> |
use lj, ONLY: module_useGeometricMixing => useGeometricMixing |
332 |
|
|
333 |
< |
end module lj |
333 |
> |
call module_useGeometricMixing() |
334 |
> |
return |
335 |
> |
end subroutine useGeometricMixing |