1 |
+ |
!! |
2 |
+ |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
+ |
!! |
4 |
+ |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
+ |
!! non-exclusive, royalty free, license to use, modify and |
6 |
+ |
!! redistribute this software in source and binary code form, provided |
7 |
+ |
!! that the following conditions are met: |
8 |
+ |
!! |
9 |
+ |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
+ |
!! publication of scientific results based in part on use of the |
11 |
+ |
!! program. An acceptable form of acknowledgement is citation of |
12 |
+ |
!! the article in which the program was described (Matthew |
13 |
+ |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
+ |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
+ |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
+ |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
+ |
!! |
18 |
+ |
!! 2. Redistributions of source code must retain the above copyright |
19 |
+ |
!! notice, this list of conditions and the following disclaimer. |
20 |
+ |
!! |
21 |
+ |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
+ |
!! notice, this list of conditions and the following disclaimer in the |
23 |
+ |
!! documentation and/or other materials provided with the |
24 |
+ |
!! distribution. |
25 |
+ |
!! |
26 |
+ |
!! This software is provided "AS IS," without a warranty of any |
27 |
+ |
!! kind. All express or implied conditions, representations and |
28 |
+ |
!! warranties, including any implied warranty of merchantability, |
29 |
+ |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
+ |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
+ |
!! be liable for any damages suffered by licensee as a result of |
32 |
+ |
!! using, modifying or distributing the software or its |
33 |
+ |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
+ |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
+ |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
+ |
!! damages, however caused and regardless of the theory of liability, |
37 |
+ |
!! arising out of the use of or inability to use software, even if the |
38 |
+ |
!! University of Notre Dame has been advised of the possibility of |
39 |
+ |
!! such damages. |
40 |
+ |
!! |
41 |
+ |
|
42 |
+ |
|
43 |
|
!! Calculates Long Range forces Lennard-Jones interactions. |
2 |
– |
!! Corresponds to the force field defined in lj_FF.cpp |
44 |
|
!! @author Charles F. Vardeman II |
45 |
|
!! @author Matthew Meineke |
46 |
< |
!! @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 $ |
46 |
> |
!! @version $Id: LJ.F90,v 1.13 2005-08-09 22:33:48 gezelter Exp $, $Date: 2005-08-09 22:33:48 $, $Name: not supported by cvs2svn $, $Revision: 1.13 $ |
47 |
|
|
48 |
+ |
|
49 |
|
module lj |
50 |
|
use atype_module |
9 |
– |
use switcheroo |
51 |
|
use vector_class |
52 |
|
use simulation |
53 |
+ |
use status |
54 |
|
#ifdef IS_MPI |
55 |
|
use mpiSimulation |
56 |
|
#endif |
58 |
|
|
59 |
|
implicit none |
60 |
|
PRIVATE |
61 |
< |
|
61 |
> |
|
62 |
|
integer, parameter :: DP = selected_real_kind(15) |
63 |
|
|
64 |
< |
#define __FORTRAN90 |
65 |
< |
#include "UseTheForce/fForceField.h" |
64 |
> |
logical, save :: useGeometricDistanceMixing = .false. |
65 |
> |
logical, save :: haveMixingMap = .false. |
66 |
|
|
67 |
< |
integer, save :: LJ_Mixing_Policy |
68 |
< |
real(kind=DP), save :: LJ_rcut |
69 |
< |
logical, save :: havePolicy = .false. |
70 |
< |
logical, save :: haveCut = .false. |
71 |
< |
logical, save :: LJ_do_shift = .false. |
72 |
< |
|
73 |
< |
!! Logical has lj force field module been initialized? |
74 |
< |
|
75 |
< |
logical, save :: LJ_FF_initialized = .false. |
76 |
< |
|
77 |
< |
!! Public methods and data |
78 |
< |
public :: init_LJ_FF |
79 |
< |
public :: setCutoffLJ |
80 |
< |
public :: do_lj_pair |
67 |
> |
real(kind=DP), save :: defaultCutoff = 0.0_DP |
68 |
> |
logical, save :: defaultShift = .false. |
69 |
> |
logical, save :: haveDefaultCutoff = .false. |
70 |
> |
|
71 |
> |
|
72 |
> |
type, private :: LJtype |
73 |
> |
integer :: atid |
74 |
> |
real(kind=dp) :: sigma |
75 |
> |
real(kind=dp) :: epsilon |
76 |
> |
logical :: isSoftCore = .false. |
77 |
> |
end type LJtype |
78 |
> |
|
79 |
> |
type, private :: LJList |
80 |
> |
integer :: nLJtypes = 0 |
81 |
> |
integer :: currentLJtype = 0 |
82 |
> |
type(LJtype), pointer :: LJtypes(:) => null() |
83 |
> |
integer, pointer :: atidToLJtype(:) => null() |
84 |
> |
end type LJList |
85 |
> |
|
86 |
> |
type(LJList), save :: LJMap |
87 |
> |
|
88 |
> |
type :: MixParameters |
89 |
> |
real(kind=DP) :: sigma |
90 |
> |
real(kind=DP) :: epsilon |
91 |
> |
real(kind=dp) :: sigma6 |
92 |
> |
real(kind=dp) :: rCut |
93 |
> |
real(kind=dp) :: delta |
94 |
> |
logical :: rCutWasSet = .false. |
95 |
> |
logical :: shiftedPot |
96 |
> |
logical :: isSoftCore = .false. |
97 |
> |
end type MixParameters |
98 |
> |
|
99 |
> |
type(MixParameters), dimension(:,:), allocatable :: MixingMap |
100 |
> |
|
101 |
|
public :: newLJtype |
102 |
< |
|
103 |
< |
!! structure for lj type parameters |
104 |
< |
type, private :: ljType |
105 |
< |
integer :: lj_ident |
106 |
< |
real(kind=dp) :: lj_sigma |
107 |
< |
real(kind=dp) :: lj_epsilon |
108 |
< |
end type ljType |
109 |
< |
|
110 |
< |
!! 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 |
< |
|
102 |
> |
public :: setLJDefaultCutoff |
103 |
> |
public :: setLJUniformCutoff |
104 |
> |
public :: setLJCutoffByTypes |
105 |
> |
public :: getSigma |
106 |
> |
public :: getEpsilon |
107 |
> |
public :: useGeometricMixing |
108 |
> |
public :: do_lj_pair |
109 |
> |
public :: destroyLJtypes |
110 |
> |
|
111 |
|
contains |
112 |
|
|
113 |
< |
subroutine newLJtype(ident,lj_sigma,lj_epsilon,status) |
114 |
< |
integer,intent(in) :: ident |
115 |
< |
real(kind=dp),intent(in) :: lj_sigma |
116 |
< |
real(kind=dp),intent(in) :: lj_epsilon |
113 |
> |
subroutine newLJtype(c_ident, sigma, epsilon, isSoftCore, status) |
114 |
> |
integer,intent(in) :: c_ident |
115 |
> |
real(kind=dp),intent(in) :: sigma |
116 |
> |
real(kind=dp),intent(in) :: epsilon |
117 |
> |
integer, intent(in) :: isSoftCore |
118 |
|
integer,intent(out) :: status |
119 |
< |
|
120 |
< |
integer,pointer :: Matchlist(:) => null() |
119 |
> |
integer :: nLJTypes, ntypes, myATID |
120 |
> |
integer, pointer :: MatchList(:) => null() |
121 |
|
integer :: current |
122 |
< |
integer :: nAtypes |
122 |
> |
|
123 |
|
status = 0 |
124 |
< |
|
125 |
< |
!! Assume that atypes has already been set and get the total number of types in atypes |
90 |
< |
|
91 |
< |
|
124 |
> |
! check to see if this is the first time into this routine... |
125 |
> |
if (.not.associated(LJMap%LJtypes)) then |
126 |
|
|
127 |
< |
! check to see if this is the first time into |
128 |
< |
if (.not.associated(ljParameterList%ljParams)) then |
129 |
< |
call getMatchingElementList(atypes, "is_lj", .true., nAtypes, MatchList) |
130 |
< |
ljParameterList%n_lj_types = nAtypes |
131 |
< |
if (nAtypes == 0) then |
132 |
< |
status = -1 |
133 |
< |
return |
134 |
< |
end if |
135 |
< |
allocate(ljParameterList%ljParams(nAtypes)) |
127 |
> |
call getMatchingElementList(atypes, "is_LennardJones", .true., & |
128 |
> |
nLJTypes, MatchList) |
129 |
> |
|
130 |
> |
LJMap%nLJtypes = nLJTypes |
131 |
> |
|
132 |
> |
allocate(LJMap%LJtypes(nLJTypes)) |
133 |
> |
|
134 |
> |
ntypes = getSize(atypes) |
135 |
> |
|
136 |
> |
allocate(LJMap%atidToLJtype(ntypes)) |
137 |
|
end if |
138 |
|
|
139 |
< |
ljParameterList%currentAddition = ljParameterList%currentAddition + 1 |
140 |
< |
current = ljParameterList%currentAddition |
141 |
< |
|
142 |
< |
! set the values for ljParameterList |
143 |
< |
ljParameterList%ljParams(current)%lj_ident = ident |
144 |
< |
ljParameterList%ljParams(current)%lj_epsilon = lj_epsilon |
145 |
< |
ljParameterList%ljParams(current)%lj_sigma = lj_sigma |
146 |
< |
|
147 |
< |
end subroutine newLJtype |
148 |
< |
|
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 |
139 |
> |
LJMap%currentLJtype = LJMap%currentLJtype + 1 |
140 |
> |
current = LJMap%currentLJtype |
141 |
> |
|
142 |
> |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
143 |
> |
LJMap%atidToLJtype(myATID) = current |
144 |
> |
LJMap%LJtypes(current)%atid = myATID |
145 |
> |
LJMap%LJtypes(current)%sigma = sigma |
146 |
> |
LJMap%LJtypes(current)%epsilon = epsilon |
147 |
> |
if (isSoftCore .eq. 1) then |
148 |
> |
LJMap%LJtypes(current)%isSoftCore = .true. |
149 |
|
else |
150 |
< |
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 |
150 |
> |
LJMap%LJtypes(current)%isSoftCore = .false. |
151 |
|
endif |
152 |
+ |
end subroutine newLJtype |
153 |
|
|
154 |
< |
havePolicy = .true. |
154 |
> |
subroutine setLJDefaultCutoff(thisRcut, shiftedPot) |
155 |
> |
real(kind=dp), intent(in) :: thisRcut |
156 |
> |
logical, intent(in) :: shiftedPot |
157 |
> |
defaultCutoff = thisRcut |
158 |
> |
defaultShift = shiftedPot |
159 |
> |
haveDefaultCutoff = .true. |
160 |
> |
end subroutine setLJDefaultCutoff |
161 |
|
|
162 |
< |
if (haveCut) then |
163 |
< |
status = 0 |
164 |
< |
call createMixingList(myStatus) |
165 |
< |
if (myStatus /= 0) then |
166 |
< |
status = -1 |
167 |
< |
return |
168 |
< |
end if |
169 |
< |
|
141 |
< |
LJ_FF_initialized = .true. |
162 |
> |
subroutine setLJUniformCutoff(thisRcut, shiftedPot) |
163 |
> |
real(kind=dp), intent(in) :: thisRcut |
164 |
> |
logical, intent(in) :: shiftedPot |
165 |
> |
integer :: nLJtypes, i, j |
166 |
> |
|
167 |
> |
if (LJMap%currentLJtype == 0) then |
168 |
> |
call handleError("LJ", "No members in LJMap") |
169 |
> |
return |
170 |
|
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 |
171 |
|
|
172 |
< |
#define __FORTRAN90 |
173 |
< |
#include "UseTheForce/fSwitchingFunction.h" |
172 |
> |
nLJtypes = LJMap%nLJtypes |
173 |
> |
if (.not. allocated(MixingMap)) then |
174 |
> |
allocate(MixingMap(nLJtypes, nLJtypes)) |
175 |
> |
endif |
176 |
|
|
177 |
< |
status = 0 |
177 |
> |
do i = 1, nLJtypes |
178 |
> |
do j = 1, nLJtypes |
179 |
> |
MixingMap(i,j)%rCut = thisRcut |
180 |
> |
MixingMap(i,j)%shiftedPot = shiftedPot |
181 |
> |
MixingMap(i,j)%rCutWasSet = .true. |
182 |
> |
enddo |
183 |
> |
enddo |
184 |
> |
call createMixingMap() |
185 |
> |
end subroutine setLJUniformCutoff |
186 |
|
|
187 |
< |
LJ_rcut = rcut |
188 |
< |
LJ_do_shift = do_shift |
189 |
< |
call set_switch(LJ_SWITCH, rcut, rcut) |
190 |
< |
haveCut = .true. |
187 |
> |
subroutine setLJCutoffByTypes(atid1, atid2, thisRcut, shiftedPot) |
188 |
> |
integer, intent(in) :: atid1, atid2 |
189 |
> |
real(kind=dp), intent(in) :: thisRcut |
190 |
> |
logical, intent(in) :: shiftedPot |
191 |
> |
integer :: nLJtypes, ljt1, ljt2 |
192 |
|
|
193 |
< |
if (havePolicy) then |
194 |
< |
status = 0 |
195 |
< |
call createMixingList(myStatus) |
196 |
< |
if (myStatus /= 0) then |
197 |
< |
status = -1 |
198 |
< |
return |
199 |
< |
end if |
200 |
< |
|
201 |
< |
LJ_FF_initialized = .true. |
202 |
< |
end if |
203 |
< |
|
193 |
> |
if (LJMap%currentLJtype == 0) then |
194 |
> |
call handleError("LJ", "No members in LJMap") |
195 |
> |
return |
196 |
> |
end if |
197 |
> |
|
198 |
> |
nLJtypes = LJMap%nLJtypes |
199 |
> |
if (.not. allocated(MixingMap)) then |
200 |
> |
allocate(MixingMap(nLJtypes, nLJtypes)) |
201 |
> |
endif |
202 |
> |
|
203 |
> |
ljt1 = LJMap%atidToLJtype(atid1) |
204 |
> |
ljt2 = LJMap%atidToLJtype(atid2) |
205 |
> |
|
206 |
> |
MixingMap(ljt1,ljt2)%rCut = thisRcut |
207 |
> |
MixingMap(ljt1,ljt2)%shiftedPot = shiftedPot |
208 |
> |
MixingMap(ljt1,ljt2)%rCutWasSet = .true. |
209 |
> |
MixingMap(ljt2,ljt1)%rCut = thisRcut |
210 |
> |
MixingMap(ljt2,ljt1)%shiftedPot = shiftedPot |
211 |
> |
MixingMap(ljt2,ljt1)%rCutWasSet = .true. |
212 |
> |
|
213 |
> |
call createMixingMap() |
214 |
> |
end subroutine setLJCutoffByTypes |
215 |
> |
|
216 |
> |
function getSigma(atid) result (s) |
217 |
> |
integer, intent(in) :: atid |
218 |
> |
integer :: ljt1 |
219 |
> |
real(kind=dp) :: s |
220 |
> |
|
221 |
> |
if (LJMap%currentLJtype == 0) then |
222 |
> |
call handleError("LJ", "No members in LJMap") |
223 |
> |
return |
224 |
> |
end if |
225 |
> |
|
226 |
> |
ljt1 = LJMap%atidToLJtype(atid) |
227 |
> |
s = LJMap%LJtypes(ljt1)%sigma |
228 |
> |
|
229 |
> |
end function getSigma |
230 |
> |
|
231 |
> |
function getEpsilon(atid) result (e) |
232 |
> |
integer, intent(in) :: atid |
233 |
> |
integer :: ljt1 |
234 |
> |
real(kind=dp) :: e |
235 |
> |
|
236 |
> |
if (LJMap%currentLJtype == 0) then |
237 |
> |
call handleError("LJ", "No members in LJMap") |
238 |
> |
return |
239 |
> |
end if |
240 |
> |
|
241 |
> |
ljt1 = LJMap%atidToLJtype(atid) |
242 |
> |
e = LJMap%LJtypes(ljt1)%epsilon |
243 |
> |
|
244 |
> |
end function getEpsilon |
245 |
> |
|
246 |
> |
subroutine useGeometricMixing() |
247 |
> |
useGeometricDistanceMixing = .true. |
248 |
> |
haveMixingMap = .false. |
249 |
|
return |
250 |
< |
end subroutine setCutoffLJ |
251 |
< |
|
252 |
< |
subroutine createMixingList(status) |
253 |
< |
integer :: nAtypes |
254 |
< |
integer :: status |
255 |
< |
integer :: i |
256 |
< |
integer :: j |
257 |
< |
real ( kind = dp ) :: mySigma_i,mySigma_j |
258 |
< |
real ( kind = dp ) :: myEpsilon_i,myEpsilon_j |
259 |
< |
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 |
250 |
> |
end subroutine useGeometricMixing |
251 |
> |
|
252 |
> |
subroutine createMixingMap() |
253 |
> |
integer :: nLJtypes, i, j |
254 |
> |
real ( kind = dp ) :: s1, s2, e1, e2 |
255 |
> |
real ( kind = dp ) :: rcut6, tp6, tp12 |
256 |
> |
logical :: isSoftCore1, isSoftCore2 |
257 |
> |
|
258 |
> |
if (LJMap%currentLJtype == 0) then |
259 |
> |
call handleError("LJ", "No members in LJMap") |
260 |
|
return |
261 |
|
end if |
262 |
< |
|
263 |
< |
if (.not. associated(ljMixed)) then |
264 |
< |
allocate(ljMixed(nAtypes, nAtypes)) |
262 |
> |
|
263 |
> |
nLJtypes = LJMap%nLJtypes |
264 |
> |
|
265 |
> |
if (.not. allocated(MixingMap)) then |
266 |
> |
allocate(MixingMap(nLJtypes, nLJtypes)) |
267 |
|
endif |
268 |
|
|
269 |
< |
rcut6 = LJ_rcut**6 |
269 |
> |
do i = 1, nLJtypes |
270 |
|
|
271 |
< |
! This loops through all atypes, even those that don't support LJ forces. |
272 |
< |
do i = 1, nAtypes |
271 |
> |
s1 = LJMap%LJtypes(i)%sigma |
272 |
> |
e1 = LJMap%LJtypes(i)%epsilon |
273 |
> |
isSoftCore1 = LJMap%LJtypes(i)%isSoftCore |
274 |
|
|
275 |
< |
myEpsilon_i = ljParameterList%ljParams(i)%lj_epsilon |
203 |
< |
mySigma_i = ljParameterList%ljParams(i)%lj_sigma |
275 |
> |
do j = i, nLJtypes |
276 |
|
|
277 |
< |
! do self mixing rule |
278 |
< |
ljMixed(i,i)%sigma = mySigma_i |
277 |
> |
s2 = LJMap%LJtypes(j)%sigma |
278 |
> |
e2 = LJMap%LJtypes(j)%epsilon |
279 |
> |
isSoftCore2 = LJMap%LJtypes(j)%isSoftCore |
280 |
|
|
281 |
< |
ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6 |
281 |
> |
MixingMap(i,j)%isSoftCore = isSoftCore1 .or. isSoftCore2 |
282 |
> |
|
283 |
> |
! only the distance parameter uses different mixing policies |
284 |
> |
if (useGeometricDistanceMixing) then |
285 |
> |
MixingMap(i,j)%sigma = dsqrt(s1 * s2) |
286 |
> |
else |
287 |
> |
MixingMap(i,j)%sigma = 0.5_dp * (s1 + s2) |
288 |
> |
endif |
289 |
|
|
290 |
< |
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 |
290 |
> |
MixingMap(i,j)%epsilon = dsqrt(e1 * e2) |
291 |
|
|
292 |
< |
myEpsilon_j = ljParameterList%ljParams(j)%lj_epsilon |
223 |
< |
mySigma_j = ljParameterList%ljParams(j)%lj_sigma |
292 |
> |
MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6 |
293 |
|
|
294 |
< |
|
295 |
< |
ljMixed(i,j)%sigma = & |
296 |
< |
calcLJMix("sigma",mySigma_i, & |
297 |
< |
mySigma_j) |
298 |
< |
|
299 |
< |
ljMixed(i,j)%sigma6 = & |
300 |
< |
(ljMixed(i,j)%sigma)**6 |
301 |
< |
|
302 |
< |
|
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 |
294 |
> |
if (MixingMap(i,j)%rCutWasSet) then |
295 |
> |
rcut6 = (MixingMap(i,j)%rcut)**6 |
296 |
> |
else |
297 |
> |
if (haveDefaultCutoff) then |
298 |
> |
rcut6 = defaultCutoff**6 |
299 |
> |
else |
300 |
> |
call handleError("LJ", "No specified or default cutoff value!") |
301 |
> |
endif |
302 |
> |
endif |
303 |
|
|
304 |
< |
end do |
305 |
< |
end do |
304 |
> |
tp6 = MixingMap(i,j)%sigma6/rcut6 |
305 |
> |
tp12 = tp6**2 |
306 |
> |
|
307 |
> |
MixingMap(i,j)%delta =-4.0_DP*MixingMap(i,j)%epsilon*(tp12 - tp6) |
308 |
> |
enddo |
309 |
> |
enddo |
310 |
|
|
311 |
< |
end subroutine createMixingList |
311 |
> |
haveMixingMap = .true. |
312 |
> |
|
313 |
> |
end subroutine createMixingMap |
314 |
|
|
315 |
|
subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
316 |
|
pot, f, do_pot) |
317 |
|
|
318 |
|
integer, intent(in) :: atom1, atom2 |
319 |
+ |
integer :: atid1, atid2, ljt1, ljt2 |
320 |
|
real( kind = dp ), intent(in) :: rij, r2 |
321 |
|
real( kind = dp ) :: pot, sw, vpair |
322 |
|
real( kind = dp ), dimension(3,nLocal) :: f |
334 |
|
real( kind = dp ) :: t6 |
335 |
|
real( kind = dp ) :: t12 |
336 |
|
real( kind = dp ) :: delta |
337 |
< |
integer :: id1, id2 |
337 |
> |
logical :: isSoftCore, shiftedPot |
338 |
> |
integer :: id1, id2, localError |
339 |
|
|
340 |
+ |
if (.not.haveMixingMap) then |
341 |
+ |
call createMixingMap() |
342 |
+ |
endif |
343 |
+ |
|
344 |
|
! Look up the correct parameters in the mixing matrix |
345 |
|
#ifdef IS_MPI |
346 |
< |
sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6 |
347 |
< |
epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon |
286 |
< |
delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta |
346 |
> |
atid1 = atid_Row(atom1) |
347 |
> |
atid2 = atid_Col(atom2) |
348 |
|
#else |
349 |
< |
sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6 |
350 |
< |
epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon |
290 |
< |
delta = ljMixed(atid(atom1),atid(atom2))%delta |
349 |
> |
atid1 = atid(atom1) |
350 |
> |
atid2 = atid(atom2) |
351 |
|
#endif |
352 |
|
|
353 |
+ |
ljt1 = LJMap%atidToLJtype(atid1) |
354 |
+ |
ljt2 = LJMap%atidToLJtype(atid2) |
355 |
+ |
|
356 |
+ |
sigma6 = MixingMap(ljt1,ljt2)%sigma6 |
357 |
+ |
epsilon = MixingMap(ljt1,ljt2)%epsilon |
358 |
+ |
delta = MixingMap(ljt1,ljt2)%delta |
359 |
+ |
isSoftCore = MixingMap(ljt1,ljt2)%isSoftCore |
360 |
+ |
shiftedPot = MixingMap(ljt1,ljt2)%shiftedPot |
361 |
+ |
|
362 |
|
r6 = r2 * r2 * r2 |
363 |
< |
|
363 |
> |
|
364 |
|
t6 = sigma6/ r6 |
365 |
|
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 |
366 |
|
|
367 |
< |
vpair = vpair + pot_temp |
367 |
> |
if (isSoftCore) then |
368 |
|
|
369 |
< |
dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij |
369 |
> |
pot_temp = 4.0E0_DP * epsilon * t6 |
370 |
> |
if (shiftedPot) then |
371 |
> |
pot_temp = pot_temp + delta |
372 |
> |
endif |
373 |
|
|
374 |
+ |
vpair = vpair + pot_temp |
375 |
+ |
|
376 |
+ |
dudr = -sw * 24.0E0_DP * epsilon * t6 / rij |
377 |
+ |
|
378 |
+ |
else |
379 |
+ |
pot_temp = 4.0E0_DP * epsilon * (t12 - t6) |
380 |
+ |
if (shiftedPot) then |
381 |
+ |
pot_temp = pot_temp + delta |
382 |
+ |
endif |
383 |
+ |
|
384 |
+ |
vpair = vpair + pot_temp |
385 |
+ |
|
386 |
+ |
dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij |
387 |
+ |
endif |
388 |
+ |
|
389 |
|
drdx = d(1) / rij |
390 |
|
drdy = d(2) / rij |
391 |
|
drdz = d(3) / rij |
392 |
< |
|
392 |
> |
|
393 |
|
fx = dudr * drdx |
394 |
|
fy = dudr * drdy |
395 |
|
fz = dudr * drdz |
396 |
< |
|
397 |
< |
|
396 |
> |
|
397 |
> |
|
398 |
|
#ifdef IS_MPI |
399 |
|
if (do_pot) then |
400 |
|
pot_Row(atom1) = pot_Row(atom1) + sw*pot_temp*0.5 |
401 |
|
pot_Col(atom2) = pot_Col(atom2) + sw*pot_temp*0.5 |
402 |
|
endif |
403 |
< |
|
403 |
> |
|
404 |
|
f_Row(1,atom1) = f_Row(1,atom1) + fx |
405 |
|
f_Row(2,atom1) = f_Row(2,atom1) + fy |
406 |
|
f_Row(3,atom1) = f_Row(3,atom1) + fz |
407 |
< |
|
407 |
> |
|
408 |
|
f_Col(1,atom2) = f_Col(1,atom2) - fx |
409 |
|
f_Col(2,atom2) = f_Col(2,atom2) - fy |
410 |
|
f_Col(3,atom2) = f_Col(3,atom2) - fz |
411 |
< |
|
411 |
> |
|
412 |
|
#else |
413 |
|
if (do_pot) pot = pot + sw*pot_temp |
414 |
|
|
415 |
|
f(1,atom1) = f(1,atom1) + fx |
416 |
|
f(2,atom1) = f(2,atom1) + fy |
417 |
|
f(3,atom1) = f(3,atom1) + fz |
418 |
< |
|
418 |
> |
|
419 |
|
f(1,atom2) = f(1,atom2) - fx |
420 |
|
f(2,atom2) = f(2,atom2) - fy |
421 |
|
f(3,atom2) = f(3,atom2) - fz |
422 |
|
#endif |
423 |
< |
|
423 |
> |
|
424 |
|
#ifdef IS_MPI |
425 |
|
id1 = AtomRowToGlobal(atom1) |
426 |
|
id2 = AtomColToGlobal(atom2) |
430 |
|
#endif |
431 |
|
|
432 |
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
433 |
< |
|
433 |
> |
|
434 |
|
fpair(1) = fpair(1) + fx |
435 |
|
fpair(2) = fpair(2) + fy |
436 |
|
fpair(3) = fpair(3) + fz |
438 |
|
endif |
439 |
|
|
440 |
|
return |
441 |
< |
|
441 |
> |
|
442 |
|
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 |
443 |
|
|
444 |
< |
integer, optional :: status |
444 |
> |
subroutine destroyLJTypes() |
445 |
|
|
446 |
< |
myMixParam = 0.0_dp |
446 |
> |
LJMap%nLJtypes = 0 |
447 |
> |
LJMap%currentLJtype = 0 |
448 |
|
|
449 |
< |
if (present(status)) status = 0 |
450 |
< |
select case (LJ_Mixing_Policy) |
451 |
< |
case (1) |
452 |
< |
select case (thisParam) |
453 |
< |
case ("sigma") |
454 |
< |
myMixParam = 0.5_dp * (param1 + param2) |
455 |
< |
case ("epsilon") |
456 |
< |
myMixParam = sqrt(param1 * param2) |
457 |
< |
case default |
458 |
< |
status = -1 |
459 |
< |
end select |
460 |
< |
case default |
387 |
< |
status = -1 |
388 |
< |
end select |
389 |
< |
end function calcLJMix |
390 |
< |
|
391 |
< |
end module lj |
449 |
> |
if (associated(LJMap%LJtypes)) then |
450 |
> |
deallocate(LJMap%LJtypes) |
451 |
> |
LJMap%LJtypes => null() |
452 |
> |
end if |
453 |
> |
|
454 |
> |
if (associated(LJMap%atidToLJtype)) then |
455 |
> |
deallocate(LJMap%atidToLJtype) |
456 |
> |
LJMap%atidToLJtype => null() |
457 |
> |
end if |
458 |
> |
|
459 |
> |
haveMixingMap = .false. |
460 |
> |
end subroutine destroyLJTypes |
461 |
|
|
462 |
< |
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 |
< |
|
462 |
> |
end module lj |