ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90
Revision: 2319
Committed: Wed Sep 21 23:45:48 2005 UTC (18 years, 9 months ago) by chuckv
File size: 13978 byte(s)
Log Message:
Bug fix: If we are not using LJ (say we are using EAM), we probably shouldn't rebuild the LJ mixing map.

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