ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/LJ.F90
Revision: 1286
Committed: Wed Sep 10 17:57:55 2008 UTC (16 years, 10 months ago) by gezelter
File size: 13620 byte(s)
Log Message:
Added support for scaled 1-2, 1-3 and 1-4 interactions.

File Contents

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