ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90
Revision: 3129
Committed: Fri Apr 20 18:15:48 2007 UTC (17 years, 3 months ago) by chrisfen
File size: 13581 byte(s)
Log Message:
SF Lennard-Jones was added for everyones' enjoyment.  The behavior is tethered to the electrostaticSummationMethod keyword.

File Contents

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