ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90
Revision: 2461
Committed: Mon Nov 21 22:59:02 2005 UTC (18 years, 9 months ago) by gezelter
File size: 12562 byte(s)
Log Message:
Cutoff mixing fixes have been made.

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.19 2005-11-21 22:59:01 gezelter Exp $, $Date: 2005-11-21 22:59:01 $, $Name: not supported by cvs2svn $, $Revision: 1.19 $
47
48
49 module lj
50 use atype_module
51 use vector_class
52 use simulation
53 use status
54 #ifdef IS_MPI
55 use mpiSimulation
56 #endif
57 use force_globals
58
59 implicit none
60 PRIVATE
61 #define __FORTRAN90
62 #include "UseTheForce/DarkSide/fInteractionMap.h"
63
64 integer, parameter :: DP = selected_real_kind(15)
65
66 logical, save :: useGeometricDistanceMixing = .false.
67 logical, save :: haveMixingMap = .false.
68
69 real(kind=DP), save :: defaultCutoff = 0.0_DP
70 logical, save :: defaultShift = .false.
71 logical, save :: haveDefaultCutoff = .false.
72
73
74 type, private :: LJtype
75 integer :: atid
76 real(kind=dp) :: sigma
77 real(kind=dp) :: epsilon
78 logical :: isSoftCore = .false.
79 end type LJtype
80
81 type, private :: LJList
82 integer :: Nljtypes = 0
83 integer :: currentLJtype = 0
84 type(LJtype), pointer :: LJtypes(:) => null()
85 integer, pointer :: atidToLJtype(:) => null()
86 end type LJList
87
88 type(LJList), save :: LJMap
89
90 type :: MixParameters
91 real(kind=DP) :: sigma
92 real(kind=DP) :: epsilon
93 real(kind=dp) :: sigma6
94 real(kind=dp) :: rCut
95 real(kind=dp) :: delta
96 logical :: rCutWasSet = .false.
97 logical :: shiftedPot
98 logical :: isSoftCore = .false.
99 end type MixParameters
100
101 type(MixParameters), dimension(:,:), allocatable :: MixingMap
102
103 public :: newLJtype
104 public :: setLJDefaultCutoff
105 public :: getSigma
106 public :: getEpsilon
107 public :: useGeometricMixing
108 public :: do_lj_pair
109 public :: destroyLJtypes
110
111 contains
112
113 subroutine newLJtype(c_ident, sigma, epsilon, isSoftCore, status)
114 integer,intent(in) :: c_ident
115 real(kind=dp),intent(in) :: sigma
116 real(kind=dp),intent(in) :: epsilon
117 integer, intent(in) :: isSoftCore
118 integer,intent(out) :: status
119 integer :: nLJTypes, ntypes, myATID
120 integer, pointer :: MatchList(:) => null()
121 integer :: current
122
123 status = 0
124 ! check to see if this is the first time into this routine...
125 if (.not.associated(LJMap%LJtypes)) then
126
127 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 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
143 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
154 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 !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 end subroutine setLJDefaultCutoff
165
166 function getSigma(atid) result (s)
167 integer, intent(in) :: atid
168 integer :: ljt1
169 real(kind=dp) :: s
170
171 if (LJMap%currentLJtype == 0) then
172 call handleError("LJ", "No members in LJMap")
173 return
174 end if
175
176 ljt1 = LJMap%atidToLJtype(atid)
177 s = LJMap%LJtypes(ljt1)%sigma
178
179 end function getSigma
180
181 function getEpsilon(atid) result (e)
182 integer, intent(in) :: atid
183 integer :: ljt1
184 real(kind=dp) :: e
185
186 if (LJMap%currentLJtype == 0) then
187 call handleError("LJ", "No members in LJMap")
188 return
189 end if
190
191 ljt1 = LJMap%atidToLJtype(atid)
192 e = LJMap%LJtypes(ljt1)%epsilon
193
194 end function getEpsilon
195
196 subroutine useGeometricMixing()
197 useGeometricDistanceMixing = .true.
198 haveMixingMap = .false.
199 return
200 end subroutine useGeometricMixing
201
202 subroutine createMixingMap()
203 integer :: nLJtypes, i, j
204 real ( kind = dp ) :: s1, s2, e1, e2
205 real ( kind = dp ) :: rcut6, tp6, tp12
206 logical :: isSoftCore1, isSoftCore2, doShift
207
208 if (LJMap%currentLJtype == 0) then
209 call handleError("LJ", "No members in LJMap")
210 return
211 end if
212
213 nLJtypes = LJMap%nLJtypes
214
215 if (.not. allocated(MixingMap)) then
216 allocate(MixingMap(nLJtypes, nLJtypes))
217 endif
218
219 do i = 1, nLJtypes
220
221 s1 = LJMap%LJtypes(i)%sigma
222 e1 = LJMap%LJtypes(i)%epsilon
223 isSoftCore1 = LJMap%LJtypes(i)%isSoftCore
224
225 do j = i, nLJtypes
226
227 s2 = LJMap%LJtypes(j)%sigma
228 e2 = LJMap%LJtypes(j)%epsilon
229 isSoftCore2 = LJMap%LJtypes(j)%isSoftCore
230
231 MixingMap(i,j)%isSoftCore = isSoftCore1 .or. isSoftCore2
232
233 ! only the distance parameter uses different mixing policies
234 if (useGeometricDistanceMixing) then
235 MixingMap(i,j)%sigma = dsqrt(s1 * s2)
236 else
237 MixingMap(i,j)%sigma = 0.5_dp * (s1 + s2)
238 endif
239
240 MixingMap(i,j)%epsilon = dsqrt(e1 * e2)
241
242 MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6
243
244 if (haveDefaultCutoff) then
245 rcut6 = defaultCutoff**6
246 tp6 = MixingMap(i,j)%sigma6/rcut6
247 tp12 = tp6**2
248 MixingMap(i,j)%delta =-4.0_DP*MixingMap(i,j)%epsilon*(tp12 - tp6)
249 MixingMap(i,j)%shiftedPot = defaultShift
250 else
251 MixingMap(i,j)%delta = 0.0_DP
252 MixingMap(i,j)%shiftedPot = defaultShift
253 endif
254
255 if (i.ne.j) then
256 MixingMap(j,i)%sigma = MixingMap(i,j)%sigma
257 MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
258 MixingMap(j,i)%sigma6 = MixingMap(i,j)%sigma6
259 MixingMap(j,i)%rCut = MixingMap(i,j)%rCut
260 MixingMap(j,i)%delta = MixingMap(i,j)%delta
261 MixingMap(j,i)%rCutWasSet = MixingMap(i,j)%rCutWasSet
262 MixingMap(j,i)%shiftedPot = MixingMap(i,j)%shiftedPot
263 MixingMap(j,i)%isSoftCore = MixingMap(i,j)%isSoftCore
264 endif
265
266 enddo
267 enddo
268
269 haveMixingMap = .true.
270
271 end subroutine createMixingMap
272
273 subroutine do_lj_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
274 pot, f, do_pot)
275
276 integer, intent(in) :: atom1, atom2
277 integer :: atid1, atid2, ljt1, ljt2
278 real( kind = dp ), intent(in) :: rij, r2, rcut
279 real( kind = dp ) :: pot, sw, vpair
280 real( kind = dp ), dimension(3,nLocal) :: f
281 real( kind = dp ), intent(in), dimension(3) :: d
282 real( kind = dp ), intent(inout), dimension(3) :: fpair
283 logical, intent(in) :: do_pot
284
285 ! local Variables
286 real( kind = dp ) :: drdx, drdy, drdz
287 real( kind = dp ) :: fx, fy, fz
288 real( kind = dp ) :: pot_temp, dudr
289 real( kind = dp ) :: sigma6
290 real( kind = dp ) :: epsilon
291 real( kind = dp ) :: r6, rc6
292 real( kind = dp ) :: t6, tp6
293 real( kind = dp ) :: t12, tp12
294 real( kind = dp ) :: delta
295 logical :: isSoftCore, shiftedPot
296 integer :: id1, id2, localError
297
298 if (.not.haveMixingMap) then
299 call createMixingMap()
300 endif
301
302 ! Look up the correct parameters in the mixing matrix
303 #ifdef IS_MPI
304 atid1 = atid_Row(atom1)
305 atid2 = atid_Col(atom2)
306 #else
307 atid1 = atid(atom1)
308 atid2 = atid(atom2)
309 #endif
310
311 ljt1 = LJMap%atidToLJtype(atid1)
312 ljt2 = LJMap%atidToLJtype(atid2)
313
314 sigma6 = MixingMap(ljt1,ljt2)%sigma6
315 epsilon = MixingMap(ljt1,ljt2)%epsilon
316 delta = MixingMap(ljt1,ljt2)%delta
317 isSoftCore = MixingMap(ljt1,ljt2)%isSoftCore
318 shiftedPot = MixingMap(ljt1,ljt2)%shiftedPot
319
320 r6 = r2 * r2 * r2
321
322 t6 = sigma6/ r6
323 t12 = t6 * t6
324
325 if (isSoftCore) then
326
327 pot_temp = 4.0E0_DP * epsilon * t6
328 if (shiftedPot) then
329 rc6 = rcut**6
330 tp6 = sigma6 / rc6
331 delta =-4.0_DP*epsilon*(tp6)
332 pot_temp = pot_temp + delta
333 endif
334
335 vpair = vpair + pot_temp
336
337 dudr = -sw * 24.0E0_DP * epsilon * t6 / rij
338
339 else
340 pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
341 if (shiftedPot) then
342 rc6 = rcut**6
343 tp6 = sigma6 / rc6
344 tp12 = tp6*tp6
345 delta =-4.0_DP*epsilon*(tp12 - tp6)
346 pot_temp = pot_temp + delta
347 endif
348
349 vpair = vpair + pot_temp
350
351 dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
352 endif
353
354 drdx = d(1) / rij
355 drdy = d(2) / rij
356 drdz = d(3) / rij
357
358 fx = dudr * drdx
359 fy = dudr * drdy
360 fz = dudr * drdz
361
362 #ifdef IS_MPI
363 if (do_pot) then
364 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 endif
367
368 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
372 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
376 #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
383 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
388 #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
398 fpair(1) = fpair(1) + fx
399 fpair(2) = fpair(2) + fy
400 fpair(3) = fpair(3) + fz
401
402 endif
403
404 return
405
406 end subroutine do_lj_pair
407
408 subroutine destroyLJTypes()
409
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 haveMixingMap = .false.
424 end subroutine destroyLJTypes
425
426 end module lj