ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/eam.F90
Revision: 1442
Committed: Mon May 10 17:28:26 2010 UTC (15 years, 2 months ago) by gezelter
File size: 14711 byte(s)
Log Message:
Adding property set to svn entries

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(atom1, atom2, atid1, atid2, d, r, rijsq, rho_i_at_j, rho_j_at_i)
220 integer :: atom1, atom2, 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 integer :: eam_err
229
230 integer :: myid_atom1 ! EAM atid
231 integer :: myid_atom2
232
233 ! check to see if we need to be cleaned at the start of a force loop
234
235 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
236 Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
237
238 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
239
240 call lookupEAMSpline(EAMList%EAMParams(myid_atom1)%rho, r, &
241 rho_i_at_j)
242
243 call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
244 rho_j_at_i)
245 endif
246 end subroutine calc_eam_prepair_rho
247
248
249 !! Calculate the functional F(rho) for all local atoms
250 subroutine calc_eam_preforce_Frho(nlocal, pot, particle_pot)
251 integer :: nlocal
252 real(kind=dp) :: pot
253 integer :: i, j
254 integer :: atom
255 real(kind=dp) :: U,U1
256 real( kind = dp ), dimension(nlocal) :: particle_pot
257 integer :: atype1
258 integer :: me, atid1
259
260 ! cleanme = .true.
261 !! Calculate F(rho) and derivative
262 do atom = 1, nlocal
263 atid1 = atid(atom)
264 me = eamList%atidToEAMtype(atid1)
265
266 call lookupEAMSpline1d(EAMList%EAMParams(me)%F, rho(atom), &
267 u, u1)
268
269 frho(atom) = u
270 dfrhodrho(atom) = u1
271 pot = pot + u
272 particle_pot(atom) = particle_pot(atom) + u
273
274 enddo
275
276 end subroutine calc_eam_preforce_Frho
277
278 !! Does EAM pairwise Force calculation.
279 subroutine do_eam_pair(atom1, atom2, atid1, atid2, d, rij, r2, sw, vpair, &
280 fpair, pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, &
281 fshift_i, fshift_j, do_pot)
282 !Arguments
283 integer, intent(in) :: atom1, atom2, atid1, atid2
284 real( kind = dp ), intent(in) :: rij, r2
285 real( kind = dp ) :: pot, sw, vpair
286 real( kind = dp ), dimension(3) :: f1
287 real( kind = dp ), intent(in), dimension(3) :: d
288 real( kind = dp ), intent(inout), dimension(3) :: fpair
289 real( kind = dp ), intent(inout) :: dfrhodrho_i, dfrhodrho_j
290 real( kind = dp ), intent(inout) :: rho_i, rho_j
291 real( kind = dp ), intent(inout):: fshift_i, fshift_j
292
293 logical, intent(in) :: do_pot
294
295 real( kind = dp ) :: drdx, drdy, drdz
296 real( kind = dp ) :: phab, pha, dvpdr
297 real( kind = dp ) :: rha, drha, dpha
298 real( kind = dp ) :: rhb, drhb, dphb
299 real( kind = dp ) :: dudr
300 real( kind = dp ) :: rci, rcj
301 real( kind = dp ) :: drhoidr, drhojdr
302 real( kind = dp ) :: Fx, Fy, Fz
303 real( kind = dp ) :: r, phb
304 real( kind = dp ) :: u1, u2
305
306 integer :: id1, id2
307 integer :: mytype_atom1
308 integer :: mytype_atom2
309
310 phab = 0.0E0_DP
311 dvpdr = 0.0E0_DP
312
313 if (rij .lt. EAM_rcut) then
314
315 mytype_atom1 = EAMList%atidToEAMType(atid1)
316 mytype_atom2 = EAMList%atidTOEAMType(atid2)
317
318
319 ! get cutoff for atom 1
320 rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
321 ! get type specific cutoff for atom 2
322 rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
323
324 drdx = d(1)/rij
325 drdy = d(2)/rij
326 drdz = d(3)/rij
327
328 if (rij.lt.rci) then
329
330 ! Calculate rho and drho for atom1
331
332 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
333 rij, rha, drha)
334
335 ! Calculate Phi(r) for atom1.
336
337 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
338 rij, pha, dpha)
339
340 endif
341
342 if (rij.lt.rcj) then
343
344 ! Calculate rho and drho for atom2
345
346 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
347 rij, rhb, drhb)
348
349 ! Calculate Phi(r) for atom2.
350
351 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
352 rij, phb, dphb)
353
354 endif
355
356 if (rij.lt.rci) then
357 phab = phab + 0.5E0_DP*(rhb/rha)*pha
358 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
359 pha*((drhb/rha) - (rhb*drha/rha/rha)))
360 endif
361
362 if (rij.lt.rcj) then
363 phab = phab + 0.5E0_DP*(rha/rhb)*phb
364 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
365 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
366 endif
367
368 drhoidr = drha
369 drhojdr = drhb
370
371 dudr = drhojdr*dfrhodrho_i + drhoidr*dfrhodrho_j + dvpdr
372
373
374 fx = dudr * drdx
375 fy = dudr * drdy
376 fz = dudr * drdz
377
378
379 if (do_pot) then
380 ! particle_pot is the difference between the full potential
381 ! and the full potential without the presence of a particular
382 ! particle (atom1).
383 !
384 ! This reduces the density at other particle locations, so
385 ! we need to recompute the density at atom2 assuming atom1
386 ! didn't contribute. This then requires recomputing the
387 ! density functional for atom2 as well.
388 !
389 ! Most of the particle_pot heavy lifting comes from the
390 ! pair interaction, and will be handled by vpair.
391
392 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, &
393 rho_i-rhb, &
394 fshift_i, u1)
395 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, &
396 rho_j-rha, &
397 fshift_j, u2)
398
399
400 end if
401
402 pot = pot + phab
403
404 f1(1) = f1(1) + fx
405 f1(2) = f1(2) + fy
406 f1(3) = f1(3) + fz
407
408 vpair = vpair + phab
409
410 endif
411 end subroutine do_eam_pair
412
413 subroutine lookupEAMSpline(cs, xval, yval)
414
415 implicit none
416
417 type (cubicSpline), intent(in) :: cs
418 real( kind = DP ), intent(in) :: xval
419 real( kind = DP ), intent(out) :: yval
420 real( kind = DP ) :: dx
421 integer :: i, j
422 !
423 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
424 ! or is nearest to xval.
425
426 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
427
428 dx = xval - cs%x(j)
429 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
430
431 return
432 end subroutine lookupEAMSpline
433
434 subroutine lookupEAMSpline1d(cs, xval, yval, dydx)
435
436 implicit none
437
438 type (cubicSpline), intent(in) :: cs
439 real( kind = DP ), intent(in) :: xval
440 real( kind = DP ), intent(out) :: yval, dydx
441 real( kind = DP ) :: dx
442 integer :: i, j
443
444 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
445 ! or is nearest to xval.
446
447
448 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
449
450 dx = xval - cs%x(j)
451 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
452
453 dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j))
454
455 return
456 end subroutine lookupEAMSpline1d
457
458 end module eam

Properties

Name Value
svn:keywords Author Id Revision Date