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