ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/eam.F90
Revision: 2291
Committed: Wed Sep 14 20:31:38 2005 UTC (18 years, 9 months ago) by chuckv
File size: 32952 byte(s)
Log Message:
EAM now uses eamlist to lookup eamAtypes instead of assuming a 1-1 correspondence to atypes.

File Contents

# User Rev Content
1 gezelter 1930 !!
2     !! Copyright (c) 2005 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. Acknowledgement of the program authors must be made in any
10     !! publication of scientific results based in part on use of the
11     !! program. An acceptable form of acknowledgement is citation of
12     !! the article in which the program was described (Matthew
13     !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42 gezelter 1608 module eam
43     use simulation
44     use force_globals
45     use status
46     use atype_module
47 chrisfen 2277 use vector_class
48 gezelter 1608 #ifdef IS_MPI
49     use mpiSimulation
50     #endif
51     implicit none
52     PRIVATE
53    
54 chuckv 1624 INTEGER, PARAMETER :: DP = selected_real_kind(15)
55 gezelter 1608
56     logical, save :: EAM_FF_initialized = .false.
57     integer, save :: EAM_Mixing_Policy
58     real(kind = dp), save :: EAM_rcut
59     logical, save :: haveRcut = .false.
60    
61     character(len = statusMsgSize) :: errMesg
62     integer :: eam_err
63    
64     character(len = 200) :: errMsg
65     character(len=*), parameter :: RoutineName = "EAM MODULE"
66 gezelter 2204 !! Logical that determines if eam arrays should be zeroed
67 gezelter 1608 logical :: cleanme = .true.
68     logical :: nmflag = .false.
69    
70 gezelter 2204
71 gezelter 1608 type, private :: EAMtype
72     integer :: eam_atype
73     real( kind = DP ) :: eam_dr
74     integer :: eam_nr
75     integer :: eam_nrho
76     real( kind = DP ) :: eam_lattice
77     real( kind = DP ) :: eam_drho
78     real( kind = DP ) :: eam_rcut
79     integer :: eam_atype_map
80 gezelter 2204
81 gezelter 1608 real( kind = DP ), pointer, dimension(:) :: eam_rvals => null()
82     real( kind = DP ), pointer, dimension(:) :: eam_rhovals => null()
83     real( kind = DP ), pointer, dimension(:) :: eam_F_rho => null()
84     real( kind = DP ), pointer, dimension(:) :: eam_Z_r => null()
85     real( kind = DP ), pointer, dimension(:) :: eam_rho_r => null()
86     real( kind = DP ), pointer, dimension(:) :: eam_phi_r => null()
87     real( kind = DP ), pointer, dimension(:) :: eam_F_rho_pp => null()
88     real( kind = DP ), pointer, dimension(:) :: eam_Z_r_pp => null()
89     real( kind = DP ), pointer, dimension(:) :: eam_rho_r_pp => null()
90     real( kind = DP ), pointer, dimension(:) :: eam_phi_r_pp => null()
91     end type EAMtype
92    
93    
94     !! Arrays for derivatives used in force calculation
95     real( kind = dp), dimension(:), allocatable :: frho
96     real( kind = dp), dimension(:), allocatable :: rho
97    
98     real( kind = dp), dimension(:), allocatable :: dfrhodrho
99     real( kind = dp), dimension(:), allocatable :: d2frhodrhodrho
100    
101    
102 gezelter 2204 !! Arrays for MPI storage
103 gezelter 1608 #ifdef IS_MPI
104     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
105     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
106     real( kind = dp),save, dimension(:), allocatable :: frho_row
107     real( kind = dp),save, dimension(:), allocatable :: frho_col
108     real( kind = dp),save, dimension(:), allocatable :: rho_row
109     real( kind = dp),save, dimension(:), allocatable :: rho_col
110     real( kind = dp),save, dimension(:), allocatable :: rho_tmp
111     real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_col
112     real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_row
113     #endif
114    
115     type, private :: EAMTypeList
116     integer :: n_eam_types = 0
117     integer :: currentAddition = 0
118 gezelter 2204
119 gezelter 1608 type (EAMtype), pointer :: EAMParams(:) => null()
120 chuckv 2276 integer, pointer :: atidToEAMType(:) => null()
121 gezelter 1608 end type EAMTypeList
122    
123    
124     type (eamTypeList), save :: EAMList
125    
126     !! standard eam stuff
127    
128    
129     public :: init_EAM_FF
130     public :: setCutoffEAM
131     public :: do_eam_pair
132     public :: newEAMtype
133     public :: calc_eam_prepair_rho
134     public :: calc_eam_preforce_Frho
135     public :: clean_EAM
136 chuckv 2169 public :: destroyEAMTypes
137 chrisfen 2277 public :: getEAMCut
138 gezelter 1608
139     contains
140    
141    
142     subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
143     eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho,&
144 chrisfen 2278 c_ident,status)
145 gezelter 1608 real (kind = dp ) :: lattice_constant
146     integer :: eam_nrho
147     real (kind = dp ) :: eam_drho
148     integer :: eam_nr
149     real (kind = dp ) :: eam_dr
150     real (kind = dp ) :: rcut
151     real (kind = dp ), dimension(eam_nr) :: eam_Z_r
152     real (kind = dp ), dimension(eam_nr) :: eam_rho_r
153     real (kind = dp ), dimension(eam_nrho) :: eam_F_rho
154 chrisfen 2278 integer :: c_ident
155 gezelter 1608 integer :: status
156    
157 chuckv 2276 integer :: nAtypes,nEAMTypes,myATID
158 gezelter 1608 integer :: maxVals
159     integer :: alloc_stat
160     integer :: current
161     integer,pointer :: Matchlist(:) => null()
162    
163     status = 0
164    
165    
166     !! Assume that atypes has already been set and get the total number of types in atypes
167     !! Also assume that every member of atypes is a EAM model.
168    
169 gezelter 2204
170 gezelter 1608 ! check to see if this is the first time into
171     if (.not.associated(EAMList%EAMParams)) then
172 chuckv 2276 call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList)
173     EAMList%n_eam_types = nEAMtypes
174     allocate(EAMList%EAMParams(nEAMTypes))
175     nAtypes = getSize(atypes)
176     allocate(EAMList%atidToEAMType(nAtypes))
177 gezelter 1608 end if
178    
179     EAMList%currentAddition = EAMList%currentAddition + 1
180     current = EAMList%currentAddition
181    
182 chrisfen 2278 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
183 chuckv 2276 EAMList%atidToEAMType(myATID) = current
184 gezelter 2204
185 gezelter 1608 call allocate_EAMType(eam_nrho,eam_nr,EAMList%EAMParams(current),stat=alloc_stat)
186     if (alloc_stat /= 0) then
187     status = -1
188     return
189     end if
190    
191 chuckv 2276
192 chrisfen 2278 EAMList%EAMParams(current)%eam_atype = c_ident
193 gezelter 1608 EAMList%EAMParams(current)%eam_lattice = lattice_constant
194     EAMList%EAMParams(current)%eam_nrho = eam_nrho
195     EAMList%EAMParams(current)%eam_drho = eam_drho
196     EAMList%EAMParams(current)%eam_nr = eam_nr
197     EAMList%EAMParams(current)%eam_dr = eam_dr
198     EAMList%EAMParams(current)%eam_rcut = rcut
199     EAMList%EAMParams(current)%eam_Z_r = eam_Z_r
200     EAMList%EAMParams(current)%eam_rho_r = eam_rho_r
201     EAMList%EAMParams(current)%eam_F_rho = eam_F_rho
202    
203     end subroutine newEAMtype
204    
205    
206 chuckv 2169 ! kills all eam types entered and sets simulation to uninitalized
207     subroutine destroyEAMtypes()
208     integer :: i
209     type(EAMType), pointer :: tempEAMType=>null()
210 gezelter 1608
211 chuckv 2169 do i = 1, EAMList%n_eam_types
212     tempEAMType => eamList%EAMParams(i)
213     call deallocate_EAMType(tempEAMType)
214     end do
215     if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
216     eamList%EAMParams => null()
217 gezelter 2204
218 chuckv 2169 eamList%n_eam_types = 0
219     eamList%currentAddition = 0
220 gezelter 2204
221 chuckv 2169 end subroutine destroyEAMtypes
222    
223 chrisfen 2277 function getEAMCut(atomID) result(cutValue)
224     integer, intent(in) :: atomID
225 chrisfen 2278 integer :: eamID
226 chrisfen 2277 real(kind=dp) :: cutValue
227    
228 chrisfen 2278 eamID = EAMList%atidToEAMType(atomID)
229     cutValue = EAMList%EAMParams(eamID)%eam_rcut
230 chrisfen 2277 end function getEAMCut
231    
232 gezelter 1608 subroutine init_EAM_FF(status)
233     integer :: status
234     integer :: i,j
235     real(kind=dp) :: current_rcut_max
236     integer :: alloc_stat
237     integer :: number_r, number_rho
238    
239    
240     status = 0
241     if (EAMList%currentAddition == 0) then
242     call handleError("init_EAM_FF","No members in EAMList")
243     status = -1
244     return
245     end if
246    
247    
248 gezelter 2204 do i = 1, EAMList%currentAddition
249 gezelter 1608
250 gezelter 2204 ! Build array of r values
251 gezelter 1608
252 gezelter 2204 do j = 1,EAMList%EAMParams(i)%eam_nr
253     EAMList%EAMParams(i)%eam_rvals(j) = &
254     real(j-1,kind=dp)* &
255     EAMList%EAMParams(i)%eam_dr
256 gezelter 1608 end do
257 gezelter 2204 ! Build array of rho values
258     do j = 1,EAMList%EAMParams(i)%eam_nrho
259     EAMList%EAMParams(i)%eam_rhovals(j) = &
260     real(j-1,kind=dp)* &
261     EAMList%EAMParams(i)%eam_drho
262     end do
263     ! convert from eV to kcal / mol:
264     EAMList%EAMParams(i)%eam_F_rho = EAMList%EAMParams(i)%eam_F_rho * 23.06054E0_DP
265 gezelter 1608
266 gezelter 2204 ! precompute the pair potential and get it into kcal / mol:
267     EAMList%EAMParams(i)%eam_phi_r(1) = 0.0E0_DP
268     do j = 2, EAMList%EAMParams(i)%eam_nr
269     EAMList%EAMParams(i)%eam_phi_r(j) = (EAMList%EAMParams(i)%eam_Z_r(j)**2)/EAMList%EAMParams(i)%eam_rvals(j)
270     EAMList%EAMParams(i)%eam_phi_r(j) = EAMList%EAMParams(i)%eam_phi_r(j)*331.999296E0_DP
271 gezelter 1608 enddo
272 gezelter 2204 end do
273 gezelter 1608
274    
275 gezelter 2204 do i = 1, EAMList%currentAddition
276     number_r = EAMList%EAMParams(i)%eam_nr
277     number_rho = EAMList%EAMParams(i)%eam_nrho
278    
279     call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
280     EAMList%EAMParams(i)%eam_rho_r, &
281     EAMList%EAMParams(i)%eam_rho_r_pp, &
282     0.0E0_DP, 0.0E0_DP, 'N')
283     call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
284     EAMList%EAMParams(i)%eam_Z_r, &
285     EAMList%EAMParams(i)%eam_Z_r_pp, &
286     0.0E0_DP, 0.0E0_DP, 'N')
287     call eam_spline(number_rho, EAMList%EAMParams(i)%eam_rhovals, &
288     EAMList%EAMParams(i)%eam_F_rho, &
289     EAMList%EAMParams(i)%eam_F_rho_pp, &
290     0.0E0_DP, 0.0E0_DP, 'N')
291     call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
292     EAMList%EAMParams(i)%eam_phi_r, &
293     EAMList%EAMParams(i)%eam_phi_r_pp, &
294     0.0E0_DP, 0.0E0_DP, 'N')
295     enddo
296    
297     ! current_rcut_max = EAMList%EAMParams(1)%eam_rcut
298     !! find the smallest rcut for any eam atype
299     ! do i = 2, EAMList%currentAddition
300     ! current_rcut_max =max(current_rcut_max,EAMList%EAMParams(i)%eam_rcut)
301     ! end do
302    
303     ! EAM_rcut = current_rcut_max
304     ! EAM_rcut_orig = current_rcut_max
305     ! do i = 1, EAMList%currentAddition
306     ! EAMList%EAMParam(i)s%eam_atype_map(eam_atype(i)) = i
307     ! end do
308     !! Allocate arrays for force calculation
309    
310     call allocateEAM(alloc_stat)
311     if (alloc_stat /= 0 ) then
312     write(*,*) "allocateEAM failed"
313     status = -1
314     return
315     endif
316    
317 gezelter 1608 end subroutine init_EAM_FF
318    
319 gezelter 2204 !! routine checks to see if array is allocated, deallocates array if allocated
320     !! and then creates the array to the required size
321 gezelter 1608 subroutine allocateEAM(status)
322     integer, intent(out) :: status
323    
324     #ifdef IS_MPI
325     integer :: nAtomsInRow
326     integer :: nAtomsInCol
327     #endif
328     integer :: alloc_stat
329    
330    
331     status = 0
332     #ifdef IS_MPI
333     nAtomsInRow = getNatomsInRow(plan_atom_row)
334     nAtomsInCol = getNatomsInCol(plan_atom_col)
335     #endif
336    
337     if (allocated(frho)) deallocate(frho)
338     allocate(frho(nlocal),stat=alloc_stat)
339     if (alloc_stat /= 0) then
340     status = -1
341     return
342     end if
343     if (allocated(rho)) deallocate(rho)
344     allocate(rho(nlocal),stat=alloc_stat)
345     if (alloc_stat /= 0) then
346     status = -1
347     return
348     end if
349    
350     if (allocated(dfrhodrho)) deallocate(dfrhodrho)
351     allocate(dfrhodrho(nlocal),stat=alloc_stat)
352     if (alloc_stat /= 0) then
353     status = -1
354     return
355     end if
356    
357     if (allocated(d2frhodrhodrho)) deallocate(d2frhodrhodrho)
358     allocate(d2frhodrhodrho(nlocal),stat=alloc_stat)
359     if (alloc_stat /= 0) then
360     status = -1
361     return
362     end if
363 gezelter 2204
364 gezelter 1608 #ifdef IS_MPI
365    
366     if (allocated(rho_tmp)) deallocate(rho_tmp)
367     allocate(rho_tmp(nlocal),stat=alloc_stat)
368     if (alloc_stat /= 0) then
369     status = -1
370     return
371     end if
372    
373    
374     if (allocated(frho_row)) deallocate(frho_row)
375     allocate(frho_row(nAtomsInRow),stat=alloc_stat)
376     if (alloc_stat /= 0) then
377     status = -1
378     return
379     end if
380     if (allocated(rho_row)) deallocate(rho_row)
381     allocate(rho_row(nAtomsInRow),stat=alloc_stat)
382     if (alloc_stat /= 0) then
383     status = -1
384     return
385     end if
386     if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
387     allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
388     if (alloc_stat /= 0) then
389     status = -1
390     return
391     end if
392     if (allocated(d2frhodrhodrho_row)) deallocate(d2frhodrhodrho_row)
393     allocate(d2frhodrhodrho_row(nAtomsInRow),stat=alloc_stat)
394     if (alloc_stat /= 0) then
395     status = -1
396     return
397     end if
398    
399    
400 gezelter 2204 ! Now do column arrays
401 gezelter 1608
402     if (allocated(frho_col)) deallocate(frho_col)
403     allocate(frho_col(nAtomsInCol),stat=alloc_stat)
404     if (alloc_stat /= 0) then
405     status = -1
406     return
407     end if
408     if (allocated(rho_col)) deallocate(rho_col)
409     allocate(rho_col(nAtomsInCol),stat=alloc_stat)
410     if (alloc_stat /= 0) then
411     status = -1
412     return
413     end if
414     if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
415     allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
416     if (alloc_stat /= 0) then
417     status = -1
418     return
419     end if
420     if (allocated(d2frhodrhodrho_col)) deallocate(d2frhodrhodrho_col)
421     allocate(d2frhodrhodrho_col(nAtomsInCol),stat=alloc_stat)
422     if (alloc_stat /= 0) then
423     status = -1
424     return
425     end if
426 gezelter 2204
427 gezelter 1608 #endif
428    
429     end subroutine allocateEAM
430    
431 gezelter 2204 !! C sets rcut to be the largest cutoff of any atype
432     !! present in this simulation. Doesn't include all atypes
433     !! sim knows about, just those in the simulation.
434 gezelter 1608 subroutine setCutoffEAM(rcut, status)
435     real(kind=dp) :: rcut
436     integer :: status
437     status = 0
438    
439     EAM_rcut = rcut
440    
441     end subroutine setCutoffEAM
442    
443    
444    
445     subroutine clean_EAM()
446 gezelter 2204
447     ! clean non-IS_MPI first
448 gezelter 1608 frho = 0.0_dp
449     rho = 0.0_dp
450     dfrhodrho = 0.0_dp
451 gezelter 2204 ! clean MPI if needed
452 gezelter 1608 #ifdef IS_MPI
453     frho_row = 0.0_dp
454     frho_col = 0.0_dp
455     rho_row = 0.0_dp
456     rho_col = 0.0_dp
457     rho_tmp = 0.0_dp
458     dfrhodrho_row = 0.0_dp
459     dfrhodrho_col = 0.0_dp
460     #endif
461     end subroutine clean_EAM
462    
463    
464    
465     subroutine allocate_EAMType(eam_n_rho,eam_n_r,thisEAMType,stat)
466     integer, intent(in) :: eam_n_rho
467     integer, intent(in) :: eam_n_r
468     type (EAMType) :: thisEAMType
469     integer, optional :: stat
470     integer :: alloc_stat
471    
472    
473    
474     if (present(stat)) stat = 0
475 gezelter 2204
476 gezelter 1608 allocate(thisEAMType%eam_rvals(eam_n_r),stat=alloc_stat)
477     if (alloc_stat /= 0 ) then
478     if (present(stat)) stat = -1
479     return
480     end if
481     allocate(thisEAMType%eam_rhovals(eam_n_rho),stat=alloc_stat)
482     if (alloc_stat /= 0 ) then
483     if (present(stat)) stat = -1
484     return
485     end if
486     allocate(thisEAMType%eam_F_rho(eam_n_rho),stat=alloc_stat)
487     if (alloc_stat /= 0 ) then
488     if (present(stat)) stat = -1
489     return
490     end if
491     allocate(thisEAMType%eam_Z_r(eam_n_r),stat=alloc_stat)
492     if (alloc_stat /= 0 ) then
493     if (present(stat)) stat = -1
494     return
495     end if
496     allocate(thisEAMType%eam_rho_r(eam_n_r),stat=alloc_stat)
497     if (alloc_stat /= 0 ) then
498     if (present(stat)) stat = -1
499     return
500     end if
501     allocate(thisEAMType%eam_phi_r(eam_n_r),stat=alloc_stat)
502     if (alloc_stat /= 0 ) then
503     if (present(stat)) stat = -1
504     return
505     end if
506     allocate(thisEAMType%eam_F_rho_pp(eam_n_rho),stat=alloc_stat)
507     if (alloc_stat /= 0 ) then
508     if (present(stat)) stat = -1
509     return
510     end if
511     allocate(thisEAMType%eam_Z_r_pp(eam_n_r),stat=alloc_stat)
512     if (alloc_stat /= 0 ) then
513     if (present(stat)) stat = -1
514     return
515     end if
516     allocate(thisEAMType%eam_rho_r_pp(eam_n_r),stat=alloc_stat)
517     if (alloc_stat /= 0 ) then
518     if (present(stat)) stat = -1
519     return
520     end if
521     allocate(thisEAMType%eam_phi_r_pp(eam_n_r),stat=alloc_stat)
522     if (alloc_stat /= 0 ) then
523     if (present(stat)) stat = -1
524     return
525     end if
526    
527 gezelter 2204
528 gezelter 1608 end subroutine allocate_EAMType
529    
530    
531     subroutine deallocate_EAMType(thisEAMType)
532     type (EAMtype), pointer :: thisEAMType
533    
534     ! free Arrays in reverse order of allocation...
535 chuckv 2169 if(associated(thisEAMType%eam_phi_r_pp)) deallocate(thisEAMType%eam_phi_r_pp)
536     if(associated(thisEAMType%eam_rho_r_pp)) deallocate(thisEAMType%eam_rho_r_pp)
537     if(associated(thisEAMType%eam_Z_r_pp)) deallocate(thisEAMType%eam_Z_r_pp)
538     if(associated(thisEAMType%eam_F_rho_pp)) deallocate(thisEAMType%eam_F_rho_pp)
539     if(associated(thisEAMType%eam_phi_r)) deallocate(thisEAMType%eam_phi_r)
540     if(associated(thisEAMType%eam_rho_r)) deallocate(thisEAMType%eam_rho_r)
541     if(associated(thisEAMType%eam_Z_r)) deallocate(thisEAMType%eam_Z_r)
542     if(associated(thisEAMType%eam_F_rho)) deallocate(thisEAMType%eam_F_rho)
543     if(associated(thisEAMType%eam_rhovals)) deallocate(thisEAMType%eam_rhovals)
544     if(associated(thisEAMType%eam_rvals)) deallocate(thisEAMType%eam_rvals)
545 gezelter 2204
546 gezelter 1608 end subroutine deallocate_EAMType
547    
548 gezelter 2204 !! Calculates rho_r
549 gezelter 1608 subroutine calc_eam_prepair_rho(atom1,atom2,d,r,rijsq)
550     integer :: atom1,atom2
551     real(kind = dp), dimension(3) :: d
552     real(kind = dp), intent(inout) :: r
553     real(kind = dp), intent(inout) :: rijsq
554     ! value of electron density rho do to atom i at atom j
555     real(kind = dp) :: rho_i_at_j
556     ! value of electron density rho do to atom j at atom i
557     real(kind = dp) :: rho_j_at_i
558    
559     ! we don't use the derivatives, dummy variables
560     real( kind = dp) :: drho,d2rho
561     integer :: eam_err
562 chuckv 2291
563     integer :: atid1,atid2 ! Global atid
564     integer :: myid_atom1 ! EAM atid
565     integer :: myid_atom2
566 gezelter 2204
567 gezelter 1608
568 gezelter 2204 ! check to see if we need to be cleaned at the start of a force loop
569 gezelter 1608
570    
571    
572 gezelter 2204
573 gezelter 1608 #ifdef IS_MPI
574 chuckv 2291 Atid1 = Atid_row(Atom1)
575     Atid2 = Atid_col(Atom2)
576 gezelter 1608 #else
577 chuckv 2291 Atid1 = Atid(Atom1)
578     Atid2 = Atid(Atom2)
579 gezelter 1608 #endif
580    
581 chuckv 2291 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
582     Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
583    
584 gezelter 1608 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
585    
586    
587    
588     call eam_splint(EAMList%EAMParams(myid_atom1)%eam_nr, &
589     EAMList%EAMParams(myid_atom1)%eam_rvals, &
590     EAMList%EAMParams(myid_atom1)%eam_rho_r, &
591     EAMList%EAMParams(myid_atom1)%eam_rho_r_pp, &
592     r, rho_i_at_j,drho,d2rho)
593    
594    
595 gezelter 2204
596 gezelter 1608 #ifdef IS_MPI
597     rho_col(atom2) = rho_col(atom2) + rho_i_at_j
598     #else
599     rho(atom2) = rho(atom2) + rho_i_at_j
600     #endif
601 chuckv 2291 ! write(*,*) atom1,atom2,r,rho_i_at_j
602 gezelter 2204 endif
603 gezelter 1608
604 gezelter 2204 if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then
605     call eam_splint(EAMList%EAMParams(myid_atom2)%eam_nr, &
606     EAMList%EAMParams(myid_atom2)%eam_rvals, &
607     EAMList%EAMParams(myid_atom2)%eam_rho_r, &
608     EAMList%EAMParams(myid_atom2)%eam_rho_r_pp, &
609     r, rho_j_at_i,drho,d2rho)
610 gezelter 1608
611 gezelter 2204
612    
613    
614 gezelter 1608 #ifdef IS_MPI
615 gezelter 2204 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
616 gezelter 1608 #else
617 gezelter 2204 rho(atom1) = rho(atom1) + rho_j_at_i
618 gezelter 1608 #endif
619 gezelter 2204 endif
620 gezelter 1608
621    
622    
623    
624    
625    
626     end subroutine calc_eam_prepair_rho
627    
628    
629    
630    
631     !! Calculate the functional F(rho) for all local atoms
632     subroutine calc_eam_preforce_Frho(nlocal,pot)
633     integer :: nlocal
634     real(kind=dp) :: pot
635     integer :: i,j
636     integer :: atom
637     real(kind=dp) :: U,U1,U2
638     integer :: atype1
639 chuckv 2291 integer :: me,atid1
640 gezelter 1608 integer :: n_rho_points
641    
642 gezelter 2204
643 gezelter 1608 cleanme = .true.
644     !! Scatter the electron density from pre-pair calculation back to local atoms
645     #ifdef IS_MPI
646     call scatter(rho_row,rho,plan_atom_row,eam_err)
647     if (eam_err /= 0 ) then
648 gezelter 2204 write(errMsg,*) " Error scattering rho_row into rho"
649     call handleError(RoutineName,errMesg)
650     endif
651 gezelter 1608 call scatter(rho_col,rho_tmp,plan_atom_col,eam_err)
652     if (eam_err /= 0 ) then
653 gezelter 2204 write(errMsg,*) " Error scattering rho_col into rho"
654     call handleError(RoutineName,errMesg)
655     endif
656 gezelter 1608
657 gezelter 2204 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
658 gezelter 1608 #endif
659    
660    
661    
662 gezelter 2204 !! Calculate F(rho) and derivative
663 gezelter 1608 do atom = 1, nlocal
664 chuckv 2291 atid1 = atid(atom)
665     me = eamList%atidToEAMtype(atid1)
666 gezelter 1608 n_rho_points = EAMList%EAMParams(me)%eam_nrho
667     ! Check to see that the density is not greater than the larges rho we have calculated
668     if (rho(atom) < EAMList%EAMParams(me)%eam_rhovals(n_rho_points)) then
669     call eam_splint(n_rho_points, &
670     EAMList%EAMParams(me)%eam_rhovals, &
671     EAMList%EAMParams(me)%eam_f_rho, &
672     EAMList%EAMParams(me)%eam_f_rho_pp, &
673     rho(atom), & ! Actual Rho
674     u, u1, u2)
675     else
676     ! Calculate F(rho with the largest available rho value
677     call eam_splint(n_rho_points, &
678     EAMList%EAMParams(me)%eam_rhovals, &
679     EAMList%EAMParams(me)%eam_f_rho, &
680     EAMList%EAMParams(me)%eam_f_rho_pp, &
681     EAMList%EAMParams(me)%eam_rhovals(n_rho_points), & ! Largest rho
682     u,u1,u2)
683     end if
684    
685    
686     frho(atom) = u
687     dfrhodrho(atom) = u1
688     d2frhodrhodrho(atom) = u2
689     pot = pot + u
690    
691     enddo
692    
693    
694 gezelter 2204
695 gezelter 1608 #ifdef IS_MPI
696     !! communicate f(rho) and derivatives back into row and column arrays
697     call gather(frho,frho_row,plan_atom_row, eam_err)
698     if (eam_err /= 0) then
699     call handleError("cal_eam_forces()","MPI gather frho_row failure")
700     endif
701     call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, eam_err)
702     if (eam_err /= 0) then
703     call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
704     endif
705     call gather(frho,frho_col,plan_atom_col, eam_err)
706     if (eam_err /= 0) then
707     call handleError("cal_eam_forces()","MPI gather frho_col failure")
708     endif
709     call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, eam_err)
710     if (eam_err /= 0) then
711     call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
712     endif
713    
714    
715    
716    
717    
718     if (nmflag) then
719     call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_atom_row)
720     call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_atom_col)
721     endif
722     #endif
723    
724 gezelter 2204
725 gezelter 1608 end subroutine calc_eam_preforce_Frho
726    
727    
728    
729    
730     !! Does EAM pairwise Force calculation.
731     subroutine do_eam_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
732     pot, f, do_pot)
733     !Arguments
734     integer, intent(in) :: atom1, atom2
735     real( kind = dp ), intent(in) :: rij, r2
736     real( kind = dp ) :: pot, sw, vpair
737     real( kind = dp ), dimension(3,nLocal) :: f
738     real( kind = dp ), intent(in), dimension(3) :: d
739     real( kind = dp ), intent(inout), dimension(3) :: fpair
740    
741     logical, intent(in) :: do_pot
742 gezelter 2204
743 gezelter 1608 real( kind = dp ) :: drdx,drdy,drdz
744     real( kind = dp ) :: d2
745     real( kind = dp ) :: phab,pha,dvpdr,d2vpdrdr
746     real( kind = dp ) :: rha,drha,d2rha, dpha
747     real( kind = dp ) :: rhb,drhb,d2rhb, dphb
748     real( kind = dp ) :: dudr
749     real( kind = dp ) :: rci,rcj
750     real( kind = dp ) :: drhoidr,drhojdr
751     real( kind = dp ) :: d2rhoidrdr
752     real( kind = dp ) :: d2rhojdrdr
753     real( kind = dp ) :: Fx,Fy,Fz
754     real( kind = dp ) :: r,d2pha,phb,d2phb
755    
756     integer :: id1,id2
757     integer :: mytype_atom1
758     integer :: mytype_atom2
759 chuckv 2291 integer :: atid1,atid2
760 gezelter 2204 !Local Variables
761 gezelter 1608
762 gezelter 2204 ! write(*,*) "Frho: ", Frho(atom1)
763    
764 gezelter 1608 phab = 0.0E0_DP
765     dvpdr = 0.0E0_DP
766     d2vpdrdr = 0.0E0_DP
767    
768     if (rij .lt. EAM_rcut) then
769    
770     #ifdef IS_MPI
771 chuckv 2291 atid1 = atid_row(atom1)
772     atid2 = atid_col(atom2)
773 gezelter 1608 #else
774 chuckv 2291 atid1 = atid(atom1)
775     atid2 = atid(atom2)
776 gezelter 1608 #endif
777 chuckv 2291
778     mytype_atom1 = EAMList%atidToEAMType(atid1)
779     mytype_atom2 = EAMList%atidTOEAMType(atid2)
780    
781    
782 gezelter 1608 ! get cutoff for atom 1
783     rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
784     ! get type specific cutoff for atom 2
785     rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
786 gezelter 2204
787 gezelter 1608 drdx = d(1)/rij
788     drdy = d(2)/rij
789     drdz = d(3)/rij
790 gezelter 2204
791 gezelter 1608 if (rij.lt.rci) then
792     call eam_splint(EAMList%EAMParams(mytype_atom1)%eam_nr, &
793     EAMList%EAMParams(mytype_atom1)%eam_rvals, &
794     EAMList%EAMParams(mytype_atom1)%eam_rho_r, &
795     EAMList%EAMParams(mytype_atom1)%eam_rho_r_pp, &
796     rij, rha,drha,d2rha)
797     !! Calculate Phi(r) for atom1.
798     call eam_splint(EAMList%EAMParams(mytype_atom1)%eam_nr, &
799     EAMList%EAMParams(mytype_atom1)%eam_rvals, &
800     EAMList%EAMParams(mytype_atom1)%eam_phi_r, &
801     EAMList%EAMParams(mytype_atom1)%eam_phi_r_pp, &
802     rij, pha,dpha,d2pha)
803     endif
804    
805     if (rij.lt.rcj) then
806     ! Calculate rho,drho and d2rho for atom1
807     call eam_splint(EAMList%EAMParams(mytype_atom2)%eam_nr, &
808     EAMList%EAMParams(mytype_atom2)%eam_rvals, &
809     EAMList%EAMParams(mytype_atom2)%eam_rho_r, &
810     EAMList%EAMParams(mytype_atom2)%eam_rho_r_pp, &
811     rij, rhb,drhb,d2rhb)
812 gezelter 2204
813 gezelter 1608 !! Calculate Phi(r) for atom2.
814     call eam_splint(EAMList%EAMParams(mytype_atom2)%eam_nr, &
815     EAMList%EAMParams(mytype_atom2)%eam_rvals, &
816     EAMList%EAMParams(mytype_atom2)%eam_phi_r, &
817     EAMList%EAMParams(mytype_atom2)%eam_phi_r_pp, &
818     rij, phb,dphb,d2phb)
819     endif
820    
821     if (rij.lt.rci) then
822     phab = phab + 0.5E0_DP*(rhb/rha)*pha
823     dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
824     pha*((drhb/rha) - (rhb*drha/rha/rha)))
825     d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rhb/rha)*d2pha + &
826     2.0E0_DP*dpha*((drhb/rha) - (rhb*drha/rha/rha)) + &
827     pha*((d2rhb/rha) - 2.0E0_DP*(drhb*drha/rha/rha) + &
828     (2.0E0_DP*rhb*drha*drha/rha/rha/rha) - (rhb*d2rha/rha/rha)))
829     endif
830 gezelter 2204
831 gezelter 1608 if (rij.lt.rcj) then
832     phab = phab + 0.5E0_DP*(rha/rhb)*phb
833     dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
834     phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
835     d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rha/rhb)*d2phb + &
836     2.0E0_DP*dphb*((drha/rhb) - (rha*drhb/rhb/rhb)) + &
837     phb*((d2rha/rhb) - 2.0E0_DP*(drha*drhb/rhb/rhb) + &
838     (2.0E0_DP*rha*drhb*drhb/rhb/rhb/rhb) - (rha*d2rhb/rhb/rhb)))
839     endif
840 gezelter 2204
841 gezelter 1608 drhoidr = drha
842     drhojdr = drhb
843    
844     d2rhoidrdr = d2rha
845     d2rhojdrdr = d2rhb
846    
847    
848     #ifdef IS_MPI
849     dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) &
850     + dvpdr
851    
852     #else
853     dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) &
854     + dvpdr
855 gezelter 2204 ! write(*,*) "Atom1,Atom2, dfrhodrho(atom1) dfrhodrho(atom2): ", atom1,atom2,dfrhodrho(atom1),dfrhodrho(atom2)
856 gezelter 1608 #endif
857    
858     fx = dudr * drdx
859     fy = dudr * drdy
860     fz = dudr * drdz
861    
862    
863     #ifdef IS_MPI
864     if (do_pot) then
865     pot_Row(atom1) = pot_Row(atom1) + phab*0.5
866     pot_Col(atom2) = pot_Col(atom2) + phab*0.5
867     end if
868    
869     f_Row(1,atom1) = f_Row(1,atom1) + fx
870     f_Row(2,atom1) = f_Row(2,atom1) + fy
871     f_Row(3,atom1) = f_Row(3,atom1) + fz
872 gezelter 2204
873 gezelter 1608 f_Col(1,atom2) = f_Col(1,atom2) - fx
874     f_Col(2,atom2) = f_Col(2,atom2) - fy
875     f_Col(3,atom2) = f_Col(3,atom2) - fz
876     #else
877    
878     if(do_pot) then
879     pot = pot + phab
880     end if
881    
882     f(1,atom1) = f(1,atom1) + fx
883     f(2,atom1) = f(2,atom1) + fy
884     f(3,atom1) = f(3,atom1) + fz
885 gezelter 2204
886 gezelter 1608 f(1,atom2) = f(1,atom2) - fx
887     f(2,atom2) = f(2,atom2) - fy
888     f(3,atom2) = f(3,atom2) - fz
889     #endif
890 gezelter 2204
891 gezelter 1608 vpair = vpair + phab
892     #ifdef IS_MPI
893     id1 = AtomRowToGlobal(atom1)
894     id2 = AtomColToGlobal(atom2)
895     #else
896     id1 = atom1
897     id2 = atom2
898     #endif
899 gezelter 2204
900 gezelter 1608 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
901 gezelter 2204
902 gezelter 1608 fpair(1) = fpair(1) + fx
903     fpair(2) = fpair(2) + fy
904     fpair(3) = fpair(3) + fz
905 gezelter 2204
906 gezelter 1608 endif
907 gezelter 2204
908 gezelter 1608 if (nmflag) then
909    
910     drhoidr = drha
911     drhojdr = drhb
912     d2rhoidrdr = d2rha
913     d2rhojdrdr = d2rhb
914    
915     #ifdef IS_MPI
916     d2 = d2vpdrdr + &
917     d2rhoidrdr*dfrhodrho_col(atom2) + &
918     d2rhojdrdr*dfrhodrho_row(atom1) + &
919     drhoidr*drhoidr*d2frhodrhodrho_col(atom2) + &
920     drhojdr*drhojdr*d2frhodrhodrho_row(atom1)
921 gezelter 2204
922 gezelter 1608 #else
923    
924     d2 = d2vpdrdr + &
925     d2rhoidrdr*dfrhodrho(atom2) + &
926     d2rhojdrdr*dfrhodrho(atom1) + &
927     drhoidr*drhoidr*d2frhodrhodrho(atom2) + &
928     drhojdr*drhojdr*d2frhodrhodrho(atom1)
929     #endif
930     end if
931 gezelter 2204
932     endif
933 gezelter 1608 end subroutine do_eam_pair
934    
935    
936     subroutine eam_splint(nx, xa, ya, yppa, x, y, dy, d2y)
937    
938     integer :: atype, nx, j
939     real( kind = DP ), dimension(:) :: xa
940     real( kind = DP ), dimension(:) :: ya
941     real( kind = DP ), dimension(:) :: yppa
942     real( kind = DP ) :: x, y
943     real( kind = DP ) :: dy, d2y
944     real( kind = DP ) :: del, h, a, b, c, d
945     integer :: pp_arraySize
946    
947 gezelter 2204
948 gezelter 1608 ! this spline code assumes that the x points are equally spaced
949     ! do not attempt to use this code if they are not.
950 gezelter 2204
951    
952 gezelter 1608 ! find the closest point with a value below our own:
953     j = FLOOR(real((nx-1),kind=dp) * (x - xa(1)) / (xa(nx) - xa(1))) + 1
954    
955     ! check to make sure we're inside the spline range:
956     if ((j.gt.nx).or.(j.lt.1)) then
957     write(errMSG,*) "EAM_splint: x is outside bounds of spline: ",x,j
958     call handleError(routineName,errMSG)
959     endif
960     ! check to make sure we haven't screwed up the calculation of j:
961     if ((x.lt.xa(j)).or.(x.gt.xa(j+1))) then
962     if (j.ne.nx) then
963 gezelter 2204 write(errMSG,*) "EAM_splint:",x," x is outside bounding range"
964     call handleError(routineName,errMSG)
965 gezelter 1608 endif
966     endif
967    
968     del = xa(j+1) - x
969     h = xa(j+1) - xa(j)
970 gezelter 2204
971 gezelter 1608 a = del / h
972     b = 1.0E0_DP - a
973     c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP
974     d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP
975 gezelter 2204
976 gezelter 1608 y = a*ya(j) + b*ya(j+1) + c*yppa(j) + d*yppa(j+1)
977    
978 gezelter 2204 dy = (ya(j+1)-ya(j))/h &
979     - (3.0E0_DP*a*a - 1.0E0_DP)*h*yppa(j)/6.0E0_DP &
980     + (3.0E0_DP*b*b - 1.0E0_DP)*h*yppa(j+1)/6.0E0_DP
981    
982    
983     d2y = a*yppa(j) + b*yppa(j+1)
984    
985    
986 gezelter 1608 end subroutine eam_splint
987    
988    
989     subroutine eam_spline(nx, xa, ya, yppa, yp1, ypn, boundary)
990    
991    
992     ! yp1 and ypn are the first derivatives of y at the two endpoints
993     ! if boundary is 'L' the lower derivative is used
994     ! if boundary is 'U' the upper derivative is used
995     ! if boundary is 'B' then both derivatives are used
996     ! if boundary is anything else, then both derivatives are assumed to be 0
997 gezelter 2204
998 gezelter 1608 integer :: nx, i, k, max_array_size
999 gezelter 2204
1000 gezelter 1608 real( kind = DP ), dimension(:) :: xa
1001     real( kind = DP ), dimension(:) :: ya
1002     real( kind = DP ), dimension(:) :: yppa
1003     real( kind = DP ), dimension(size(xa)) :: u
1004     real( kind = DP ) :: yp1,ypn,un,qn,sig,p
1005     character(len=*) :: boundary
1006 gezelter 2204
1007 gezelter 1608 ! make sure the sizes match
1008     if ((nx /= size(xa)) .or. (nx /= size(ya))) then
1009     call handleWarning("EAM_SPLINE","Array size mismatch")
1010     end if
1011    
1012     if ((boundary.eq.'l').or.(boundary.eq.'L').or. &
1013     (boundary.eq.'b').or.(boundary.eq.'B')) then
1014     yppa(1) = -0.5E0_DP
1015     u(1) = (3.0E0_DP/(xa(2)-xa(1)))*((ya(2)-&
1016     ya(1))/(xa(2)-xa(1))-yp1)
1017     else
1018     yppa(1) = 0.0E0_DP
1019     u(1) = 0.0E0_DP
1020     endif
1021 gezelter 2204
1022 gezelter 1608 do i = 2, nx - 1
1023     sig = (xa(i) - xa(i-1)) / (xa(i+1) - xa(i-1))
1024     p = sig * yppa(i-1) + 2.0E0_DP
1025     yppa(i) = (sig - 1.0E0_DP) / p
1026     u(i) = (6.0E0_DP*((ya(i+1)-ya(i))/(xa(i+1)-xa(i)) - &
1027     (ya(i)-ya(i-1))/(xa(i)-xa(i-1)))/ &
1028     (xa(i+1)-xa(i-1)) - sig * u(i-1))/p
1029     enddo
1030 gezelter 2204
1031 gezelter 1608 if ((boundary.eq.'u').or.(boundary.eq.'U').or. &
1032     (boundary.eq.'b').or.(boundary.eq.'B')) then
1033     qn = 0.5E0_DP
1034     un = (3.0E0_DP/(xa(nx)-xa(nx-1)))* &
1035     (ypn-(ya(nx)-ya(nx-1))/(xa(nx)-xa(nx-1)))
1036     else
1037     qn = 0.0E0_DP
1038     un = 0.0E0_DP
1039     endif
1040    
1041     yppa(nx)=(un-qn*u(nx-1))/(qn*yppa(nx-1)+1.0E0_DP)
1042 gezelter 2204
1043 gezelter 1608 do k = nx-1, 1, -1
1044     yppa(k)=yppa(k)*yppa(k+1)+u(k)
1045     enddo
1046    
1047     end subroutine eam_spline
1048    
1049     end module eam