ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/LJ.F90
Revision: 1386
Committed: Fri Oct 23 18:41:09 2009 UTC (15 years, 9 months ago) by gezelter
File size: 12485 byte(s)
Log Message:
removing MPI responsibilities from the lowest level force routines.  This is
in preparation for migrating these routines (LJ, electrostatic, eam, suttonchen,
gay-berne, sticky) to C++

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