ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/eam.F90
Revision: 1578
Committed: Thu Jun 9 17:42:11 2011 UTC (14 years, 1 month ago) by gezelter
File size: 15118 byte(s)
Log Message:
be careful with non-EAM types in density calc

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 EAMList%atidToEAMType = -1
137 end if
138
139 EAMList%currentAddition = EAMList%currentAddition + 1
140 current = EAMList%currentAddition
141
142 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
143 EAMList%atidToEAMType(myATID) = current
144
145 EAMList%EAMParams(current)%eam_atype = c_ident
146 EAMList%EAMParams(current)%eam_lattice = lattice_constant
147 EAMList%EAMParams(current)%eam_rcut = rcut
148
149 ! Build array of r values
150 do j = 1, eam_nr
151 rvals(j) = real(j-1,kind=dp) * eam_dr
152 end do
153
154 ! Build array of rho values
155 do j = 1, eam_nrho
156 rhovals(j) = real(j-1,kind=dp) * eam_drho
157 end do
158
159 ! convert from eV to kcal / mol:
160 do j = 1, eam_nrho
161 eam_F_rho(j) = eam_F_rho(j) * 23.06054E0_DP
162 end do
163
164 ! precompute the pair potential and get it into kcal / mol:
165 eam_phi_r(1) = 0.0E0_DP
166 do j = 2, eam_nr
167 eam_phi_r(j) = 331.999296E0_DP * (eam_Z_r(j)**2) / rvals(j)
168 enddo
169
170 call newSpline(EAMList%EAMParams(current)%rho, rvals, eam_rho_r, .true.)
171 call newSpline(EAMList%EAMParams(current)%Z, rvals, eam_Z_r, .true.)
172 call newSpline(EAMList%EAMParams(current)%F, rhovals, eam_F_rho, .true.)
173 call newSpline(EAMList%EAMParams(current)%phi, rvals, eam_phi_r, .true.)
174 end subroutine newEAMtype
175
176
177 ! kills all eam types entered and sets simulation to uninitalized
178 subroutine destroyEAMtypes()
179 integer :: i
180 type(EAMType), pointer :: tempEAMType=>null()
181
182 do i = 1, EAMList%n_eam_types
183 tempEAMType => eamList%EAMParams(i)
184 call deallocate_EAMType(tempEAMType)
185 end do
186 if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
187 eamList%EAMParams => null()
188
189 eamList%n_eam_types = 0
190 eamList%currentAddition = 0
191 end subroutine destroyEAMtypes
192
193 function getEAMCut(atomID) result(cutValue)
194 integer, intent(in) :: atomID
195 integer :: eamID
196 real(kind=dp) :: cutValue
197
198 eamID = EAMList%atidToEAMType(atomID)
199 cutValue = EAMList%EAMParams(eamID)%eam_rcut
200 end function getEAMCut
201
202
203 subroutine setCutoffEAM(rcut)
204 real(kind=dp) :: rcut
205 EAM_rcut = rcut
206 end subroutine setCutoffEAM
207
208
209 subroutine deallocate_EAMType(thisEAMType)
210 type (EAMtype), pointer :: thisEAMType
211
212 call deleteSpline(thisEAMType%F)
213 call deleteSpline(thisEAMType%rho)
214 call deleteSpline(thisEAMType%phi)
215 call deleteSpline(thisEAMType%Z)
216
217 end subroutine deallocate_EAMType
218
219 !! Calculates rho_r
220 subroutine calc_eam_prepair_rho(atid1, atid2, d, r, rijsq, rho_i_at_j, rho_j_at_i)
221 integer :: Atid1, Atid2
222 real(kind = dp), dimension(3) :: d
223 real(kind = dp), intent(inout) :: r
224 real(kind = dp), intent(inout) :: rijsq
225 ! value of electron density rho do to atom i at atom j
226 real(kind = dp), intent(inout) :: rho_i_at_j
227 ! value of electron density rho do to atom j at atom i
228 real(kind = dp), intent(inout) :: rho_j_at_i
229 real(kind = dp) :: rci, rcj
230 integer :: eam_err
231
232 integer :: myid_atom1 ! EAM atid
233 integer :: myid_atom2
234
235 ! check to see if we need to be cleaned at the start of a force loop
236
237 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
238 Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
239
240 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
241
242 ! get cutoff for atom 1
243 rci = EAMList%EAMParams(myid_atom1)%eam_rcut
244 ! get type specific cutoff for atom 2
245 rcj = EAMList%EAMParams(myid_atom2)%eam_rcut
246
247 if (r.lt.rci) then
248 call lookupEAMSpline(EAMList%EAMParams(myid_atom1)%rho, r, &
249 rho_i_at_j)
250 endif
251
252 if (r.lt.rcj) then
253 call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
254 rho_j_at_i)
255 endif
256 endif
257 end subroutine calc_eam_prepair_rho
258
259
260 !! Calculate the functional F(rho) for all local atoms
261 subroutine calc_eam_preforce_Frho(nlocal, pot, particle_pot)
262 integer :: nlocal
263 real(kind=dp) :: pot
264 integer :: i, j
265 integer :: atom
266 real(kind=dp) :: U,U1
267 real( kind = dp ), dimension(nlocal) :: particle_pot
268 integer :: atype1
269 integer :: me, atid1
270
271 ! cleanme = .true.
272 !! Calculate F(rho) and derivative
273 do atom = 1, nlocal
274 atid1 = atid(atom)
275 me = EAMList%atidtoEAMtype(atid1)
276 ! me is set to -1 for non EAM atoms.
277 ! Punt if we are a non-EAM atom type.
278 if (me == -1) then
279 frho(atom) = 0.0_dp
280 dfrhodrho(atom) = 0.0_dp
281 else
282
283 call lookupEAMSpline1d(EAMList%EAMParams(me)%F, rho(atom), &
284 u, u1)
285
286 frho(atom) = u
287 dfrhodrho(atom) = u1
288 endif
289 pot = pot + u
290 particle_pot(atom) = particle_pot(atom) + u
291
292 enddo
293
294 end subroutine calc_eam_preforce_Frho
295
296 !! Does EAM pairwise Force calculation.
297 subroutine do_eam_pair(atid1, atid2, d, rij, r2, sw, vpair, &
298 fpair, pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, &
299 fshift_i, fshift_j)
300 !Arguments
301 integer, intent(in) :: atid1, atid2
302 real( kind = dp ), intent(in) :: rij, r2
303 real( kind = dp ) :: pot, sw, vpair
304 real( kind = dp ), dimension(3) :: f1
305 real( kind = dp ), intent(in), dimension(3) :: d
306 real( kind = dp ), intent(inout), dimension(3) :: fpair
307 real( kind = dp ), intent(inout) :: dfrhodrho_i, dfrhodrho_j
308 real( kind = dp ), intent(inout) :: rho_i, rho_j
309 real( kind = dp ), intent(inout):: fshift_i, fshift_j
310
311 real( kind = dp ) :: drdx, drdy, drdz
312 real( kind = dp ) :: phab, pha, dvpdr
313 real( kind = dp ) :: rha, drha, dpha
314 real( kind = dp ) :: rhb, drhb, dphb
315 real( kind = dp ) :: dudr
316 real( kind = dp ) :: rci, rcj
317 real( kind = dp ) :: drhoidr, drhojdr
318 real( kind = dp ) :: Fx, Fy, Fz
319 real( kind = dp ) :: r, phb
320 real( kind = dp ) :: u1, u2
321
322 integer :: id1, id2
323 integer :: mytype_atom1
324 integer :: mytype_atom2
325
326 phab = 0.0E0_DP
327 dvpdr = 0.0E0_DP
328 drha = 0.0
329 drhb = 0.0
330
331 if (rij .lt. EAM_rcut) then
332
333 mytype_atom1 = EAMList%atidToEAMType(atid1)
334 mytype_atom2 = EAMList%atidTOEAMType(atid2)
335
336
337 ! get cutoff for atom 1
338 rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
339 ! get type specific cutoff for atom 2
340 rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
341
342 drdx = d(1)/rij
343 drdy = d(2)/rij
344 drdz = d(3)/rij
345
346 if (rij.lt.rci) then
347
348 ! Calculate rho and drho for atom1
349
350 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
351 rij, rha, drha)
352
353 ! Calculate Phi(r) for atom1.
354
355 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
356 rij, pha, dpha)
357
358 endif
359
360 if (rij.lt.rcj) then
361
362 ! Calculate rho and drho for atom2
363
364 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
365 rij, rhb, drhb)
366
367 ! Calculate Phi(r) for atom2.
368
369 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
370 rij, phb, dphb)
371
372 endif
373
374 if (rij.lt.rci) then
375 phab = phab + 0.5E0_DP*(rhb/rha)*pha
376 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
377 pha*((drhb/rha) - (rhb*drha/rha/rha)))
378 endif
379
380 if (rij.lt.rcj) then
381 phab = phab + 0.5E0_DP*(rha/rhb)*phb
382 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
383 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
384 endif
385
386 drhoidr = drha
387 drhojdr = drhb
388
389 dudr = drhojdr*dfrhodrho_i + drhoidr*dfrhodrho_j + dvpdr
390
391 fx = dudr * drdx
392 fy = dudr * drdy
393 fz = dudr * drdz
394
395 ! particle_pot is the difference between the full potential
396 ! and the full potential without the presence of a particular
397 ! particle (atom1).
398 !
399 ! This reduces the density at other particle locations, so
400 ! we need to recompute the density at atom2 assuming atom1
401 ! didn't contribute. This then requires recomputing the
402 ! density functional for atom2 as well.
403 !
404 ! Most of the particle_pot heavy lifting comes from the
405 ! pair interaction, and will be handled by vpair.
406
407 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, &
408 rho_i-rhb, fshift_i, u1)
409 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, &
410 rho_j-rha, fshift_j, u2)
411
412 pot = pot + phab
413
414 f1(1) = f1(1) + fx
415 f1(2) = f1(2) + fy
416 f1(3) = f1(3) + fz
417
418 vpair = vpair + phab
419
420 endif
421 end subroutine do_eam_pair
422
423 subroutine lookupEAMSpline(cs, xval, yval)
424
425 implicit none
426
427 type (cubicSpline), intent(in) :: cs
428 real( kind = DP ), intent(in) :: xval
429 real( kind = DP ), intent(out) :: yval
430 real( kind = DP ) :: dx
431 integer :: i, j
432 !
433 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
434 ! or is nearest to xval.
435
436 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
437
438 dx = xval - cs%x(j)
439 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
440
441 return
442 end subroutine lookupEAMSpline
443
444 subroutine lookupEAMSpline1d(cs, xval, yval, dydx)
445
446 implicit none
447
448 type (cubicSpline), intent(in) :: cs
449 real( kind = DP ), intent(in) :: xval
450 real( kind = DP ), intent(out) :: yval, dydx
451 real( kind = DP ) :: dx
452 integer :: i, j
453
454 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
455 ! or is nearest to xval.
456
457
458 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
459
460 dx = xval - cs%x(j)
461 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
462
463 dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j))
464
465 return
466 end subroutine lookupEAMSpline1d
467
468 end module eam

Properties

Name Value
svn:keywords Author Id Revision Date