ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/eam.F90
Revision: 2355
Committed: Wed Oct 12 18:59:16 2005 UTC (18 years, 10 months ago) by chuckv
File size: 33054 byte(s)
Log Message:
Breaky Breaky: Add Support for seperating potential contributions.

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