ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_eam.F90
Revision: 882
Committed: Wed Dec 17 20:13:33 2003 UTC (20 years, 6 months ago) by chuckv
File size: 29812 byte(s)
Log Message:
Fixed bug in parallel EAM. rho_row and rho_col were scattered into the same array. Unfortunately, MPI zeros the array between scatters so half of the sum was being lost. Fixed by added a temp array for column scatter, then sum loop over nlocal.

File Contents

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