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 (18 years, 4 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

# User Rev Content
1 gezelter 1930 !!
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 gezelter 1608 !! Calculates Long Range forces Lennard-Jones interactions.
44     !! @author Charles F. Vardeman II
45     !! @author Matthew Meineke
46 chrisfen 3129 !! @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 gezelter 1608
48 gezelter 1930
49 gezelter 1608 module lj
50 gezelter 2756 use definitions
51 gezelter 1608 use atype_module
52     use vector_class
53     use simulation
54 gezelter 1628 use status
55 chuckv 2497 use fForceOptions
56 gezelter 1608 #ifdef IS_MPI
57     use mpiSimulation
58     #endif
59     use force_globals
60    
61     implicit none
62     PRIVATE
63 chuckv 2355 #define __FORTRAN90
64     #include "UseTheForce/DarkSide/fInteractionMap.h"
65 gezelter 2204
66 gezelter 2271 logical, save :: useGeometricDistanceMixing = .false.
67     logical, save :: haveMixingMap = .false.
68    
69     real(kind=DP), save :: defaultCutoff = 0.0_DP
70 chrisfen 3129 logical, save :: defaultShiftPot = .false.
71     logical, save :: defaultShiftFrc = .false.
72 gezelter 2271 logical, save :: haveDefaultCutoff = .false.
73    
74     type, private :: LJtype
75     integer :: atid
76 gezelter 1628 real(kind=dp) :: sigma
77     real(kind=dp) :: epsilon
78 gezelter 2271 logical :: isSoftCore = .false.
79     end type LJtype
80 gezelter 2204
81 gezelter 2271 type, private :: LJList
82 chuckv 2319 integer :: Nljtypes = 0
83 gezelter 2271 integer :: currentLJtype = 0
84     type(LJtype), pointer :: LJtypes(:) => null()
85     integer, pointer :: atidToLJtype(:) => null()
86     end type LJList
87 gezelter 2204
88 gezelter 2271 type(LJList), save :: LJMap
89 gezelter 2204
90 gezelter 1628 type :: MixParameters
91     real(kind=DP) :: sigma
92     real(kind=DP) :: epsilon
93 gezelter 2717 real(kind=dp) :: sigmai
94 gezelter 2271 real(kind=dp) :: rCut
95     logical :: rCutWasSet = .false.
96     logical :: shiftedPot
97     logical :: isSoftCore = .false.
98 chrisfen 3129 logical :: shiftedFrc
99 gezelter 1628 end type MixParameters
100 gezelter 2204
101 gezelter 1628 type(MixParameters), dimension(:,:), allocatable :: MixingMap
102 gezelter 2204
103 gezelter 2271 public :: newLJtype
104     public :: setLJDefaultCutoff
105     public :: getSigma
106     public :: getEpsilon
107 gezelter 1628 public :: do_lj_pair
108 gezelter 2271 public :: destroyLJtypes
109 gezelter 2204
110 gezelter 1608 contains
111    
112 gezelter 2271 subroutine newLJtype(c_ident, sigma, epsilon, isSoftCore, status)
113 gezelter 1930 integer,intent(in) :: c_ident
114 gezelter 1628 real(kind=dp),intent(in) :: sigma
115     real(kind=dp),intent(in) :: epsilon
116 gezelter 2271 integer, intent(in) :: isSoftCore
117 chuckv 1624 integer,intent(out) :: status
118 gezelter 2271 integer :: nLJTypes, ntypes, myATID
119 gezelter 1930 integer, pointer :: MatchList(:) => null()
120 gezelter 2271 integer :: current
121 gezelter 1628
122 chuckv 1624 status = 0
123 gezelter 2271 ! check to see if this is the first time into this routine...
124     if (.not.associated(LJMap%LJtypes)) then
125 gezelter 2204
126 gezelter 2271 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 gezelter 1930 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
142 gezelter 2271 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 chuckv 1624
153 chrisfen 3129 subroutine setLJDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
154 gezelter 2271 real(kind=dp), intent(in) :: thisRcut
155     logical, intent(in) :: shiftedPot
156 chrisfen 3129 logical, intent(in) :: shiftedFrc
157 gezelter 2271 defaultCutoff = thisRcut
158 chrisfen 3129 defaultShiftPot = shiftedPot
159     defaultShiftFrc = shiftedFrc
160 gezelter 2271 haveDefaultCutoff = .true.
161 chrisfen 3129
162 chrisfen 2727 !we only want to build LJ Mixing map if LJ is being used.
163 chuckv 2319 if(LJMap%nLJTypes /= 0) then
164     call createMixingMap()
165     end if
166 gezelter 2717
167 gezelter 2271 end subroutine setLJDefaultCutoff
168 gezelter 2204
169 gezelter 1633 function getSigma(atid) result (s)
170     integer, intent(in) :: atid
171 gezelter 2271 integer :: ljt1
172 gezelter 1633 real(kind=dp) :: s
173 gezelter 2204
174 gezelter 2271 if (LJMap%currentLJtype == 0) then
175     call handleError("LJ", "No members in LJMap")
176 gezelter 1633 return
177     end if
178 gezelter 2204
179 gezelter 2271 ljt1 = LJMap%atidToLJtype(atid)
180     s = LJMap%LJtypes(ljt1)%sigma
181    
182 gezelter 1633 end function getSigma
183    
184     function getEpsilon(atid) result (e)
185     integer, intent(in) :: atid
186 gezelter 2271 integer :: ljt1
187 gezelter 1633 real(kind=dp) :: e
188 gezelter 2204
189 gezelter 2271 if (LJMap%currentLJtype == 0) then
190     call handleError("LJ", "No members in LJMap")
191 gezelter 1633 return
192     end if
193 gezelter 2204
194 gezelter 2271 ljt1 = LJMap%atidToLJtype(atid)
195     e = LJMap%LJtypes(ljt1)%epsilon
196    
197 gezelter 1633 end function getEpsilon
198    
199 gezelter 2271 subroutine createMixingMap()
200     integer :: nLJtypes, i, j
201     real ( kind = dp ) :: s1, s2, e1, e2
202     real ( kind = dp ) :: rcut6, tp6, tp12
203 gezelter 2289 logical :: isSoftCore1, isSoftCore2, doShift
204 gezelter 1628
205 gezelter 2271 if (LJMap%currentLJtype == 0) then
206     call handleError("LJ", "No members in LJMap")
207 gezelter 2086 return
208 gezelter 2271 end if
209 gezelter 2204
210 gezelter 2271 nLJtypes = LJMap%nLJtypes
211 gezelter 2204
212 gezelter 2086 if (.not. allocated(MixingMap)) then
213 gezelter 2271 allocate(MixingMap(nLJtypes, nLJtypes))
214 gezelter 2086 endif
215    
216 chuckv 2497 useGeometricDistanceMixing = usesGeometricDistanceMixing()
217 gezelter 2271 do i = 1, nLJtypes
218 gezelter 2204
219 gezelter 2271 s1 = LJMap%LJtypes(i)%sigma
220     e1 = LJMap%LJtypes(i)%epsilon
221     isSoftCore1 = LJMap%LJtypes(i)%isSoftCore
222 gezelter 2204
223 gezelter 2271 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 gezelter 2204
231 gezelter 2271 ! only the distance parameter uses different mixing policies
232     if (useGeometricDistanceMixing) then
233 gezelter 2756 MixingMap(i,j)%sigma = sqrt(s1 * s2)
234 gezelter 2271 else
235     MixingMap(i,j)%sigma = 0.5_dp * (s1 + s2)
236     endif
237    
238 gezelter 2756 MixingMap(i,j)%epsilon = sqrt(e1 * e2)
239 gezelter 2204
240 gezelter 2717 MixingMap(i,j)%sigmai = 1.0_DP / (MixingMap(i,j)%sigma)
241 gezelter 2204
242 gezelter 2461 if (haveDefaultCutoff) then
243 chrisfen 3129 MixingMap(i,j)%shiftedPot = defaultShiftPot
244     MixingMap(i,j)%shiftedFrc = defaultShiftFrc
245 gezelter 2271 else
246 chrisfen 3129 MixingMap(i,j)%shiftedPot = defaultShiftPot
247     MixingMap(i,j)%shiftedFrc = defaultShiftFrc
248 gezelter 2461 endif
249 gezelter 2204
250 gezelter 2289 if (i.ne.j) then
251     MixingMap(j,i)%sigma = MixingMap(i,j)%sigma
252     MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
253 gezelter 2717 MixingMap(j,i)%sigmai = MixingMap(i,j)%sigmai
254 gezelter 2289 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 chrisfen 3129 MixingMap(j,i)%shiftedFrc = MixingMap(i,j)%shiftedFrc
259 gezelter 2289 endif
260    
261 gezelter 2271 enddo
262     enddo
263    
264 chrisfen 1754 haveMixingMap = .true.
265 chrisfen 3129
266 gezelter 1628 end subroutine createMixingMap
267 chrisfen 2727
268 gezelter 2461 subroutine do_lj_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
269 gezelter 1608 pot, f, do_pot)
270 gezelter 2461
271 gezelter 1608 integer, intent(in) :: atom1, atom2
272 gezelter 2271 integer :: atid1, atid2, ljt1, ljt2
273 gezelter 2461 real( kind = dp ), intent(in) :: rij, r2, rcut
274 gezelter 1608 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 gezelter 2717 real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos
284 gezelter 1608 real( kind = dp ) :: pot_temp, dudr
285 gezelter 2717 real( kind = dp ) :: sigmai
286 gezelter 1608 real( kind = dp ) :: epsilon
287 chrisfen 3129 logical :: isSoftCore, shiftedPot, shiftedFrc
288 gezelter 1628 integer :: id1, id2, localError
289 gezelter 1608
290 gezelter 1628 if (.not.haveMixingMap) then
291 gezelter 2271 call createMixingMap()
292 gezelter 1628 endif
293    
294 gezelter 1608 ! Look up the correct parameters in the mixing matrix
295     #ifdef IS_MPI
296 gezelter 2271 atid1 = atid_Row(atom1)
297     atid2 = atid_Col(atom2)
298 gezelter 1608 #else
299 gezelter 2271 atid1 = atid(atom1)
300     atid2 = atid(atom2)
301 gezelter 1608 #endif
302    
303 gezelter 2271 ljt1 = LJMap%atidToLJtype(atid1)
304     ljt2 = LJMap%atidToLJtype(atid2)
305    
306 gezelter 2717 sigmai = MixingMap(ljt1,ljt2)%sigmai
307 gezelter 2271 epsilon = MixingMap(ljt1,ljt2)%epsilon
308     isSoftCore = MixingMap(ljt1,ljt2)%isSoftCore
309     shiftedPot = MixingMap(ljt1,ljt2)%shiftedPot
310 chrisfen 3129 shiftedFrc = MixingMap(ljt1,ljt2)%shiftedFrc
311 gezelter 2271
312 gezelter 2717 ros = rij * sigmai
313 gezelter 2204
314 gezelter 2717 if (isSoftCore) then
315 tim 2062
316 gezelter 2717 call getSoftFunc(ros, myPot, myDeriv)
317    
318 gezelter 2271 if (shiftedPot) then
319 gezelter 2717 rcos = rcut * sigmai
320     call getSoftFunc(rcos, myPotC, myDerivC)
321 chrisfen 3129 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 gezelter 2271 endif
330 gezelter 2717
331 tim 2062 else
332 gezelter 2717
333     call getLJfunc(ros, myPot, myDeriv)
334    
335 gezelter 2271 if (shiftedPot) then
336 gezelter 2717 rcos = rcut * sigmai
337 gezelter 3126 call getLJfunc(rcos, myPotC, myDerivC)
338 chrisfen 3129 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 gezelter 2204 endif
347 gezelter 2271
348 gezelter 1608 endif
349 gezelter 2204
350 chrisfen 3129 pot_temp = epsilon * (myPot - myPotC)
351     vpair = vpair + pot_temp
352 gezelter 3127 dudr = sw * epsilon * (myDeriv - myDerivC) * sigmai
353 gezelter 2717
354 gezelter 1608 drdx = d(1) / rij
355     drdy = d(2) / rij
356     drdz = d(3) / rij
357 gezelter 2204
358 gezelter 1608 fx = dudr * drdx
359     fy = dudr * drdy
360     fz = dudr * drdz
361 gezelter 2204
362 gezelter 1608 #ifdef IS_MPI
363     if (do_pot) then
364 gezelter 2361 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 gezelter 1608 endif
367 gezelter 2204
368 gezelter 1608 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 gezelter 2204
372 gezelter 1608 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 gezelter 2204
376 gezelter 1608 #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 gezelter 2204
383 gezelter 1608 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 gezelter 2204
388 gezelter 1608 #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 gezelter 2204
398 gezelter 1608 fpair(1) = fpair(1) + fx
399     fpair(2) = fpair(2) + fy
400     fpair(3) = fpair(3) + fz
401    
402     endif
403    
404     return
405 gezelter 2204
406 gezelter 1608 end subroutine do_lj_pair
407 gezelter 2204
408 chuckv 2189 subroutine destroyLJTypes()
409 gezelter 2271
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 chuckv 2189 haveMixingMap = .false.
424 gezelter 2717
425 chuckv 2189 end subroutine destroyLJTypes
426    
427 gezelter 2717 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 chrisfen 2727 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 gezelter 2717 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 chrisfen 2727 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 gezelter 2717 return
464     end subroutine getSoftFunc
465    
466 gezelter 1608 end module lj