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. |
44 |
!! @author Charles F. Vardeman II |
45 |
!! @author Matthew Meineke |
46 |
!! @version $Id: LJ.F90,v 1.27 2007-04-09 18:24:00 gezelter Exp $, $Date: 2007-04-09 18:24:00 $, $Name: not supported by cvs2svn $, $Revision: 1.27 $ |
47 |
|
48 |
|
49 |
module lj |
50 |
use definitions |
51 |
use atype_module |
52 |
use vector_class |
53 |
use simulation |
54 |
use status |
55 |
use fForceOptions |
56 |
#ifdef IS_MPI |
57 |
use mpiSimulation |
58 |
#endif |
59 |
use force_globals |
60 |
|
61 |
implicit none |
62 |
PRIVATE |
63 |
#define __FORTRAN90 |
64 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
65 |
|
66 |
logical, save :: useGeometricDistanceMixing = .false. |
67 |
logical, save :: haveMixingMap = .false. |
68 |
|
69 |
real(kind=DP), save :: defaultCutoff = 0.0_DP |
70 |
logical, save :: defaultShift = .false. |
71 |
logical, save :: haveDefaultCutoff = .false. |
72 |
|
73 |
type, private :: LJtype |
74 |
integer :: atid |
75 |
real(kind=dp) :: sigma |
76 |
real(kind=dp) :: epsilon |
77 |
logical :: isSoftCore = .false. |
78 |
end type LJtype |
79 |
|
80 |
type, private :: LJList |
81 |
integer :: Nljtypes = 0 |
82 |
integer :: currentLJtype = 0 |
83 |
type(LJtype), pointer :: LJtypes(:) => null() |
84 |
integer, pointer :: atidToLJtype(:) => null() |
85 |
end type LJList |
86 |
|
87 |
type(LJList), save :: LJMap |
88 |
|
89 |
type :: MixParameters |
90 |
real(kind=DP) :: sigma |
91 |
real(kind=DP) :: epsilon |
92 |
real(kind=dp) :: sigmai |
93 |
real(kind=dp) :: rCut |
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 |
public :: setLJDefaultCutoff |
103 |
public :: getSigma |
104 |
public :: getEpsilon |
105 |
public :: do_lj_pair |
106 |
public :: destroyLJtypes |
107 |
|
108 |
contains |
109 |
|
110 |
subroutine newLJtype(c_ident, sigma, epsilon, isSoftCore, status) |
111 |
integer,intent(in) :: c_ident |
112 |
real(kind=dp),intent(in) :: sigma |
113 |
real(kind=dp),intent(in) :: epsilon |
114 |
integer, intent(in) :: isSoftCore |
115 |
integer,intent(out) :: status |
116 |
integer :: nLJTypes, ntypes, myATID |
117 |
integer, pointer :: MatchList(:) => null() |
118 |
integer :: current |
119 |
|
120 |
status = 0 |
121 |
! check to see if this is the first time into this routine... |
122 |
if (.not.associated(LJMap%LJtypes)) then |
123 |
|
124 |
call getMatchingElementList(atypes, "is_LennardJones", .true., & |
125 |
nLJTypes, MatchList) |
126 |
|
127 |
LJMap%nLJtypes = nLJTypes |
128 |
|
129 |
allocate(LJMap%LJtypes(nLJTypes)) |
130 |
|
131 |
ntypes = getSize(atypes) |
132 |
|
133 |
allocate(LJMap%atidToLJtype(ntypes)) |
134 |
end if |
135 |
|
136 |
LJMap%currentLJtype = LJMap%currentLJtype + 1 |
137 |
current = LJMap%currentLJtype |
138 |
|
139 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
140 |
LJMap%atidToLJtype(myATID) = current |
141 |
LJMap%LJtypes(current)%atid = myATID |
142 |
LJMap%LJtypes(current)%sigma = sigma |
143 |
LJMap%LJtypes(current)%epsilon = epsilon |
144 |
if (isSoftCore .eq. 1) then |
145 |
LJMap%LJtypes(current)%isSoftCore = .true. |
146 |
else |
147 |
LJMap%LJtypes(current)%isSoftCore = .false. |
148 |
endif |
149 |
end subroutine newLJtype |
150 |
|
151 |
subroutine setLJDefaultCutoff(thisRcut, shiftedPot) |
152 |
real(kind=dp), intent(in) :: thisRcut |
153 |
logical, intent(in) :: shiftedPot |
154 |
defaultCutoff = thisRcut |
155 |
defaultShift = shiftedPot |
156 |
haveDefaultCutoff = .true. |
157 |
!we only want to build LJ Mixing map if LJ is being used. |
158 |
if(LJMap%nLJTypes /= 0) then |
159 |
call createMixingMap() |
160 |
end if |
161 |
|
162 |
end subroutine setLJDefaultCutoff |
163 |
|
164 |
function getSigma(atid) result (s) |
165 |
integer, intent(in) :: atid |
166 |
integer :: ljt1 |
167 |
real(kind=dp) :: s |
168 |
|
169 |
if (LJMap%currentLJtype == 0) then |
170 |
call handleError("LJ", "No members in LJMap") |
171 |
return |
172 |
end if |
173 |
|
174 |
ljt1 = LJMap%atidToLJtype(atid) |
175 |
s = LJMap%LJtypes(ljt1)%sigma |
176 |
|
177 |
end function getSigma |
178 |
|
179 |
function getEpsilon(atid) result (e) |
180 |
integer, intent(in) :: atid |
181 |
integer :: ljt1 |
182 |
real(kind=dp) :: e |
183 |
|
184 |
if (LJMap%currentLJtype == 0) then |
185 |
call handleError("LJ", "No members in LJMap") |
186 |
return |
187 |
end if |
188 |
|
189 |
ljt1 = LJMap%atidToLJtype(atid) |
190 |
e = LJMap%LJtypes(ljt1)%epsilon |
191 |
|
192 |
end function getEpsilon |
193 |
|
194 |
subroutine createMixingMap() |
195 |
integer :: nLJtypes, i, j |
196 |
real ( kind = dp ) :: s1, s2, e1, e2 |
197 |
real ( kind = dp ) :: rcut6, tp6, tp12 |
198 |
logical :: isSoftCore1, isSoftCore2, doShift |
199 |
|
200 |
if (LJMap%currentLJtype == 0) then |
201 |
call handleError("LJ", "No members in LJMap") |
202 |
return |
203 |
end if |
204 |
|
205 |
nLJtypes = LJMap%nLJtypes |
206 |
|
207 |
if (.not. allocated(MixingMap)) then |
208 |
allocate(MixingMap(nLJtypes, nLJtypes)) |
209 |
endif |
210 |
|
211 |
useGeometricDistanceMixing = usesGeometricDistanceMixing() |
212 |
do i = 1, nLJtypes |
213 |
|
214 |
s1 = LJMap%LJtypes(i)%sigma |
215 |
e1 = LJMap%LJtypes(i)%epsilon |
216 |
isSoftCore1 = LJMap%LJtypes(i)%isSoftCore |
217 |
|
218 |
do j = i, nLJtypes |
219 |
|
220 |
s2 = LJMap%LJtypes(j)%sigma |
221 |
e2 = LJMap%LJtypes(j)%epsilon |
222 |
isSoftCore2 = LJMap%LJtypes(j)%isSoftCore |
223 |
|
224 |
MixingMap(i,j)%isSoftCore = isSoftCore1 .or. isSoftCore2 |
225 |
|
226 |
! only the distance parameter uses different mixing policies |
227 |
if (useGeometricDistanceMixing) then |
228 |
MixingMap(i,j)%sigma = sqrt(s1 * s2) |
229 |
else |
230 |
MixingMap(i,j)%sigma = 0.5_dp * (s1 + s2) |
231 |
endif |
232 |
|
233 |
MixingMap(i,j)%epsilon = sqrt(e1 * e2) |
234 |
|
235 |
MixingMap(i,j)%sigmai = 1.0_DP / (MixingMap(i,j)%sigma) |
236 |
|
237 |
if (haveDefaultCutoff) then |
238 |
MixingMap(i,j)%shiftedPot = defaultShift |
239 |
else |
240 |
MixingMap(i,j)%shiftedPot = defaultShift |
241 |
endif |
242 |
|
243 |
if (i.ne.j) then |
244 |
MixingMap(j,i)%sigma = MixingMap(i,j)%sigma |
245 |
MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon |
246 |
MixingMap(j,i)%sigmai = MixingMap(i,j)%sigmai |
247 |
MixingMap(j,i)%rCut = MixingMap(i,j)%rCut |
248 |
MixingMap(j,i)%rCutWasSet = MixingMap(i,j)%rCutWasSet |
249 |
MixingMap(j,i)%shiftedPot = MixingMap(i,j)%shiftedPot |
250 |
MixingMap(j,i)%isSoftCore = MixingMap(i,j)%isSoftCore |
251 |
endif |
252 |
|
253 |
enddo |
254 |
enddo |
255 |
|
256 |
haveMixingMap = .true. |
257 |
|
258 |
end subroutine createMixingMap |
259 |
|
260 |
subroutine do_lj_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, & |
261 |
pot, f, do_pot) |
262 |
|
263 |
integer, intent(in) :: atom1, atom2 |
264 |
integer :: atid1, atid2, ljt1, ljt2 |
265 |
real( kind = dp ), intent(in) :: rij, r2, rcut |
266 |
real( kind = dp ) :: pot, sw, vpair |
267 |
real( kind = dp ), dimension(3,nLocal) :: f |
268 |
real( kind = dp ), intent(in), dimension(3) :: d |
269 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
270 |
logical, intent(in) :: do_pot |
271 |
|
272 |
! local Variables |
273 |
real( kind = dp ) :: drdx, drdy, drdz |
274 |
real( kind = dp ) :: fx, fy, fz |
275 |
real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos |
276 |
real( kind = dp ) :: pot_temp, dudr |
277 |
real( kind = dp ) :: sigmai |
278 |
real( kind = dp ) :: epsilon |
279 |
logical :: isSoftCore, shiftedPot |
280 |
integer :: id1, id2, localError |
281 |
|
282 |
if (.not.haveMixingMap) then |
283 |
call createMixingMap() |
284 |
endif |
285 |
|
286 |
! Look up the correct parameters in the mixing matrix |
287 |
#ifdef IS_MPI |
288 |
atid1 = atid_Row(atom1) |
289 |
atid2 = atid_Col(atom2) |
290 |
#else |
291 |
atid1 = atid(atom1) |
292 |
atid2 = atid(atom2) |
293 |
#endif |
294 |
|
295 |
ljt1 = LJMap%atidToLJtype(atid1) |
296 |
ljt2 = LJMap%atidToLJtype(atid2) |
297 |
|
298 |
sigmai = MixingMap(ljt1,ljt2)%sigmai |
299 |
epsilon = MixingMap(ljt1,ljt2)%epsilon |
300 |
isSoftCore = MixingMap(ljt1,ljt2)%isSoftCore |
301 |
shiftedPot = MixingMap(ljt1,ljt2)%shiftedPot |
302 |
|
303 |
ros = rij * sigmai |
304 |
myPotC = 0.0_DP |
305 |
myDerivC = 0.0_DP |
306 |
|
307 |
if (isSoftCore) then |
308 |
|
309 |
call getSoftFunc(ros, myPot, myDeriv) |
310 |
|
311 |
if (shiftedPot) then |
312 |
rcos = rcut * sigmai |
313 |
call getSoftFunc(rcos, myPotC, myDerivC) |
314 |
endif |
315 |
|
316 |
else |
317 |
|
318 |
call getLJfunc(ros, myPot, myDeriv) |
319 |
|
320 |
if (shiftedPot) then |
321 |
rcos = rcut * sigmai |
322 |
call getLJfunc(rcos, myPotC, myDerivC) |
323 |
endif |
324 |
|
325 |
endif |
326 |
|
327 |
!! these are the shifted POTENTIAL variants. |
328 |
! pot_temp = epsilon * (myPot - myPotC) |
329 |
! dudr = sw * epsilon * myDeriv * sigmai |
330 |
|
331 |
!! these are the shifted FORCE variants. |
332 |
|
333 |
pot_temp = epsilon * (myPot - myPotC - myDerivC * (rij - rcut) * sigmai) |
334 |
dudr = sw * epsilon * (myDeriv - myDerivC) * sigmai |
335 |
|
336 |
vpair = vpair + pot_temp |
337 |
|
338 |
drdx = d(1) / rij |
339 |
drdy = d(2) / rij |
340 |
drdz = d(3) / rij |
341 |
|
342 |
fx = dudr * drdx |
343 |
fy = dudr * drdy |
344 |
fz = dudr * drdz |
345 |
|
346 |
#ifdef IS_MPI |
347 |
if (do_pot) then |
348 |
pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5 |
349 |
pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5 |
350 |
endif |
351 |
|
352 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
353 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
354 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
355 |
|
356 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
357 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
358 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
359 |
|
360 |
#else |
361 |
if (do_pot) pot = pot + sw*pot_temp |
362 |
|
363 |
f(1,atom1) = f(1,atom1) + fx |
364 |
f(2,atom1) = f(2,atom1) + fy |
365 |
f(3,atom1) = f(3,atom1) + fz |
366 |
|
367 |
f(1,atom2) = f(1,atom2) - fx |
368 |
f(2,atom2) = f(2,atom2) - fy |
369 |
f(3,atom2) = f(3,atom2) - fz |
370 |
#endif |
371 |
|
372 |
#ifdef IS_MPI |
373 |
id1 = AtomRowToGlobal(atom1) |
374 |
id2 = AtomColToGlobal(atom2) |
375 |
#else |
376 |
id1 = atom1 |
377 |
id2 = atom2 |
378 |
#endif |
379 |
|
380 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
381 |
|
382 |
fpair(1) = fpair(1) + fx |
383 |
fpair(2) = fpair(2) + fy |
384 |
fpair(3) = fpair(3) + fz |
385 |
|
386 |
endif |
387 |
|
388 |
return |
389 |
|
390 |
end subroutine do_lj_pair |
391 |
|
392 |
subroutine destroyLJTypes() |
393 |
|
394 |
LJMap%nLJtypes = 0 |
395 |
LJMap%currentLJtype = 0 |
396 |
|
397 |
if (associated(LJMap%LJtypes)) then |
398 |
deallocate(LJMap%LJtypes) |
399 |
LJMap%LJtypes => null() |
400 |
end if |
401 |
|
402 |
if (associated(LJMap%atidToLJtype)) then |
403 |
deallocate(LJMap%atidToLJtype) |
404 |
LJMap%atidToLJtype => null() |
405 |
end if |
406 |
|
407 |
haveMixingMap = .false. |
408 |
|
409 |
end subroutine destroyLJTypes |
410 |
|
411 |
subroutine getLJfunc(r, myPot, myDeriv) |
412 |
|
413 |
real(kind=dp), intent(in) :: r |
414 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
415 |
real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13 |
416 |
real(kind=dp) :: a, b, c, d, dx |
417 |
integer :: j |
418 |
|
419 |
ri = 1.0_DP / r |
420 |
ri2 = ri*ri |
421 |
ri6 = ri2*ri2*ri2 |
422 |
ri7 = ri6*ri |
423 |
ri12 = ri6*ri6 |
424 |
ri13 = ri12*ri |
425 |
|
426 |
myPot = 4.0_DP * (ri12 - ri6) |
427 |
myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13) |
428 |
|
429 |
return |
430 |
end subroutine getLJfunc |
431 |
|
432 |
subroutine getSoftFunc(r, myPot, myDeriv) |
433 |
|
434 |
real(kind=dp), intent(in) :: r |
435 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
436 |
real(kind=dp) :: ri, ri2, ri6, ri7 |
437 |
real(kind=dp) :: a, b, c, d, dx |
438 |
integer :: j |
439 |
|
440 |
ri = 1.0_DP / r |
441 |
ri2 = ri*ri |
442 |
ri6 = ri2*ri2*ri2 |
443 |
ri7 = ri6*ri |
444 |
myPot = 4.0_DP * (ri6) |
445 |
myDeriv = - 24.0_DP * ri7 |
446 |
|
447 |
return |
448 |
end subroutine getSoftFunc |
449 |
|
450 |
end module lj |