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, 3 months ago) by gezelter
File size: 15118 byte(s)
Log Message:
be careful with non-EAM types in density calc

File Contents

# User Rev Content
1 gezelter 246 !!
2 chuckv 1388 !! Copyright (c) 2005, 2009 The University of Notre Dame. All Rights Reserved.
3 gezelter 246 !!
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 gezelter 1390 !! 1. Redistributions of source code must retain the above copyright
10 gezelter 246 !! notice, this list of conditions and the following disclaimer.
11     !!
12 gezelter 1390 !! 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 !! 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 gezelter 1390 !! 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 gezelter 115 module eam
42 gezelter 960 use definitions
43 gezelter 115 use simulation
44     use force_globals
45     use status
46     use atype_module
47 chrisfen 578 use vector_class
48 gezelter 938 use interpolation
49 gezelter 115 implicit none
50     PRIVATE
51 chuckv 656 #define __FORTRAN90
52     #include "UseTheForce/DarkSide/fInteractionMap.h"
53 gezelter 115
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 gezelter 507 !! Logical that determines if eam arrays should be zeroed
65 gezelter 115 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 gezelter 938 type(cubicSpline) :: rho
73     type(cubicSpline) :: Z
74     type(cubicSpline) :: F
75     type(cubicSpline) :: phi
76 gezelter 115 end type EAMtype
77    
78    
79     type, private :: EAMTypeList
80     integer :: n_eam_types = 0
81     integer :: currentAddition = 0
82 gezelter 507
83 gezelter 115 type (EAMtype), pointer :: EAMParams(:) => null()
84 chuckv 577 integer, pointer :: atidToEAMType(:) => null()
85 gezelter 115 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 chuckv 472 public :: destroyEAMTypes
98 chrisfen 578 public :: getEAMCut
99 chrisfen 943 public :: lookupEAMSpline
100     public :: lookupEAMSpline1d
101 gezelter 115
102     contains
103    
104     subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
105 gezelter 938 eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho, c_ident, status)
106 gezelter 115 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 gezelter 938 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 chrisfen 579 integer :: c_ident
116 gezelter 115 integer :: status
117    
118 chuckv 577 integer :: nAtypes,nEAMTypes,myATID
119 gezelter 115 integer :: maxVals
120     integer :: alloc_stat
121 gezelter 938 integer :: current, j
122 gezelter 115 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 chuckv 577 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 gezelter 1578 EAMList%atidToEAMType = -1
137 gezelter 115 end if
138    
139     EAMList%currentAddition = EAMList%currentAddition + 1
140     current = EAMList%currentAddition
141    
142 chrisfen 579 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
143 chuckv 577 EAMList%atidToEAMType(myATID) = current
144 gezelter 507
145 chrisfen 579 EAMList%EAMParams(current)%eam_atype = c_ident
146 gezelter 115 EAMList%EAMParams(current)%eam_lattice = lattice_constant
147     EAMList%EAMParams(current)%eam_rcut = rcut
148    
149 gezelter 938 ! 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 gezelter 939 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 gezelter 115 end subroutine newEAMtype
175    
176    
177 chuckv 472 ! kills all eam types entered and sets simulation to uninitalized
178     subroutine destroyEAMtypes()
179     integer :: i
180     type(EAMType), pointer :: tempEAMType=>null()
181 gezelter 115
182 chuckv 472 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 gezelter 507
189 chuckv 472 eamList%n_eam_types = 0
190     eamList%currentAddition = 0
191     end subroutine destroyEAMtypes
192    
193 chrisfen 578 function getEAMCut(atomID) result(cutValue)
194     integer, intent(in) :: atomID
195 chrisfen 579 integer :: eamID
196 chrisfen 578 real(kind=dp) :: cutValue
197    
198 chrisfen 579 eamID = EAMList%atidToEAMType(atomID)
199     cutValue = EAMList%EAMParams(eamID)%eam_rcut
200 chrisfen 578 end function getEAMCut
201    
202 gezelter 115
203 gezelter 938 subroutine setCutoffEAM(rcut)
204 gezelter 115 real(kind=dp) :: rcut
205     EAM_rcut = rcut
206     end subroutine setCutoffEAM
207    
208 gezelter 507
209 gezelter 115 subroutine deallocate_EAMType(thisEAMType)
210     type (EAMtype), pointer :: thisEAMType
211    
212 gezelter 938 call deleteSpline(thisEAMType%F)
213     call deleteSpline(thisEAMType%rho)
214     call deleteSpline(thisEAMType%phi)
215     call deleteSpline(thisEAMType%Z)
216 gezelter 507
217 gezelter 115 end subroutine deallocate_EAMType
218    
219 gezelter 507 !! Calculates rho_r
220 gezelter 1464 subroutine calc_eam_prepair_rho(atid1, atid2, d, r, rijsq, rho_i_at_j, rho_j_at_i)
221     integer :: Atid1, Atid2
222 gezelter 115 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 chuckv 1388 real(kind = dp), intent(inout) :: rho_i_at_j
227 gezelter 115 ! value of electron density rho do to atom j at atom i
228 chuckv 1388 real(kind = dp), intent(inout) :: rho_j_at_i
229 gezelter 1573 real(kind = dp) :: rci, rcj
230 gezelter 115 integer :: eam_err
231 gezelter 1386
232 chuckv 592 integer :: myid_atom1 ! EAM atid
233     integer :: myid_atom2
234 gezelter 507
235     ! check to see if we need to be cleaned at the start of a force loop
236 gezelter 115
237 chuckv 592 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
238     Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
239 gezelter 1386
240 gezelter 115 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
241 gezelter 1386
242 gezelter 1573 ! 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 gezelter 115
252 gezelter 1573 if (r.lt.rcj) then
253     call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
254     rho_j_at_i)
255     endif
256 gezelter 507 endif
257 gezelter 115 end subroutine calc_eam_prepair_rho
258    
259    
260     !! Calculate the functional F(rho) for all local atoms
261 gezelter 1285 subroutine calc_eam_preforce_Frho(nlocal, pot, particle_pot)
262 gezelter 115 integer :: nlocal
263     real(kind=dp) :: pot
264 gezelter 938 integer :: i, j
265 gezelter 115 integer :: atom
266 gezelter 938 real(kind=dp) :: U,U1
267 gezelter 1285 real( kind = dp ), dimension(nlocal) :: particle_pot
268 gezelter 115 integer :: atype1
269 gezelter 938 integer :: me, atid1
270 gezelter 115
271 chuckv 1388 ! cleanme = .true.
272 gezelter 507 !! Calculate F(rho) and derivative
273 gezelter 115 do atom = 1, nlocal
274 chuckv 592 atid1 = atid(atom)
275 gezelter 1578 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 gezelter 938
286 gezelter 1578 frho(atom) = u
287     dfrhodrho(atom) = u1
288     endif
289 gezelter 115 pot = pot + u
290 gezelter 1285 particle_pot(atom) = particle_pot(atom) + u
291 gezelter 939
292 gezelter 115 enddo
293    
294     end subroutine calc_eam_preforce_Frho
295 gezelter 938
296 gezelter 115 !! Does EAM pairwise Force calculation.
297 gezelter 1464 subroutine do_eam_pair(atid1, atid2, d, rij, r2, sw, vpair, &
298 gezelter 1419 fpair, pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, &
299 gezelter 1464 fshift_i, fshift_j)
300 gezelter 115 !Arguments
301 gezelter 1464 integer, intent(in) :: atid1, atid2
302 gezelter 115 real( kind = dp ), intent(in) :: rij, r2
303     real( kind = dp ) :: pot, sw, vpair
304 gezelter 1386 real( kind = dp ), dimension(3) :: f1
305 gezelter 115 real( kind = dp ), intent(in), dimension(3) :: d
306     real( kind = dp ), intent(inout), dimension(3) :: fpair
307 chuckv 1388 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 gezelter 115
311 gezelter 938 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 gezelter 115 real( kind = dp ) :: dudr
316 gezelter 938 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 chuckv 1388 real( kind = dp ) :: u1, u2
321 gezelter 115
322 gezelter 938 integer :: id1, id2
323 gezelter 115 integer :: mytype_atom1
324     integer :: mytype_atom2
325    
326     phab = 0.0E0_DP
327     dvpdr = 0.0E0_DP
328 gezelter 1574 drha = 0.0
329     drhb = 0.0
330 gezelter 115
331     if (rij .lt. EAM_rcut) then
332    
333 chuckv 592 mytype_atom1 = EAMList%atidToEAMType(atid1)
334     mytype_atom2 = EAMList%atidTOEAMType(atid2)
335    
336    
337 gezelter 115 ! 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 gezelter 507
342 gezelter 115 drdx = d(1)/rij
343     drdy = d(2)/rij
344     drdz = d(3)/rij
345 gezelter 507
346 gezelter 115 if (rij.lt.rci) then
347 gezelter 938
348     ! Calculate rho and drho for atom1
349    
350 chrisfen 943 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
351 gezelter 938 rij, rha, drha)
352    
353     ! Calculate Phi(r) for atom1.
354    
355 chrisfen 943 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
356 gezelter 938 rij, pha, dpha)
357    
358 gezelter 115 endif
359    
360     if (rij.lt.rcj) then
361 gezelter 507
362 gezelter 938 ! Calculate rho and drho for atom2
363    
364 chrisfen 943 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
365 gezelter 938 rij, rhb, drhb)
366    
367     ! Calculate Phi(r) for atom2.
368    
369 chrisfen 943 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
370 gezelter 938 rij, phb, dphb)
371    
372 gezelter 115 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 gezelter 507
380 gezelter 115 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 gezelter 507
386 gezelter 115 drhoidr = drha
387     drhojdr = drhb
388    
389 gezelter 1574 dudr = drhojdr*dfrhodrho_i + drhoidr*dfrhodrho_j + dvpdr
390 gezelter 115
391     fx = dudr * drdx
392     fy = dudr * drdy
393     fz = dudr * drdz
394    
395 gezelter 1464 ! 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 gezelter 115
412 gezelter 1386 pot = pot + phab
413 gezelter 938
414 gezelter 1386 f1(1) = f1(1) + fx
415     f1(2) = f1(2) + fy
416     f1(3) = f1(3) + fz
417 gezelter 507
418 gezelter 1386 vpair = vpair + phab
419 gezelter 507
420     endif
421 gezelter 115 end subroutine do_eam_pair
422    
423 chrisfen 943 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 gezelter 960 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
437 chrisfen 943
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 gezelter 960 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
459 chrisfen 943
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 gezelter 115 end module eam

Properties

Name Value
svn:keywords Author Id Revision Date