ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/LJ.F90
Revision: 2062
Committed: Fri Feb 25 21:22:00 2005 UTC (19 years, 6 months ago) by tim
File size: 12416 byte(s)
Log Message:
adding soft potential to LJ Module

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.8 2005-02-25 21:22:00 tim Exp $, $Date: 2005-02-25 21:22:00 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $
47
48
49 module lj
50 use atype_module
51 use switcheroo
52 use vector_class
53 use simulation
54 use status
55 #ifdef IS_MPI
56 use mpiSimulation
57 #endif
58 use force_globals
59
60 implicit none
61 PRIVATE
62
63 integer, parameter :: DP = selected_real_kind(15)
64
65 type, private :: LjType
66 integer :: c_ident
67 integer :: atid
68 real(kind=dp) :: sigma
69 real(kind=dp) :: epsilon
70 logical :: soft_pot
71 end type LjType
72
73 type(LjType), dimension(:), allocatable :: ParameterMap
74
75 logical, save :: haveMixingMap = .false.
76
77 type :: MixParameters
78 real(kind=DP) :: sigma
79 real(kind=DP) :: epsilon
80 real(kind=dp) :: sigma6
81 real(kind=dp) :: tp6
82 real(kind=dp) :: tp12
83 real(kind=dp) :: delta
84 logical :: soft_pot
85 end type MixParameters
86
87 type(MixParameters), dimension(:,:), allocatable :: MixingMap
88
89 real(kind=DP), save :: LJ_rcut
90 logical, save :: have_rcut = .false.
91 logical, save :: LJ_do_shift = .false.
92 logical, save :: useGeometricDistanceMixing = .false.
93
94 !! Public methods and data
95
96 public :: setCutoffLJ
97 public :: useGeometricMixing
98 public :: do_lj_pair
99 public :: newLJtype
100 public :: getSigma
101 public :: getEpsilon
102
103 contains
104
105 subroutine newLJtype(c_ident, sigma, epsilon, soft_pot, status)
106 integer,intent(in) :: c_ident
107 real(kind=dp),intent(in) :: sigma
108 real(kind=dp),intent(in) :: epsilon
109 integer, intent(in) :: soft_pot
110 integer,intent(out) :: status
111 integer :: nATypes, myATID
112 integer, pointer :: MatchList(:) => null()
113
114 status = 0
115
116 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
117
118 if (.not.allocated(ParameterMap)) then
119
120 !call getMatchingElementList(atypes, "is_LennardJones", .true., &
121 ! nLJTypes, MatchList)
122 nAtypes = getSize(atypes)
123 if (nAtypes == 0) then
124 status = -1
125 return
126 end if
127
128 if (.not. allocated(ParameterMap)) then
129 allocate(ParameterMap(nAtypes))
130 endif
131
132 end if
133
134 if (myATID .gt. size(ParameterMap)) then
135 status = -1
136 return
137 endif
138
139 ! set the values for ParameterMap for this atom type:
140
141 ParameterMap(myATID)%c_ident = c_ident
142 ParameterMap(myATID)%atid = myATID
143 ParameterMap(myATID)%epsilon = epsilon
144 ParameterMap(myATID)%sigma = sigma
145 if (soft_pot == 1) then
146 ParameterMap(myATID)%soft_pot = .true.
147 else
148 ParameterMap(myATID)%soft_pot = .false.
149 endif
150 end subroutine newLJtype
151
152 function getSigma(atid) result (s)
153 integer, intent(in) :: atid
154 integer :: localError
155 real(kind=dp) :: s
156
157 if (.not.allocated(ParameterMap)) then
158 call handleError("LJ", "no ParameterMap was present before first call of getSigma!")
159 return
160 end if
161
162 s = ParameterMap(atid)%sigma
163 end function getSigma
164
165 function getEpsilon(atid) result (e)
166 integer, intent(in) :: atid
167 integer :: localError
168 real(kind=dp) :: e
169
170 if (.not.allocated(ParameterMap)) then
171 call handleError("LJ", "no ParameterMap was present before first call of getEpsilon!")
172 return
173 end if
174
175 e = ParameterMap(atid)%epsilon
176 end function getEpsilon
177
178
179 subroutine setCutoffLJ(rcut, do_shift, status)
180 logical, intent(in):: do_shift
181 integer :: status, myStatus
182 real(kind=dp) :: rcut
183
184 #define __FORTRAN90
185 #include "UseTheForce/fSwitchingFunction.h"
186
187 status = 0
188
189 LJ_rcut = rcut
190 LJ_do_shift = do_shift
191 call set_switch(LJ_SWITCH, rcut, rcut)
192 have_rcut = .true.
193
194 return
195 end subroutine setCutoffLJ
196
197 subroutine useGeometricMixing()
198 useGeometricDistanceMixing = .true.
199 haveMixingMap = .false.
200 return
201 end subroutine useGeometricMixing
202
203 subroutine createMixingMap(status)
204 integer :: nATIDs
205 integer :: status
206 integer :: i
207 integer :: j
208 real ( kind = dp ) :: Sigma_i, Sigma_j
209 real ( kind = dp ) :: Epsilon_i, Epsilon_j
210 real ( kind = dp ) :: rcut6
211
212 status = 0
213
214 nATIDs = size(ParameterMap)
215
216 if (nATIDs == 0) then
217 status = -1
218 return
219 end if
220
221 if (.not.have_rcut) then
222 status = -1
223 return
224 endif
225
226 if (.not. allocated(MixingMap)) then
227 allocate(MixingMap(nATIDs, nATIDs))
228 endif
229
230 rcut6 = LJ_rcut**6
231
232 ! This loops through all atypes, even those that don't support LJ forces.
233 do i = 1, nATIDs
234
235 Epsilon_i = ParameterMap(i)%epsilon
236 Sigma_i = ParameterMap(i)%sigma
237
238 ! do self mixing rule
239 MixingMap(i,i)%sigma = Sigma_i
240 MixingMap(i,i)%sigma6 = Sigma_i ** 6
241 MixingMap(i,i)%tp6 = (MixingMap(i,i)%sigma6)/rcut6
242 MixingMap(i,i)%tp12 = (MixingMap(i,i)%tp6) ** 2
243 MixingMap(i,i)%epsilon = Epsilon_i
244 MixingMap(i,i)%delta = -4.0_DP * MixingMap(i,i)%epsilon * &
245 (MixingMap(i,i)%tp12 - MixingMap(i,i)%tp6)
246 MixingMap(i,i)%soft_pot = ParameterMap(i)%soft_pot
247
248 do j = i + 1, nATIDs
249
250 Epsilon_j = ParameterMap(j)%epsilon
251 Sigma_j = ParameterMap(j)%sigma
252
253 ! only the distance parameter uses different mixing policies
254 if (useGeometricDistanceMixing) then
255 ! only for OPLS as far as we can tell
256 MixingMap(i,j)%sigma = dsqrt(Sigma_i * Sigma_j)
257 else
258 ! everyone else
259 MixingMap(i,j)%sigma = 0.5_dp * (Sigma_i + Sigma_j)
260 endif
261
262 ! energy parameter is always geometric mean:
263 MixingMap(i,j)%epsilon = dsqrt(Epsilon_i * Epsilon_j)
264
265 MixingMap(i,j)%sigma6 = (MixingMap(i,j)%sigma)**6
266 MixingMap(i,j)%tp6 = MixingMap(i,j)%sigma6/rcut6
267 MixingMap(i,j)%tp12 = (MixingMap(i,j)%tp6) ** 2
268
269 MixingMap(i,j)%delta = -4.0_DP * MixingMap(i,j)%epsilon * &
270 (MixingMap(i,j)%tp12 - MixingMap(i,j)%tp6)
271
272 MixingMap(i,j)%soft_pot = ParameterMap(i)%soft_pot .or. ParameterMap(j)%soft_pot
273
274
275 MixingMap(j,i)%sigma = MixingMap(i,j)%sigma
276 MixingMap(j,i)%sigma6 = MixingMap(i,j)%sigma6
277 MixingMap(j,i)%tp6 = MixingMap(i,j)%tp6
278 MixingMap(j,i)%tp12 = MixingMap(i,j)%tp12
279 MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
280 MixingMap(j,i)%delta = MixingMap(i,j)%delta
281 MixingMap(j,i)%soft_pot = MixingMap(i,j)%soft_pot
282
283 end do
284 end do
285
286 haveMixingMap = .true.
287
288 end subroutine createMixingMap
289
290 subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
291 pot, f, do_pot)
292
293 integer, intent(in) :: atom1, atom2
294 real( kind = dp ), intent(in) :: rij, r2
295 real( kind = dp ) :: pot, sw, vpair
296 real( kind = dp ), dimension(3,nLocal) :: f
297 real( kind = dp ), intent(in), dimension(3) :: d
298 real( kind = dp ), intent(inout), dimension(3) :: fpair
299 logical, intent(in) :: do_pot
300
301 ! local Variables
302 real( kind = dp ) :: drdx, drdy, drdz
303 real( kind = dp ) :: fx, fy, fz
304 real( kind = dp ) :: pot_temp, dudr
305 real( kind = dp ) :: sigma6
306 real( kind = dp ) :: epsilon
307 real( kind = dp ) :: r6
308 real( kind = dp ) :: t6
309 real( kind = dp ) :: t12
310 real( kind = dp ) :: delta
311 logical :: soft_pot
312 integer :: id1, id2, localError
313
314 if (.not.haveMixingMap) then
315 localError = 0
316 call createMixingMap(localError)
317 if ( localError .ne. 0 ) then
318 call handleError("LJ", "MixingMap creation failed!")
319 return
320 end if
321 endif
322
323 ! Look up the correct parameters in the mixing matrix
324 #ifdef IS_MPI
325 sigma6 = MixingMap(atid_Row(atom1),atid_Col(atom2))%sigma6
326 epsilon = MixingMap(atid_Row(atom1),atid_Col(atom2))%epsilon
327 delta = MixingMap(atid_Row(atom1),atid_Col(atom2))%delta
328 soft_pot = MixingMap(atid_Row(atom1),atid_Col(atom2))%soft_pot
329 #else
330 sigma6 = MixingMap(atid(atom1),atid(atom2))%sigma6
331 epsilon = MixingMap(atid(atom1),atid(atom2))%epsilon
332 delta = MixingMap(atid(atom1),atid(atom2))%delta
333 soft_pot = MixingMap(atid(atom1),atid(atom2))%soft_pot
334 #endif
335
336 r6 = r2 * r2 * r2
337
338 t6 = sigma6/ r6
339 t12 = t6 * t6
340
341 if (soft_pot) then
342 pot_temp = 4.0E0_DP * epsilon * t12
343 if (LJ_do_shift) then
344 pot_temp = pot_temp + delta
345 endif
346
347 vpair = vpair + pot_temp
348
349 dudr = sw * 24.0E0_DP * epsilon * (-2.0E0_DP)*t12 / rij
350
351 else
352 pot_temp = 4.0E0_DP * epsilon * (t12 - t6)
353 if (LJ_do_shift) then
354 pot_temp = pot_temp + delta
355 endif
356
357 vpair = vpair + pot_temp
358
359 dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
360 endif
361
362 drdx = d(1) / rij
363 drdy = d(2) / rij
364 drdz = d(3) / rij
365
366 fx = dudr * drdx
367 fy = dudr * drdy
368 fz = dudr * drdz
369
370
371 #ifdef IS_MPI
372 if (do_pot) then
373 pot_Row(atom1) = pot_Row(atom1) + sw*pot_temp*0.5
374 pot_Col(atom2) = pot_Col(atom2) + sw*pot_temp*0.5
375 endif
376
377 f_Row(1,atom1) = f_Row(1,atom1) + fx
378 f_Row(2,atom1) = f_Row(2,atom1) + fy
379 f_Row(3,atom1) = f_Row(3,atom1) + fz
380
381 f_Col(1,atom2) = f_Col(1,atom2) - fx
382 f_Col(2,atom2) = f_Col(2,atom2) - fy
383 f_Col(3,atom2) = f_Col(3,atom2) - fz
384
385 #else
386 if (do_pot) pot = pot + sw*pot_temp
387
388 f(1,atom1) = f(1,atom1) + fx
389 f(2,atom1) = f(2,atom1) + fy
390 f(3,atom1) = f(3,atom1) + fz
391
392 f(1,atom2) = f(1,atom2) - fx
393 f(2,atom2) = f(2,atom2) - fy
394 f(3,atom2) = f(3,atom2) - fz
395 #endif
396
397 #ifdef IS_MPI
398 id1 = AtomRowToGlobal(atom1)
399 id2 = AtomColToGlobal(atom2)
400 #else
401 id1 = atom1
402 id2 = atom2
403 #endif
404
405 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
406
407 fpair(1) = fpair(1) + fx
408 fpair(2) = fpair(2) + fy
409 fpair(3) = fpair(3) + fz
410
411 endif
412
413 return
414
415 end subroutine do_lj_pair
416
417
418 !! Calculates the mixing for sigma or epslon
419
420 end module lj