ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/eam.F90
Revision: 1573
Committed: Mon May 30 18:59:36 2011 UTC (14 years, 1 month ago) by gezelter
File size: 14811 byte(s)
Log Message:
fixed a cutoff bug (i think) in calc_eam_prepair_rho

File Contents

# Content
1 !!
2 !! Copyright (c) 2005, 2009 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. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41 module eam
42 use definitions
43 use simulation
44 use force_globals
45 use status
46 use atype_module
47 use vector_class
48 use interpolation
49 implicit none
50 PRIVATE
51 #define __FORTRAN90
52 #include "UseTheForce/DarkSide/fInteractionMap.h"
53
54 logical, save :: EAM_FF_initialized = .false.
55 integer, save :: EAM_Mixing_Policy
56 real(kind = dp), save :: EAM_rcut
57 logical, save :: haveRcut = .false.
58
59 character(len = statusMsgSize) :: errMesg
60 integer :: eam_err
61
62 character(len = 200) :: errMsg
63 character(len=*), parameter :: RoutineName = "EAM MODULE"
64 !! Logical that determines if eam arrays should be zeroed
65 logical :: cleanme = .true.
66
67 type, private :: EAMtype
68 integer :: eam_atype
69 real( kind = DP ) :: eam_lattice
70 real( kind = DP ) :: eam_rcut
71 integer :: eam_atype_map
72 type(cubicSpline) :: rho
73 type(cubicSpline) :: Z
74 type(cubicSpline) :: F
75 type(cubicSpline) :: phi
76 end type EAMtype
77
78
79 type, private :: EAMTypeList
80 integer :: n_eam_types = 0
81 integer :: currentAddition = 0
82
83 type (EAMtype), pointer :: EAMParams(:) => null()
84 integer, pointer :: atidToEAMType(:) => null()
85 end type EAMTypeList
86
87 type (eamTypeList), save :: EAMList
88
89 !! standard eam stuff
90
91
92 public :: setCutoffEAM
93 public :: do_eam_pair
94 public :: newEAMtype
95 public :: calc_eam_prepair_rho
96 public :: calc_eam_preforce_Frho
97 public :: destroyEAMTypes
98 public :: getEAMCut
99 public :: lookupEAMSpline
100 public :: lookupEAMSpline1d
101
102 contains
103
104 subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
105 eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho, c_ident, status)
106 real (kind = dp ) :: lattice_constant
107 integer :: eam_nrho
108 real (kind = dp ) :: eam_drho
109 integer :: eam_nr
110 real (kind = dp ) :: eam_dr
111 real (kind = dp ) :: rcut
112 real (kind = dp ), dimension(eam_nr) :: eam_Z_r, rvals
113 real (kind = dp ), dimension(eam_nr) :: eam_rho_r, eam_phi_r
114 real (kind = dp ), dimension(eam_nrho) :: eam_F_rho, rhovals
115 integer :: c_ident
116 integer :: status
117
118 integer :: nAtypes,nEAMTypes,myATID
119 integer :: maxVals
120 integer :: alloc_stat
121 integer :: current, j
122 integer,pointer :: Matchlist(:) => null()
123
124 status = 0
125
126 !! Assume that atypes has already been set and get the total number of types in atypes
127 !! Also assume that every member of atypes is a EAM model.
128
129 ! check to see if this is the first time into
130 if (.not.associated(EAMList%EAMParams)) then
131 call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList)
132 EAMList%n_eam_types = nEAMtypes
133 allocate(EAMList%EAMParams(nEAMTypes))
134 nAtypes = getSize(atypes)
135 allocate(EAMList%atidToEAMType(nAtypes))
136 end if
137
138 EAMList%currentAddition = EAMList%currentAddition + 1
139 current = EAMList%currentAddition
140
141 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
142 EAMList%atidToEAMType(myATID) = current
143
144 EAMList%EAMParams(current)%eam_atype = c_ident
145 EAMList%EAMParams(current)%eam_lattice = lattice_constant
146 EAMList%EAMParams(current)%eam_rcut = rcut
147
148 ! Build array of r values
149 do j = 1, eam_nr
150 rvals(j) = real(j-1,kind=dp) * eam_dr
151 end do
152
153 ! Build array of rho values
154 do j = 1, eam_nrho
155 rhovals(j) = real(j-1,kind=dp) * eam_drho
156 end do
157
158 ! convert from eV to kcal / mol:
159 do j = 1, eam_nrho
160 eam_F_rho(j) = eam_F_rho(j) * 23.06054E0_DP
161 end do
162
163 ! precompute the pair potential and get it into kcal / mol:
164 eam_phi_r(1) = 0.0E0_DP
165 do j = 2, eam_nr
166 eam_phi_r(j) = 331.999296E0_DP * (eam_Z_r(j)**2) / rvals(j)
167 enddo
168
169 call newSpline(EAMList%EAMParams(current)%rho, rvals, eam_rho_r, .true.)
170 call newSpline(EAMList%EAMParams(current)%Z, rvals, eam_Z_r, .true.)
171 call newSpline(EAMList%EAMParams(current)%F, rhovals, eam_F_rho, .true.)
172 call newSpline(EAMList%EAMParams(current)%phi, rvals, eam_phi_r, .true.)
173 end subroutine newEAMtype
174
175
176 ! kills all eam types entered and sets simulation to uninitalized
177 subroutine destroyEAMtypes()
178 integer :: i
179 type(EAMType), pointer :: tempEAMType=>null()
180
181 do i = 1, EAMList%n_eam_types
182 tempEAMType => eamList%EAMParams(i)
183 call deallocate_EAMType(tempEAMType)
184 end do
185 if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
186 eamList%EAMParams => null()
187
188 eamList%n_eam_types = 0
189 eamList%currentAddition = 0
190 end subroutine destroyEAMtypes
191
192 function getEAMCut(atomID) result(cutValue)
193 integer, intent(in) :: atomID
194 integer :: eamID
195 real(kind=dp) :: cutValue
196
197 eamID = EAMList%atidToEAMType(atomID)
198 cutValue = EAMList%EAMParams(eamID)%eam_rcut
199 end function getEAMCut
200
201
202 subroutine setCutoffEAM(rcut)
203 real(kind=dp) :: rcut
204 EAM_rcut = rcut
205 end subroutine setCutoffEAM
206
207
208 subroutine deallocate_EAMType(thisEAMType)
209 type (EAMtype), pointer :: thisEAMType
210
211 call deleteSpline(thisEAMType%F)
212 call deleteSpline(thisEAMType%rho)
213 call deleteSpline(thisEAMType%phi)
214 call deleteSpline(thisEAMType%Z)
215
216 end subroutine deallocate_EAMType
217
218 !! Calculates rho_r
219 subroutine calc_eam_prepair_rho(atid1, atid2, d, r, rijsq, rho_i_at_j, rho_j_at_i)
220 integer :: Atid1, Atid2
221 real(kind = dp), dimension(3) :: d
222 real(kind = dp), intent(inout) :: r
223 real(kind = dp), intent(inout) :: rijsq
224 ! value of electron density rho do to atom i at atom j
225 real(kind = dp), intent(inout) :: rho_i_at_j
226 ! value of electron density rho do to atom j at atom i
227 real(kind = dp), intent(inout) :: rho_j_at_i
228 real(kind = dp) :: rci, rcj
229 integer :: eam_err
230
231 integer :: myid_atom1 ! EAM atid
232 integer :: myid_atom2
233
234 ! check to see if we need to be cleaned at the start of a force loop
235
236 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
237 Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
238
239 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
240
241 ! get cutoff for atom 1
242 rci = EAMList%EAMParams(myid_atom1)%eam_rcut
243 ! get type specific cutoff for atom 2
244 rcj = EAMList%EAMParams(myid_atom2)%eam_rcut
245
246 if (r.lt.rci) then
247 call lookupEAMSpline(EAMList%EAMParams(myid_atom1)%rho, r, &
248 rho_i_at_j)
249 endif
250
251 if (r.lt.rcj) then
252 call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
253 rho_j_at_i)
254 endif
255 endif
256 end subroutine calc_eam_prepair_rho
257
258
259 !! Calculate the functional F(rho) for all local atoms
260 subroutine calc_eam_preforce_Frho(nlocal, pot, particle_pot)
261 integer :: nlocal
262 real(kind=dp) :: pot
263 integer :: i, j
264 integer :: atom
265 real(kind=dp) :: U,U1
266 real( kind = dp ), dimension(nlocal) :: particle_pot
267 integer :: atype1
268 integer :: me, atid1
269
270 ! cleanme = .true.
271 !! Calculate F(rho) and derivative
272 do atom = 1, nlocal
273 atid1 = atid(atom)
274 me = eamList%atidToEAMtype(atid1)
275
276 call lookupEAMSpline1d(EAMList%EAMParams(me)%F, rho(atom), &
277 u, u1)
278
279 frho(atom) = u
280 dfrhodrho(atom) = u1
281 pot = pot + u
282 particle_pot(atom) = particle_pot(atom) + u
283
284 enddo
285
286 end subroutine calc_eam_preforce_Frho
287
288 !! Does EAM pairwise Force calculation.
289 subroutine do_eam_pair(atid1, atid2, d, rij, r2, sw, vpair, &
290 fpair, pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, &
291 fshift_i, fshift_j)
292 !Arguments
293 integer, intent(in) :: atid1, atid2
294 real( kind = dp ), intent(in) :: rij, r2
295 real( kind = dp ) :: pot, sw, vpair
296 real( kind = dp ), dimension(3) :: f1
297 real( kind = dp ), intent(in), dimension(3) :: d
298 real( kind = dp ), intent(inout), dimension(3) :: fpair
299 real( kind = dp ), intent(inout) :: dfrhodrho_i, dfrhodrho_j
300 real( kind = dp ), intent(inout) :: rho_i, rho_j
301 real( kind = dp ), intent(inout):: fshift_i, fshift_j
302
303 real( kind = dp ) :: drdx, drdy, drdz
304 real( kind = dp ) :: phab, pha, dvpdr
305 real( kind = dp ) :: rha, drha, dpha
306 real( kind = dp ) :: rhb, drhb, dphb
307 real( kind = dp ) :: dudr
308 real( kind = dp ) :: rci, rcj
309 real( kind = dp ) :: drhoidr, drhojdr
310 real( kind = dp ) :: Fx, Fy, Fz
311 real( kind = dp ) :: r, phb
312 real( kind = dp ) :: u1, u2
313
314 integer :: id1, id2
315 integer :: mytype_atom1
316 integer :: mytype_atom2
317
318 phab = 0.0E0_DP
319 dvpdr = 0.0E0_DP
320
321 if (rij .lt. EAM_rcut) then
322
323 mytype_atom1 = EAMList%atidToEAMType(atid1)
324 mytype_atom2 = EAMList%atidTOEAMType(atid2)
325
326
327 ! get cutoff for atom 1
328 rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
329 ! get type specific cutoff for atom 2
330 rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
331
332 drdx = d(1)/rij
333 drdy = d(2)/rij
334 drdz = d(3)/rij
335
336 if (rij.lt.rci) then
337
338 ! Calculate rho and drho for atom1
339
340 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
341 rij, rha, drha)
342
343 ! Calculate Phi(r) for atom1.
344
345 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
346 rij, pha, dpha)
347
348 endif
349
350 if (rij.lt.rcj) then
351
352 ! Calculate rho and drho for atom2
353
354 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
355 rij, rhb, drhb)
356
357 ! Calculate Phi(r) for atom2.
358
359 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
360 rij, phb, dphb)
361
362 endif
363
364 if (rij.lt.rci) then
365 phab = phab + 0.5E0_DP*(rhb/rha)*pha
366 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
367 pha*((drhb/rha) - (rhb*drha/rha/rha)))
368 endif
369
370 if (rij.lt.rcj) then
371 phab = phab + 0.5E0_DP*(rha/rhb)*phb
372 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
373 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
374 endif
375
376 drhoidr = drha
377 drhojdr = drhb
378
379 dudr = drhojdr*dfrhodrho_i + drhoidr*dfrhodrho_j + dvpdr
380
381
382 fx = dudr * drdx
383 fy = dudr * drdy
384 fz = dudr * drdz
385
386 ! particle_pot is the difference between the full potential
387 ! and the full potential without the presence of a particular
388 ! particle (atom1).
389 !
390 ! This reduces the density at other particle locations, so
391 ! we need to recompute the density at atom2 assuming atom1
392 ! didn't contribute. This then requires recomputing the
393 ! density functional for atom2 as well.
394 !
395 ! Most of the particle_pot heavy lifting comes from the
396 ! pair interaction, and will be handled by vpair.
397
398 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, &
399 rho_i-rhb, fshift_i, u1)
400 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, &
401 rho_j-rha, fshift_j, u2)
402
403 pot = pot + phab
404
405 f1(1) = f1(1) + fx
406 f1(2) = f1(2) + fy
407 f1(3) = f1(3) + fz
408
409 vpair = vpair + phab
410
411 endif
412 end subroutine do_eam_pair
413
414 subroutine lookupEAMSpline(cs, xval, yval)
415
416 implicit none
417
418 type (cubicSpline), intent(in) :: cs
419 real( kind = DP ), intent(in) :: xval
420 real( kind = DP ), intent(out) :: yval
421 real( kind = DP ) :: dx
422 integer :: i, j
423 !
424 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
425 ! or is nearest to xval.
426
427 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
428
429 dx = xval - cs%x(j)
430 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
431
432 return
433 end subroutine lookupEAMSpline
434
435 subroutine lookupEAMSpline1d(cs, xval, yval, dydx)
436
437 implicit none
438
439 type (cubicSpline), intent(in) :: cs
440 real( kind = DP ), intent(in) :: xval
441 real( kind = DP ), intent(out) :: yval, dydx
442 real( kind = DP ) :: dx
443 integer :: i, j
444
445 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
446 ! or is nearest to xval.
447
448
449 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
450
451 dx = xval - cs%x(j)
452 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
453
454 dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j))
455
456 return
457 end subroutine lookupEAMSpline1d
458
459 end module eam

Properties

Name Value
svn:keywords Author Id Revision Date