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 |
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(atid1, atid2, d, rij, r2, sw, vpair, & |
280 |
pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, & |
281 |
fshift_i, fshift_j) |
282 |
!Arguments |
283 |
integer, intent(in) :: 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) :: dfrhodrho_i, dfrhodrho_j |
289 |
real( kind = dp ), intent(inout) :: rho_i, rho_j |
290 |
real( kind = dp ), intent(inout):: fshift_i, fshift_j |
291 |
|
292 |
real( kind = dp ) :: drdx, drdy, drdz |
293 |
real( kind = dp ) :: phab, pha, dvpdr |
294 |
real( kind = dp ) :: rha, drha, dpha |
295 |
real( kind = dp ) :: rhb, drhb, dphb |
296 |
real( kind = dp ) :: dudr |
297 |
real( kind = dp ) :: rci, rcj |
298 |
real( kind = dp ) :: drhoidr, drhojdr |
299 |
real( kind = dp ) :: Fx, Fy, Fz |
300 |
real( kind = dp ) :: r, phb |
301 |
real( kind = dp ) :: u1, u2 |
302 |
|
303 |
integer :: id1, id2 |
304 |
integer :: mytype_atom1 |
305 |
integer :: mytype_atom2 |
306 |
|
307 |
phab = 0.0E0_DP |
308 |
dvpdr = 0.0E0_DP |
309 |
|
310 |
if (rij .lt. EAM_rcut) then |
311 |
|
312 |
mytype_atom1 = EAMList%atidToEAMType(atid1) |
313 |
mytype_atom2 = EAMList%atidTOEAMType(atid2) |
314 |
|
315 |
|
316 |
! get cutoff for atom 1 |
317 |
rci = EAMList%EAMParams(mytype_atom1)%eam_rcut |
318 |
! get type specific cutoff for atom 2 |
319 |
rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut |
320 |
|
321 |
drdx = d(1)/rij |
322 |
drdy = d(2)/rij |
323 |
drdz = d(3)/rij |
324 |
|
325 |
if (rij.lt.rci) then |
326 |
|
327 |
! Calculate rho and drho for atom1 |
328 |
|
329 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, & |
330 |
rij, rha, drha) |
331 |
|
332 |
! Calculate Phi(r) for atom1. |
333 |
|
334 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, & |
335 |
rij, pha, dpha) |
336 |
|
337 |
endif |
338 |
|
339 |
if (rij.lt.rcj) then |
340 |
|
341 |
! Calculate rho and drho for atom2 |
342 |
|
343 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, & |
344 |
rij, rhb, drhb) |
345 |
|
346 |
! Calculate Phi(r) for atom2. |
347 |
|
348 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, & |
349 |
rij, phb, dphb) |
350 |
|
351 |
endif |
352 |
|
353 |
if (rij.lt.rci) then |
354 |
phab = phab + 0.5E0_DP*(rhb/rha)*pha |
355 |
dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + & |
356 |
pha*((drhb/rha) - (rhb*drha/rha/rha))) |
357 |
endif |
358 |
|
359 |
if (rij.lt.rcj) then |
360 |
phab = phab + 0.5E0_DP*(rha/rhb)*phb |
361 |
dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + & |
362 |
phb*((drha/rhb) - (rha*drhb/rhb/rhb))) |
363 |
endif |
364 |
|
365 |
drhoidr = drha |
366 |
drhojdr = drhb |
367 |
|
368 |
dudr = drhojdr*dfrhodrho_i + drhoidr*dfrhodrho_j + dvpdr |
369 |
|
370 |
|
371 |
fx = dudr * drdx |
372 |
fy = dudr * drdy |
373 |
fz = dudr * drdz |
374 |
|
375 |
! particle_pot is the difference between the full potential |
376 |
! and the full potential without the presence of a particular |
377 |
! particle (atom1). |
378 |
! |
379 |
! This reduces the density at other particle locations, so |
380 |
! we need to recompute the density at atom2 assuming atom1 |
381 |
! didn't contribute. This then requires recomputing the |
382 |
! density functional for atom2 as well. |
383 |
! |
384 |
! Most of the particle_pot heavy lifting comes from the |
385 |
! pair interaction, and will be handled by vpair. |
386 |
|
387 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, & |
388 |
rho_i-rhb, fshift_i, u1) |
389 |
call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, & |
390 |
rho_j-rha, fshift_j, u2) |
391 |
|
392 |
pot = pot + phab |
393 |
|
394 |
f1(1) = f1(1) + fx |
395 |
f1(2) = f1(2) + fy |
396 |
f1(3) = f1(3) + fz |
397 |
|
398 |
vpair = vpair + phab |
399 |
|
400 |
endif |
401 |
end subroutine do_eam_pair |
402 |
|
403 |
subroutine lookupEAMSpline(cs, xval, yval) |
404 |
|
405 |
implicit none |
406 |
|
407 |
type (cubicSpline), intent(in) :: cs |
408 |
real( kind = DP ), intent(in) :: xval |
409 |
real( kind = DP ), intent(out) :: yval |
410 |
real( kind = DP ) :: dx |
411 |
integer :: i, j |
412 |
! |
413 |
! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains |
414 |
! or is nearest to xval. |
415 |
|
416 |
j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1)) |
417 |
|
418 |
dx = xval - cs%x(j) |
419 |
yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) |
420 |
|
421 |
return |
422 |
end subroutine lookupEAMSpline |
423 |
|
424 |
subroutine lookupEAMSpline1d(cs, xval, yval, dydx) |
425 |
|
426 |
implicit none |
427 |
|
428 |
type (cubicSpline), intent(in) :: cs |
429 |
real( kind = DP ), intent(in) :: xval |
430 |
real( kind = DP ), intent(out) :: yval, dydx |
431 |
real( kind = DP ) :: dx |
432 |
integer :: i, j |
433 |
|
434 |
! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains |
435 |
! or is nearest to xval. |
436 |
|
437 |
|
438 |
j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1)) |
439 |
|
440 |
dx = xval - cs%x(j) |
441 |
yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) |
442 |
|
443 |
dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j)) |
444 |
|
445 |
return |
446 |
end subroutine lookupEAMSpline1d |
447 |
|
448 |
end module eam |