ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_eam.F90
Revision: 801
Committed: Sat Oct 4 18:46:12 2003 UTC (20 years, 9 months ago) by chuckv
File size: 29488 byte(s)
Log Message:
Fixed bug in calc_eam.

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