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