ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/LJ.F90
Revision: 939
Committed: Thu Apr 20 18:24:24 2006 UTC (19 years, 1 month ago) by gezelter
File size: 15402 byte(s)
Log Message:
Complete rewrite of spline code and everything that uses it.

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 939 !! @version $Id: LJ.F90,v 1.22 2006-04-20 18:24:24 gezelter Exp $, $Date: 2006-04-20 18:24:24 $, $Name: not supported by cvs2svn $, $Revision: 1.22 $
47 gezelter 115
48 gezelter 246
49 gezelter 115 module lj
50     use atype_module
51     use vector_class
52     use simulation
53 gezelter 135 use status
54 chuckv 798 use fForceOptions
55 gezelter 938 use interpolation
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 chuckv 131 integer, parameter :: DP = selected_real_kind(15)
67 gezelter 507
68 gezelter 572 logical, save :: useGeometricDistanceMixing = .false.
69     logical, save :: haveMixingMap = .false.
70 gezelter 938 logical, save :: useSplines = .false.
71 gezelter 572
72     real(kind=DP), save :: defaultCutoff = 0.0_DP
73     logical, save :: defaultShift = .false.
74     logical, save :: haveDefaultCutoff = .false.
75    
76     type, private :: LJtype
77     integer :: atid
78 gezelter 135 real(kind=dp) :: sigma
79     real(kind=dp) :: epsilon
80 gezelter 572 logical :: isSoftCore = .false.
81     end type LJtype
82 gezelter 507
83 gezelter 572 type, private :: LJList
84 chuckv 620 integer :: Nljtypes = 0
85 gezelter 572 integer :: currentLJtype = 0
86     type(LJtype), pointer :: LJtypes(:) => null()
87     integer, pointer :: atidToLJtype(:) => null()
88     end type LJList
89 gezelter 507
90 gezelter 572 type(LJList), save :: LJMap
91 gezelter 507
92 gezelter 135 type :: MixParameters
93     real(kind=DP) :: sigma
94     real(kind=DP) :: epsilon
95 gezelter 938 real(kind=dp) :: sigmai
96 gezelter 572 real(kind=dp) :: rCut
97     logical :: rCutWasSet = .false.
98     logical :: shiftedPot
99     logical :: isSoftCore = .false.
100 gezelter 135 end type MixParameters
101 gezelter 507
102 gezelter 135 type(MixParameters), dimension(:,:), allocatable :: MixingMap
103 gezelter 507
104 gezelter 938 type(cubicSpline), save :: vLJspline
105     type(cubicSpline), save :: vLJpspline
106     type(cubicSpline), save :: vSoftSpline
107     type(cubicSpline), save :: vSoftpSpline
108    
109 gezelter 572 public :: newLJtype
110     public :: setLJDefaultCutoff
111     public :: getSigma
112     public :: getEpsilon
113 gezelter 135 public :: do_lj_pair
114 gezelter 572 public :: destroyLJtypes
115 gezelter 938 public :: setLJsplineRmax
116 gezelter 507
117 gezelter 115 contains
118    
119 gezelter 572 subroutine newLJtype(c_ident, sigma, epsilon, isSoftCore, status)
120 gezelter 246 integer,intent(in) :: c_ident
121 gezelter 135 real(kind=dp),intent(in) :: sigma
122     real(kind=dp),intent(in) :: epsilon
123 gezelter 572 integer, intent(in) :: isSoftCore
124 chuckv 131 integer,intent(out) :: status
125 gezelter 572 integer :: nLJTypes, ntypes, myATID
126 gezelter 246 integer, pointer :: MatchList(:) => null()
127 gezelter 572 integer :: current
128 gezelter 135
129 chuckv 131 status = 0
130 gezelter 572 ! check to see if this is the first time into this routine...
131     if (.not.associated(LJMap%LJtypes)) then
132 gezelter 507
133 gezelter 572 call getMatchingElementList(atypes, "is_LennardJones", .true., &
134     nLJTypes, MatchList)
135    
136     LJMap%nLJtypes = nLJTypes
137    
138     allocate(LJMap%LJtypes(nLJTypes))
139    
140     ntypes = getSize(atypes)
141    
142     allocate(LJMap%atidToLJtype(ntypes))
143     end if
144    
145     LJMap%currentLJtype = LJMap%currentLJtype + 1
146     current = LJMap%currentLJtype
147    
148 gezelter 246 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
149 gezelter 572 LJMap%atidToLJtype(myATID) = current
150     LJMap%LJtypes(current)%atid = myATID
151     LJMap%LJtypes(current)%sigma = sigma
152     LJMap%LJtypes(current)%epsilon = epsilon
153     if (isSoftCore .eq. 1) then
154     LJMap%LJtypes(current)%isSoftCore = .true.
155     else
156     LJMap%LJtypes(current)%isSoftCore = .false.
157     endif
158     end subroutine newLJtype
159 chuckv 131
160 gezelter 572 subroutine setLJDefaultCutoff(thisRcut, shiftedPot)
161     real(kind=dp), intent(in) :: thisRcut
162     logical, intent(in) :: shiftedPot
163     defaultCutoff = thisRcut
164     defaultShift = shiftedPot
165     haveDefaultCutoff = .true.
166 gezelter 938 !we only want to build LJ Mixing map and spline if LJ is being used.
167 chuckv 620 if(LJMap%nLJTypes /= 0) then
168     call createMixingMap()
169 gezelter 938 call setLJsplineRmax(defaultCutoff)
170 chuckv 620 end if
171 gezelter 938
172 gezelter 572 end subroutine setLJDefaultCutoff
173 gezelter 507
174 gezelter 140 function getSigma(atid) result (s)
175     integer, intent(in) :: atid
176 gezelter 572 integer :: ljt1
177 gezelter 140 real(kind=dp) :: s
178 gezelter 507
179 gezelter 572 if (LJMap%currentLJtype == 0) then
180     call handleError("LJ", "No members in LJMap")
181 gezelter 140 return
182     end if
183 gezelter 507
184 gezelter 572 ljt1 = LJMap%atidToLJtype(atid)
185     s = LJMap%LJtypes(ljt1)%sigma
186    
187 gezelter 140 end function getSigma
188    
189     function getEpsilon(atid) result (e)
190     integer, intent(in) :: atid
191 gezelter 572 integer :: ljt1
192 gezelter 140 real(kind=dp) :: e
193 gezelter 507
194 gezelter 572 if (LJMap%currentLJtype == 0) then
195     call handleError("LJ", "No members in LJMap")
196 gezelter 140 return
197     end if
198 gezelter 507
199 gezelter 572 ljt1 = LJMap%atidToLJtype(atid)
200     e = LJMap%LJtypes(ljt1)%epsilon
201    
202 gezelter 140 end function getEpsilon
203    
204 gezelter 572 subroutine createMixingMap()
205     integer :: nLJtypes, i, j
206     real ( kind = dp ) :: s1, s2, e1, e2
207     real ( kind = dp ) :: rcut6, tp6, tp12
208 gezelter 590 logical :: isSoftCore1, isSoftCore2, doShift
209 gezelter 135
210 gezelter 572 if (LJMap%currentLJtype == 0) then
211     call handleError("LJ", "No members in LJMap")
212 gezelter 402 return
213 gezelter 572 end if
214 gezelter 507
215 gezelter 572 nLJtypes = LJMap%nLJtypes
216 gezelter 507
217 gezelter 402 if (.not. allocated(MixingMap)) then
218 gezelter 572 allocate(MixingMap(nLJtypes, nLJtypes))
219 gezelter 402 endif
220    
221 chuckv 798 useGeometricDistanceMixing = usesGeometricDistanceMixing()
222 gezelter 572 do i = 1, nLJtypes
223 gezelter 507
224 gezelter 572 s1 = LJMap%LJtypes(i)%sigma
225     e1 = LJMap%LJtypes(i)%epsilon
226     isSoftCore1 = LJMap%LJtypes(i)%isSoftCore
227 gezelter 507
228 gezelter 572 do j = i, nLJtypes
229    
230     s2 = LJMap%LJtypes(j)%sigma
231     e2 = LJMap%LJtypes(j)%epsilon
232     isSoftCore2 = LJMap%LJtypes(j)%isSoftCore
233    
234     MixingMap(i,j)%isSoftCore = isSoftCore1 .or. isSoftCore2
235 gezelter 507
236 gezelter 572 ! only the distance parameter uses different mixing policies
237     if (useGeometricDistanceMixing) then
238     MixingMap(i,j)%sigma = dsqrt(s1 * s2)
239     else
240     MixingMap(i,j)%sigma = 0.5_dp * (s1 + s2)
241     endif
242    
243     MixingMap(i,j)%epsilon = dsqrt(e1 * e2)
244 gezelter 507
245 gezelter 938 MixingMap(i,j)%sigmai = 1.0_DP / (MixingMap(i,j)%sigma)
246 gezelter 507
247 gezelter 762 if (haveDefaultCutoff) then
248     MixingMap(i,j)%shiftedPot = defaultShift
249 gezelter 572 else
250 gezelter 762 MixingMap(i,j)%shiftedPot = defaultShift
251     endif
252 gezelter 507
253 gezelter 590 if (i.ne.j) then
254     MixingMap(j,i)%sigma = MixingMap(i,j)%sigma
255     MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
256 gezelter 938 MixingMap(j,i)%sigmai = MixingMap(i,j)%sigmai
257 gezelter 590 MixingMap(j,i)%rCut = MixingMap(i,j)%rCut
258     MixingMap(j,i)%rCutWasSet = MixingMap(i,j)%rCutWasSet
259     MixingMap(j,i)%shiftedPot = MixingMap(i,j)%shiftedPot
260     MixingMap(j,i)%isSoftCore = MixingMap(i,j)%isSoftCore
261     endif
262    
263 gezelter 572 enddo
264     enddo
265    
266 chrisfen 217 haveMixingMap = .true.
267 gezelter 572
268 gezelter 135 end subroutine createMixingMap
269 gezelter 572
270 gezelter 938 subroutine setLJsplineRmax(largestRcut)
271     real( kind = dp ), intent(in) :: largestRcut
272     real( kind = dp ) :: s, bigS, smallS, rmax, rmin
273     integer :: np, i
274    
275     if (LJMap%nLJtypes .ne. 0) then
276    
277     !
278     ! find the largest and smallest values of sigma that we'll need
279     !
280     bigS = 0.0_DP
281     smallS = 1.0e9
282     do i = 1, LJMap%nLJtypes
283     s = LJMap%LJtypes(i)%sigma
284     if (s .gt. bigS) bigS = s
285     if (s .lt. smallS) smallS = s
286     enddo
287    
288     !
289     ! give ourselves a 20% margin just in case
290     !
291     rmax = 1.2 * largestRcut / smallS
292     !
293     ! assume atoms will never get closer than 1 angstrom
294     !
295     rmin = 1 / bigS
296     !
297     ! assume 500 points is enough
298     !
299     np = 500
300    
301     write(*,*) 'calling setupSplines with rmin = ', rmin, ' rmax = ', rmax, &
302     ' np = ', np
303    
304     call setupSplines(rmin, rmax, np)
305    
306     endif
307     return
308     end subroutine setLJsplineRmax
309    
310 gezelter 762 subroutine do_lj_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
311 gezelter 115 pot, f, do_pot)
312 gezelter 762
313 gezelter 115 integer, intent(in) :: atom1, atom2
314 gezelter 572 integer :: atid1, atid2, ljt1, ljt2
315 gezelter 762 real( kind = dp ), intent(in) :: rij, r2, rcut
316 gezelter 115 real( kind = dp ) :: pot, sw, vpair
317     real( kind = dp ), dimension(3,nLocal) :: f
318     real( kind = dp ), intent(in), dimension(3) :: d
319     real( kind = dp ), intent(inout), dimension(3) :: fpair
320     logical, intent(in) :: do_pot
321    
322     ! local Variables
323     real( kind = dp ) :: drdx, drdy, drdz
324     real( kind = dp ) :: fx, fy, fz
325 gezelter 938 real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos
326 gezelter 115 real( kind = dp ) :: pot_temp, dudr
327 gezelter 938 real( kind = dp ) :: sigmai
328 gezelter 115 real( kind = dp ) :: epsilon
329 gezelter 572 logical :: isSoftCore, shiftedPot
330 gezelter 135 integer :: id1, id2, localError
331 gezelter 115
332 gezelter 135 if (.not.haveMixingMap) then
333 gezelter 572 call createMixingMap()
334 gezelter 135 endif
335    
336 gezelter 115 ! Look up the correct parameters in the mixing matrix
337     #ifdef IS_MPI
338 gezelter 572 atid1 = atid_Row(atom1)
339     atid2 = atid_Col(atom2)
340 gezelter 115 #else
341 gezelter 572 atid1 = atid(atom1)
342     atid2 = atid(atom2)
343 gezelter 115 #endif
344    
345 gezelter 572 ljt1 = LJMap%atidToLJtype(atid1)
346     ljt2 = LJMap%atidToLJtype(atid2)
347    
348 gezelter 938 sigmai = MixingMap(ljt1,ljt2)%sigmai
349 gezelter 572 epsilon = MixingMap(ljt1,ljt2)%epsilon
350     isSoftCore = MixingMap(ljt1,ljt2)%isSoftCore
351     shiftedPot = MixingMap(ljt1,ljt2)%shiftedPot
352    
353 gezelter 938 ros = rij * sigmai
354     myPotC = 0.0_DP
355 gezelter 507
356 gezelter 938 if (isSoftCore) then
357 tim 378
358 gezelter 938 call getSoftFunc(ros, myPot, myDeriv)
359    
360 gezelter 572 if (shiftedPot) then
361 gezelter 938 rcos = rcut * sigmai
362     call getSoftFunc(rcos, myPotC, myDerivC)
363 gezelter 572 endif
364 gezelter 938
365 tim 378 else
366 gezelter 938
367     call getLJfunc(ros, myPot, myDeriv)
368    
369 gezelter 572 if (shiftedPot) then
370 gezelter 938 rcos = rcut * sigmai
371     call getLJfunc(rcos, myPotC, myDerivC)
372 gezelter 507 endif
373 gezelter 572
374 gezelter 115 endif
375 gezelter 507
376 gezelter 938 !write(*,*) rij, ros, rcos, myPot, myDeriv, myPotC
377    
378     pot_temp = epsilon * (myPot - myPotC)
379     vpair = vpair + pot_temp
380     dudr = sw * epsilon * myDeriv * sigmai
381    
382 gezelter 115 drdx = d(1) / rij
383     drdy = d(2) / rij
384     drdz = d(3) / rij
385 gezelter 507
386 gezelter 115 fx = dudr * drdx
387     fy = dudr * drdy
388     fz = dudr * drdz
389 gezelter 507
390 gezelter 115 #ifdef IS_MPI
391     if (do_pot) then
392 gezelter 662 pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5
393     pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5
394 gezelter 115 endif
395 gezelter 507
396 gezelter 115 f_Row(1,atom1) = f_Row(1,atom1) + fx
397     f_Row(2,atom1) = f_Row(2,atom1) + fy
398     f_Row(3,atom1) = f_Row(3,atom1) + fz
399 gezelter 507
400 gezelter 115 f_Col(1,atom2) = f_Col(1,atom2) - fx
401     f_Col(2,atom2) = f_Col(2,atom2) - fy
402     f_Col(3,atom2) = f_Col(3,atom2) - fz
403 gezelter 507
404 gezelter 115 #else
405     if (do_pot) pot = pot + sw*pot_temp
406    
407     f(1,atom1) = f(1,atom1) + fx
408     f(2,atom1) = f(2,atom1) + fy
409     f(3,atom1) = f(3,atom1) + fz
410 gezelter 507
411 gezelter 115 f(1,atom2) = f(1,atom2) - fx
412     f(2,atom2) = f(2,atom2) - fy
413     f(3,atom2) = f(3,atom2) - fz
414     #endif
415 gezelter 507
416 gezelter 115 #ifdef IS_MPI
417     id1 = AtomRowToGlobal(atom1)
418     id2 = AtomColToGlobal(atom2)
419     #else
420     id1 = atom1
421     id2 = atom2
422     #endif
423    
424     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
425 gezelter 507
426 gezelter 115 fpair(1) = fpair(1) + fx
427     fpair(2) = fpair(2) + fy
428     fpair(3) = fpair(3) + fz
429    
430     endif
431    
432     return
433 gezelter 507
434 gezelter 115 end subroutine do_lj_pair
435 gezelter 507
436 chuckv 492 subroutine destroyLJTypes()
437 gezelter 572
438     LJMap%nLJtypes = 0
439     LJMap%currentLJtype = 0
440    
441     if (associated(LJMap%LJtypes)) then
442     deallocate(LJMap%LJtypes)
443     LJMap%LJtypes => null()
444     end if
445    
446     if (associated(LJMap%atidToLJtype)) then
447     deallocate(LJMap%atidToLJtype)
448     LJMap%atidToLJtype => null()
449     end if
450    
451 chuckv 492 haveMixingMap = .false.
452 gezelter 938
453     call deleteSpline(vLJspline)
454     call deleteSpline(vLJpspline)
455     call deleteSpline(vSoftSpline)
456     call deleteSpline(vSoftpSpline)
457    
458 chuckv 492 end subroutine destroyLJTypes
459    
460 gezelter 938 subroutine getLJfunc(r, myPot, myDeriv)
461    
462     real(kind=dp), intent(in) :: r
463     real(kind=dp), intent(inout) :: myPot, myDeriv
464     real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13
465     real(kind=dp) :: a, b, c, d, dx
466     integer :: j
467    
468     if (useSplines) then
469    
470 gezelter 939 call lookupUniformSpline1d(vLJSpline, r, myPot, myDeriv)
471 gezelter 938
472     else
473     ri = 1.0_DP / r
474     ri2 = ri*ri
475     ri6 = ri2*ri2*ri2
476     ri7 = ri6*ri
477     ri12 = ri6*ri6
478     ri13 = ri12*ri
479    
480     myPot = 4.0_DP * (ri12 - ri6)
481     myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13)
482     endif
483    
484     return
485     end subroutine getLJfunc
486    
487     subroutine getSoftFunc(r, myPot, myDeriv)
488    
489     real(kind=dp), intent(in) :: r
490     real(kind=dp), intent(inout) :: myPot, myDeriv
491     real(kind=dp) :: ri, ri2, ri6, ri7
492     real(kind=dp) :: a, b, c, d, dx
493     integer :: j
494    
495     if (useSplines) then
496    
497 gezelter 939 call lookupUniformSpline1d(vSoftSpline, r, myPot, myDeriv)
498 gezelter 938
499     else
500     ri = 1.0_DP / r
501     ri2 = ri*ri
502     ri6 = ri2*ri2*ri2
503     ri7 = ri6*ri
504     myPot = 4.0_DP * (ri6)
505     myDeriv = - 24.0_DP * ri7
506     endif
507    
508     return
509     end subroutine getSoftFunc
510    
511     subroutine setupSplines(rmin, rmax, np)
512     real( kind = dp ), intent(in) :: rmin, rmax
513     integer, intent(in) :: np
514     real( kind = dp ) :: rvals(np), vLJ(np), vLJp(np), vSoft(np), vSoftp(np)
515     real( kind = dp ) :: dr, r, ri, ri2, ri6, ri7, ri12, ri13
516     integer :: i
517    
518     dr = (rmax-rmin) / float(np-1)
519    
520     do i = 1, np
521     r = rmin + dble(i-1)*dr
522     ri = 1.0_DP / r
523     ri2 = ri*ri
524     ri6 = ri2*ri2*ri2
525     ri7 = ri6*ri
526     ri12 = ri6*ri6
527     ri13 = ri12*ri
528    
529     rvals(i) = r
530     vLJ(i) = 4.0_DP * (ri12 - ri6)
531     vLJp(i) = 24.0_DP * (ri7 - 2.0_DP * ri13)
532    
533     vSoft(i) = 4.0_DP * (ri6)
534     vSoftp(i) = - 24.0_DP * ri7
535     enddo
536    
537 gezelter 939 call newSpline(vLJspline, rvals, vLJ, .true.)
538     call newSpline(vLJpspline, rvals, vLJp, .true.)
539     call newSpline(vSoftSpline, rvals, vSoft, .true.)
540     call newSpline(vSoftpSpline, rvals, vSoftp, .true.)
541 gezelter 938
542     return
543     end subroutine setupSplines
544 gezelter 115 end module lj