ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/LJ.F90
Revision: 2733
Committed: Tue Apr 25 02:09:01 2006 UTC (19 years, 3 months ago) by gezelter
File size: 12676 byte(s)
Log Message:
Adding spherical boundary conditions to LD integrator

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