1 |
chuckv |
492 |
module calc_eam |
2 |
|
|
use definitions, ONLY : DP |
3 |
|
|
use force_globals |
4 |
|
|
#ifdef MPI |
5 |
|
|
use mpiSimulation |
6 |
|
|
#endif |
7 |
chuckv |
497 |
PRIVATE |
8 |
chuckv |
492 |
|
9 |
|
|
|
10 |
chuckv |
497 |
logical, save :: EAM_FF_initialized = .false. |
11 |
|
|
integer, save :: EAM_Mixing_Policy |
12 |
|
|
integer, save :: EAM_rcut |
13 |
chuckv |
492 |
|
14 |
|
|
|
15 |
chuckv |
497 |
|
16 |
|
|
|
17 |
chuckv |
492 |
!! standard eam stuff |
18 |
|
|
integer :: n_eam_atypes |
19 |
|
|
integer, allocatable, dimension(:) :: eam_atype |
20 |
|
|
real( kind = DP ), allocatable, dimension(:) :: eam_dr |
21 |
|
|
integer, allocatable, dimension(:) :: eam_nr |
22 |
|
|
integer, allocatable, dimension(:) :: eam_nrho |
23 |
|
|
real( kind = DP ), allocatable, dimension(:) :: eam_lattice |
24 |
|
|
real( kind = DP ), allocatable, dimension(:) :: eam_drho |
25 |
|
|
integer , allocatable, dimension(:) :: eam_atype_map |
26 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_rvals |
27 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_rhovals |
28 |
|
|
real( kind = DP ), allocatable, dimension(:) :: eam_rcut |
29 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_F_rho |
30 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_Z_r |
31 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_rho_r |
32 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_phi_r |
33 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_F_rho_pp |
34 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_Z_r_pp |
35 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_rho_r_pp |
36 |
|
|
real( kind = DP ), allocatable, dimension(:,:) :: eam_phi_r_pp |
37 |
|
|
|
38 |
|
|
|
39 |
chuckv |
497 |
public :: init_EAM_FF |
40 |
|
|
public :: EAM_new_rcut |
41 |
|
|
public :: do_EAM_pair |
42 |
chuckv |
492 |
|
43 |
|
|
|
44 |
|
|
|
45 |
|
|
contains |
46 |
chuckv |
497 |
subroutine init_EAM_FF() |
47 |
chuckv |
492 |
|
48 |
|
|
|
49 |
|
|
|
50 |
|
|
character(len=80) :: eam_pot_file |
51 |
|
|
integer :: i, j, max_size, prev_max_size |
52 |
|
|
integer :: number_rho, number_r |
53 |
|
|
integer :: eam_unit |
54 |
|
|
integer :: this_error |
55 |
|
|
character(len=300) :: msg |
56 |
|
|
integer, external :: nfiles |
57 |
|
|
!for mpi |
58 |
|
|
|
59 |
|
|
|
60 |
|
|
#ifdef MPI |
61 |
|
|
if (node == 0) & |
62 |
|
|
n_eam_atypes = nfiles(trim(eam_pot_dir)//char(0)) |
63 |
|
|
|
64 |
|
|
call mpi_bcast(n_eam_atypes,1,mpi_integer,0,mpi_comm_world,mpi_err) |
65 |
|
|
if (n_eam_atypes == -1) then |
66 |
|
|
call error("INITIALIZE_EAM","NO EAM potentials found!") |
67 |
|
|
endif |
68 |
|
|
write(msg,'(a5,i4,a12,i5,a14)') 'Node: ',node,' reading ...', & |
69 |
|
|
n_eam_atypes, ' eam atom types' |
70 |
|
|
call info('INITIALIZE_EAM', trim(msg)) |
71 |
|
|
#else |
72 |
|
|
n_eam_atypes = nfiles(trim(eam_pot_dir)//char(0)) |
73 |
|
|
if (n_eam_atypes == -1) then |
74 |
|
|
call error("INITIALIZE_EAM","NO EAM potentials found!") |
75 |
|
|
endif |
76 |
|
|
|
77 |
|
|
write(msg,'(a12,i5,a14)') ' Reading ...', & |
78 |
|
|
n_eam_atypes, ' eam atom types' |
79 |
|
|
call info('INITIALIZE_EAM', trim(msg)) |
80 |
|
|
#endif |
81 |
|
|
|
82 |
|
|
|
83 |
|
|
call allocate_eam_atype(n_eam_atypes) |
84 |
|
|
|
85 |
|
|
|
86 |
|
|
|
87 |
|
|
!! get largest number of data points for any potential |
88 |
|
|
#ifdef MPI |
89 |
|
|
if (node == 0) then |
90 |
|
|
#endif |
91 |
|
|
prev_max_size = 0 |
92 |
|
|
do i = 1, n_eam_atypes |
93 |
|
|
call getfilename(i, eam_pot_file) |
94 |
|
|
max_size = max(get_eam_sizes( & |
95 |
|
|
trim(eam_pot_dir) // '/' // eam_pot_file), & |
96 |
|
|
prev_max_size) |
97 |
|
|
prev_max_size = max_size |
98 |
|
|
end do |
99 |
|
|
#ifdef MPI |
100 |
|
|
end if |
101 |
|
|
|
102 |
|
|
|
103 |
|
|
call mpi_bcast(max_size,1,mpi_integer,0,mpi_comm_world,mpi_err) |
104 |
|
|
#endif |
105 |
|
|
|
106 |
|
|
call allocate_eam_module(n_eam_atypes,max_size) |
107 |
|
|
allocate(eam_atype_map(get_max_atype())) |
108 |
|
|
|
109 |
|
|
#ifdef MPI |
110 |
|
|
if (node == 0) then |
111 |
|
|
#endif |
112 |
|
|
do i = 1, n_eam_atypes |
113 |
|
|
call getfilename(i, eam_pot_file) |
114 |
|
|
call read_eam_pot(i,trim(eam_pot_dir) // '/' // eam_pot_file, & |
115 |
|
|
this_error) |
116 |
|
|
|
117 |
|
|
do j = 1, eam_nr(i) |
118 |
|
|
eam_rvals(j,i) = dble(j-1)*eam_dr(i) |
119 |
|
|
enddo |
120 |
|
|
|
121 |
|
|
do j = 1, eam_nrho(i) |
122 |
|
|
eam_rhovals(j,i) = dble(j-1)*eam_drho(i) |
123 |
|
|
enddo |
124 |
|
|
|
125 |
|
|
! convert from eV to kcal / mol: |
126 |
|
|
do j = 1, eam_nrho(i) |
127 |
|
|
eam_F_rho(j,i) = eam_F_rho(j,i)*23.06054E0_DP |
128 |
|
|
enddo |
129 |
|
|
|
130 |
|
|
! precompute the pair potential and get it into kcal / mol: |
131 |
|
|
eam_phi_r(1,i) = 0.0E0_DP |
132 |
|
|
do j = 2, eam_nr(i) |
133 |
|
|
eam_phi_r(j,i) = (eam_Z_r(j,i)**2)/eam_rvals(j,i) |
134 |
|
|
eam_phi_r(j,i) = eam_phi_r(j,i)*331.999296E0_DP |
135 |
|
|
enddo |
136 |
|
|
|
137 |
|
|
end do |
138 |
|
|
#ifdef MPI |
139 |
|
|
call info('INITIALIZE_EAM','NODE 0: Distributing spline arrays') |
140 |
|
|
endif |
141 |
|
|
|
142 |
|
|
call mpi_bcast(this_error,n_eam_atypes,mpi_integer,0, & |
143 |
|
|
mpi_comm_world,mpi_err) |
144 |
|
|
if (this_error /= 0) then |
145 |
|
|
call error('INITIALIZE_EAM',"Cannot read eam files") |
146 |
|
|
endif |
147 |
|
|
|
148 |
|
|
call mpi_bcast(eam_atype,n_eam_atypes,mpi_integer,0, & |
149 |
|
|
mpi_comm_world,mpi_err) |
150 |
|
|
|
151 |
|
|
!! distribute values to cluster...... |
152 |
|
|
call mpi_bcast(eam_nr,n_eam_atypes,mpi_integer,& |
153 |
|
|
0,mpi_comm_world,mpi_err) |
154 |
|
|
call mpi_bcast(eam_nrho,n_eam_atypes,mpi_integer,& |
155 |
|
|
0,mpi_comm_world,mpi_err) |
156 |
|
|
call mpi_bcast(eam_rvals,n_eam_atypes*max_size,mpi_double_precision, & |
157 |
|
|
0,mpi_comm_world,mpi_err) |
158 |
|
|
call mpi_bcast(eam_rcut,n_eam_atypes,mpi_double_precision, & |
159 |
|
|
0,mpi_comm_world,mpi_err) |
160 |
|
|
call mpi_bcast(eam_rhovals,n_eam_atypes*max_size,mpi_double_precision, & |
161 |
|
|
0,mpi_comm_world,mpi_err) |
162 |
|
|
|
163 |
|
|
!! distribute arrays |
164 |
|
|
call mpi_bcast(eam_rho_r,n_eam_atypes*max_size,mpi_double_precision, & |
165 |
|
|
0,mpi_comm_world,mpi_err) |
166 |
|
|
call mpi_bcast(eam_Z_r,n_eam_atypes*max_size,mpi_double_precision, & |
167 |
|
|
0,mpi_comm_world,mpi_err) |
168 |
|
|
call mpi_bcast(eam_F_rho,n_eam_atypes*max_size,mpi_double_precision, & |
169 |
|
|
0,mpi_comm_world,mpi_err) |
170 |
|
|
call mpi_bcast(eam_phi_r,n_eam_atypes*max_size,mpi_double_precision, & |
171 |
|
|
0,mpi_comm_world,mpi_err) |
172 |
|
|
|
173 |
|
|
#endif |
174 |
|
|
call info('INITIALIZE_EAM', 'creating splines') |
175 |
|
|
|
176 |
|
|
do i = 1, n_eam_atypes |
177 |
|
|
number_r = eam_nr(i) |
178 |
|
|
number_rho = eam_nrho(i) |
179 |
|
|
|
180 |
|
|
call eam_spline(i, number_r, eam_rvals, eam_rho_r, eam_rho_r_pp, & |
181 |
|
|
0.0E0_DP, 0.0E0_DP, 'N') |
182 |
|
|
call eam_spline(i, number_r, eam_rvals, eam_Z_r, eam_Z_r_pp, & |
183 |
|
|
0.0E0_DP, 0.0E0_DP, 'N') |
184 |
|
|
call eam_spline(i, number_rho, eam_rhovals, eam_F_rho, eam_F_rho_pp, & |
185 |
|
|
0.0E0_DP, 0.0E0_DP, 'N') |
186 |
|
|
call eam_spline(i, number_r, eam_rvals, eam_phi_r, eam_phi_r_pp, & |
187 |
|
|
0.0E0_DP, 0.0E0_DP, 'N') |
188 |
|
|
enddo |
189 |
|
|
|
190 |
|
|
do i = 1, n_eam_atypes |
191 |
|
|
eam_atype_map(eam_atype(i)) = i |
192 |
|
|
end do |
193 |
|
|
|
194 |
|
|
|
195 |
|
|
|
196 |
|
|
call info('INITIALIZE_EAM','Done creating splines') |
197 |
|
|
|
198 |
|
|
return |
199 |
|
|
end subroutine initialize_eam |
200 |
|
|
|
201 |
|
|
|
202 |
|
|
subroutine allocate_eam_atype(n_size_atype) |
203 |
|
|
integer, intent(in) :: n_size_atype |
204 |
|
|
|
205 |
|
|
allocate(eam_atype(n_size_atype)) |
206 |
|
|
allocate(eam_drho(n_size_atype)) |
207 |
|
|
allocate(eam_dr(n_size_atype)) |
208 |
|
|
allocate(eam_nr(n_size_atype)) |
209 |
|
|
allocate(eam_nrho(n_size_atype)) |
210 |
|
|
allocate(eam_lattice(n_size_atype)) |
211 |
|
|
allocate(eam_rcut(n_size_atype)) |
212 |
|
|
|
213 |
|
|
end subroutine allocate_eam_atype |
214 |
|
|
|
215 |
|
|
subroutine allocate_eam_module(n_size_atype,n_eam_points) |
216 |
|
|
integer, intent(in) :: n_eam_points |
217 |
|
|
integer, intent(in) :: n_size_atype |
218 |
|
|
|
219 |
|
|
allocate(eam_rvals(n_eam_points,n_size_atype)) |
220 |
|
|
allocate(eam_rhovals(n_eam_points,n_size_atype)) |
221 |
|
|
allocate(eam_F_rho(n_eam_points,n_size_atype)) |
222 |
|
|
allocate(eam_Z_r(n_eam_points,n_size_atype)) |
223 |
|
|
allocate(eam_rho_r(n_eam_points,n_size_atype)) |
224 |
|
|
allocate(eam_phi_r(n_eam_points,n_size_atype)) |
225 |
|
|
allocate(eam_F_rho_pp(n_eam_points,n_size_atype)) |
226 |
|
|
allocate(eam_Z_r_pp(n_eam_points,n_size_atype)) |
227 |
|
|
allocate(eam_rho_r_pp(n_eam_points,n_size_atype)) |
228 |
|
|
allocate(eam_phi_r_pp(n_eam_points,n_size_atype)) |
229 |
|
|
|
230 |
|
|
end subroutine allocate_eam_module |
231 |
|
|
|
232 |
|
|
subroutine deallocate_eam_module() |
233 |
|
|
|
234 |
|
|
deallocate(eam_atype) |
235 |
|
|
deallocate(eam_drho) |
236 |
|
|
deallocate(eam_dr) |
237 |
|
|
deallocate(eam_nr) |
238 |
|
|
deallocate(eam_nrho) |
239 |
|
|
deallocate(eam_lattice) |
240 |
|
|
deallocate(eam_atype_map) |
241 |
|
|
deallocate(eam_rvals) |
242 |
|
|
deallocate(eam_rhovals) |
243 |
|
|
deallocate(eam_rcut) |
244 |
|
|
deallocate(eam_Z_r) |
245 |
|
|
deallocate(eam_rho_r) |
246 |
|
|
deallocate(eam_phi_r) |
247 |
|
|
deallocate(eam_F_rho_pp) |
248 |
|
|
deallocate(eam_Z_r_pp) |
249 |
|
|
deallocate(eam_rho_r_pp) |
250 |
|
|
deallocate(eam_phi_r_pp) |
251 |
|
|
|
252 |
|
|
end subroutine deallocate_eam_module |
253 |
|
|
|
254 |
|
|
|
255 |
|
|
subroutine calc_eam_pair(atom1,atom2,d,rij,r2,pot,f,do_pot,do_stress) |
256 |
|
|
!Arguments |
257 |
|
|
integer, intent(in) :: atom1, atom2 |
258 |
|
|
real( kind = dp ), intent(in) :: rij, r2 |
259 |
|
|
real( kind = dp ) :: pot |
260 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: f |
261 |
|
|
real( kind = dp ), intent(in), dimension(3) :: d |
262 |
|
|
logical, intent(in) :: do_pot, do_stress |
263 |
|
|
|
264 |
|
|
!Local Variables |
265 |
|
|
|
266 |
|
|
|
267 |
|
|
if (rij .lt. EAM_rcut) then |
268 |
|
|
|
269 |
|
|
r = dsqrt(rijsq) |
270 |
|
|
efr(1,j) = -rxij |
271 |
|
|
efr(2,j) = -ryij |
272 |
|
|
efr(3,j) = -rzij |
273 |
|
|
|
274 |
|
|
|
275 |
|
|
call calc_eam_rho(r, rha, drha, d2rha, atype1) |
276 |
|
|
call calc_eam_phi(r, pha, dpha, d2pha, atype1) |
277 |
|
|
rci = eam_rcut(eam_atype_map(atype1)) |
278 |
|
|
#ifdef MPI |
279 |
|
|
atype2 = ident_col(j) |
280 |
|
|
#else |
281 |
|
|
atype2 = ident(j) |
282 |
|
|
#endif |
283 |
|
|
|
284 |
|
|
call calc_eam_rho(r, rhb, drhb, d2rhb, atype2) |
285 |
|
|
call calc_eam_phi(r, phb, dphb, d2phb, atype2) |
286 |
|
|
rcj = eam_rcut(eam_atype_map(atype2)) |
287 |
|
|
|
288 |
|
|
phab = 0.0E0_DP |
289 |
|
|
dvpdr = 0.0E0_DP |
290 |
|
|
d2vpdrdr = 0.0E0_DP |
291 |
|
|
|
292 |
|
|
if (r.lt.rci) then |
293 |
|
|
phab = phab + 0.5E0_DP*(rhb/rha)*pha |
294 |
|
|
dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + & |
295 |
|
|
pha*((drhb/rha) - (rhb*drha/rha/rha))) |
296 |
|
|
d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rhb/rha)*d2pha + & |
297 |
|
|
2.0E0_DP*dpha*((drhb/rha) - (rhb*drha/rha/rha)) + & |
298 |
|
|
pha*((d2rhb/rha) - 2.0E0_DP*(drhb*drha/rha/rha) + & |
299 |
|
|
(2.0E0_DP*rhb*drha*drha/rha/rha/rha) - (rhb*d2rha/rha/rha))) |
300 |
|
|
endif |
301 |
|
|
|
302 |
|
|
|
303 |
|
|
if (r.lt.rcj) then |
304 |
|
|
phab = phab + 0.5E0_DP*(rha/rhb)*phb |
305 |
|
|
dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + & |
306 |
|
|
phb*((drha/rhb) - (rha*drhb/rhb/rhb))) |
307 |
|
|
d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rha/rhb)*d2phb + & |
308 |
|
|
2.0E0_DP*dphb*((drha/rhb) - (rha*drhb/rhb/rhb)) + & |
309 |
|
|
phb*((d2rha/rhb) - 2.0E0_DP*(drha*drhb/rhb/rhb) + & |
310 |
|
|
(2.0E0_DP*rha*drhb*drhb/rhb/rhb/rhb) - (rha*d2rhb/rhb/rhb))) |
311 |
|
|
endif |
312 |
|
|
|
313 |
|
|
|
314 |
|
|
#ifdef MPI |
315 |
|
|
|
316 |
|
|
e_row(i) = e_row(i) + phab*0.5 |
317 |
|
|
e_col(i) = e_col(i) + phab*0.5 |
318 |
|
|
#else |
319 |
|
|
if (do_pot) pot = pot + phab |
320 |
|
|
#endif |
321 |
|
|
|
322 |
|
|
drhoidr = drha |
323 |
|
|
drhojdr = drhb |
324 |
|
|
|
325 |
|
|
d2rhoidrdr = d2rha |
326 |
|
|
d2rhojdrdr = d2rhb |
327 |
|
|
#ifdef MPI |
328 |
|
|
dudr = drhojdr*dfrhodrho_row(i)+drhoidr*dfrhodrho_col(j) & |
329 |
|
|
+ dvpdr |
330 |
|
|
|
331 |
|
|
if (nmflag) then |
332 |
|
|
d2 = d2vpdrdr + & |
333 |
|
|
d2rhoidrdr*dfrhodrho_col(j) + & |
334 |
|
|
d2rhojdrdr*dfrhodrho_row(i) + & |
335 |
|
|
drhoidr*drhoidr*d2frhodrhodrho_col(j) + & |
336 |
|
|
drhojdr*drhojdr*d2frhodrhodrho_row(i) |
337 |
|
|
endif |
338 |
|
|
#else |
339 |
|
|
dudr = drhojdr*dfrhodrho(i)+drhoidr*dfrhodrho(j) & |
340 |
|
|
+ dvpdr |
341 |
|
|
|
342 |
|
|
d2 = d2vpdrdr + & |
343 |
|
|
d2rhoidrdr*dfrhodrho(j) + & |
344 |
|
|
d2rhojdrdr*dfrhodrho(i) + & |
345 |
|
|
drhoidr*drhoidr*d2frhodrhodrho(j) + & |
346 |
|
|
drhojdr*drhojdr*d2frhodrhodrho(i) |
347 |
|
|
#endif |
348 |
|
|
|
349 |
|
|
|
350 |
|
|
do dim = 1, 3 |
351 |
|
|
|
352 |
|
|
drdx1 = efr(dim,j) / r |
353 |
|
|
ftmp = dudr * drdx1 |
354 |
|
|
|
355 |
|
|
#ifdef MPI |
356 |
|
|
f_col(dim,j) = f_col(dim,j) - ftmp |
357 |
|
|
f_row(dim,i) = f_row(dim,i) + ftmp |
358 |
|
|
#else |
359 |
|
|
f(dim,j) = f(dim,j) - ftmp |
360 |
|
|
f(dim,i) = f(dim,i) + ftmp |
361 |
|
|
#endif |
362 |
|
|
|
363 |
|
|
if (nmflag) then |
364 |
|
|
idim = 3 * (i-1) + dim |
365 |
|
|
jdim = 3 * (j-1) + dim |
366 |
|
|
|
367 |
|
|
do dim2 = 1, 3 |
368 |
|
|
|
369 |
|
|
kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r |
370 |
|
|
kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r |
371 |
|
|
|
372 |
|
|
if (dim.eq.dim2) then |
373 |
|
|
kt3 = dudr / r |
374 |
|
|
else |
375 |
|
|
kt3 = 0.0E0_DP |
376 |
|
|
endif |
377 |
|
|
|
378 |
|
|
! The factor of 2 below is to compensate for |
379 |
|
|
! overcounting. |
380 |
|
|
! Mass weighting is done separately... |
381 |
|
|
|
382 |
|
|
ktmp = (kt1+kt2+kt3)/2.0E0_DP |
383 |
|
|
idim2 = 3 * (i-1) + dim2 |
384 |
|
|
jdim2 = 3 * (j-1) + dim2 |
385 |
|
|
|
386 |
|
|
d(idim, idim2) = d(idim,idim2) + ktmp |
387 |
|
|
d(idim2, idim) = d(idim2,idim) + ktmp |
388 |
|
|
|
389 |
|
|
d(idim, jdim2) = d(idim,jdim2) - ktmp |
390 |
|
|
d(idim2, jdim) = d(idim2,jdim) - ktmp |
391 |
|
|
|
392 |
|
|
d(jdim, idim2) = d(jdim,idim2) - ktmp |
393 |
|
|
d(jdim2, idim) = d(jdim2,idim) - ktmp |
394 |
|
|
|
395 |
|
|
d(jdim, jdim2) = d(jdim,jdim2) + ktmp |
396 |
|
|
d(jdim2, jdim) = d(jdim2,jdim) + ktmp |
397 |
|
|
|
398 |
|
|
enddo |
399 |
|
|
endif |
400 |
|
|
enddo |
401 |
|
|
|
402 |
|
|
endif |
403 |
|
|
enddo |
404 |
|
|
endif |
405 |
|
|
|
406 |
|
|
enddo |
407 |
|
|
|
408 |
|
|
|
409 |
|
|
|
410 |
|
|
|
411 |
|
|
end subroutine calc_eam_pair |
412 |
|
|
|
413 |
|
|
subroutine calc_eam_rho(r, rho, drho, d2rho, atype) |
414 |
|
|
|
415 |
|
|
! include 'headers/sizes.h' |
416 |
|
|
|
417 |
|
|
|
418 |
|
|
integer atype, etype, number_r |
419 |
|
|
real( kind = DP ) :: r, rho, drho, d2rho |
420 |
|
|
integer :: i |
421 |
|
|
|
422 |
|
|
|
423 |
|
|
etype = eam_atype_map(atype) |
424 |
|
|
|
425 |
|
|
if (r.lt.eam_rcut(etype)) then |
426 |
|
|
number_r = eam_nr(etype) |
427 |
|
|
call eam_splint(etype, number_r, eam_rvals, eam_rho_r, & |
428 |
|
|
eam_rho_r_pp, r, rho, drho, d2rho) |
429 |
|
|
else |
430 |
|
|
rho = 0.0E0_DP |
431 |
|
|
drho = 0.0E0_DP |
432 |
|
|
d2rho = 0.0E0_DP |
433 |
|
|
endif |
434 |
|
|
|
435 |
|
|
return |
436 |
|
|
end subroutine calc_eam_rho |
437 |
|
|
|
438 |
|
|
subroutine calc_eam_frho(dens, u, u1, u2, atype) |
439 |
|
|
|
440 |
|
|
! include 'headers/sizes.h' |
441 |
|
|
|
442 |
|
|
integer atype, etype, number_rho |
443 |
|
|
real( kind = DP ) :: dens, u, u1, u2 |
444 |
|
|
real( kind = DP ) :: rho_vals |
445 |
|
|
|
446 |
|
|
etype = eam_atype_map(atype) |
447 |
|
|
number_rho = eam_nrho(etype) |
448 |
|
|
if (dens.lt.eam_rhovals(number_rho, etype)) then |
449 |
|
|
call eam_splint(etype, number_rho, eam_rhovals, eam_f_rho, & |
450 |
|
|
eam_f_rho_pp, dens, u, u1, u2) |
451 |
|
|
else |
452 |
|
|
rho_vals = eam_rhovals(number_rho,etype) |
453 |
|
|
call eam_splint(etype, number_rho, eam_rhovals, eam_f_rho, & |
454 |
|
|
eam_f_rho_pp, rho_vals, u, u1, u2) |
455 |
|
|
endif |
456 |
|
|
|
457 |
|
|
return |
458 |
|
|
end subroutine calc_eam_frho |
459 |
|
|
|
460 |
|
|
subroutine calc_eam_phi(r, phi, dphi, d2phi, atype) |
461 |
|
|
|
462 |
|
|
|
463 |
|
|
|
464 |
|
|
|
465 |
|
|
integer atype, etype, number_r |
466 |
|
|
real( kind = DP ) :: r, phi, dphi, d2phi |
467 |
|
|
|
468 |
|
|
etype = eam_atype_map(atype) |
469 |
|
|
|
470 |
|
|
if (r.lt.eam_rcut(etype)) then |
471 |
|
|
number_r = eam_nr(etype) |
472 |
|
|
call eam_splint(etype, number_r, eam_rvals, eam_phi_r, & |
473 |
|
|
eam_phi_r_pp, r, phi, dphi, d2phi) |
474 |
|
|
else |
475 |
|
|
phi = 0.0E0_DP |
476 |
|
|
dphi = 0.0E0_DP |
477 |
|
|
d2phi = 0.0E0_DP |
478 |
|
|
endif |
479 |
|
|
|
480 |
|
|
return |
481 |
|
|
end subroutine calc_eam_phi |
482 |
|
|
|
483 |
|
|
|
484 |
|
|
subroutine eam_splint(atype, nx, xa, ya, yppa, x, y, dy, d2y) |
485 |
|
|
|
486 |
|
|
! include 'headers/sizes.h' |
487 |
|
|
|
488 |
|
|
real( kind = DP ), dimension(:,:) :: xa |
489 |
|
|
real( kind = DP ), dimension(:,:) :: ya |
490 |
|
|
real( kind = DP ), dimension(:,:) :: yppa |
491 |
|
|
real( kind = DP ) :: x, y, dy, d2y |
492 |
|
|
real( kind = DP ) :: del, h, a, b, c, d |
493 |
|
|
|
494 |
|
|
|
495 |
|
|
integer atype, nx, j |
496 |
|
|
|
497 |
|
|
|
498 |
|
|
! this spline code assumes that the x points are equally spaced |
499 |
|
|
! do not attempt to use this code if they are not. |
500 |
|
|
|
501 |
|
|
|
502 |
|
|
! find the closest point with a value below our own: |
503 |
|
|
j = FLOOR(dble(nx-1) * (x - xa(1,atype)) / (xa(nx,atype) - xa(1,atype))) + 1 |
504 |
|
|
|
505 |
|
|
! check to make sure we're inside the spline range: |
506 |
|
|
if ((j.gt.nx).or.(j.lt.1)) call error('eam_splint', & |
507 |
|
|
'x is outside bounds of spline') |
508 |
|
|
|
509 |
|
|
! check to make sure we haven't screwed up the calculation of j: |
510 |
|
|
if ((x.lt.xa(j,atype)).or.(x.gt.xa(j+1,atype))) then |
511 |
|
|
if (j.ne.nx) then |
512 |
|
|
call error('eam_splint', & |
513 |
|
|
'x is outside bounding range') |
514 |
|
|
endif |
515 |
|
|
endif |
516 |
|
|
|
517 |
|
|
del = xa(j+1,atype) - x |
518 |
|
|
h = xa(j+1,atype) - xa(j,atype) |
519 |
|
|
|
520 |
|
|
a = del / h |
521 |
|
|
b = 1.0E0_DP - a |
522 |
|
|
c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP |
523 |
|
|
d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP |
524 |
|
|
|
525 |
|
|
y = a*ya(j,atype) + b*ya(j+1,atype) + c*yppa(j,atype) + d*yppa(j+1,atype) |
526 |
|
|
|
527 |
|
|
dy = (ya(j+1,atype)-ya(j,atype))/h & |
528 |
|
|
- (3.0E0_DP*a*a - 1.0E0_DP)*h*yppa(j,atype)/6.0E0_DP & |
529 |
|
|
+ (3.0E0_DP*b*b - 1.0E0_DP)*h*yppa(j+1,atype)/6.0E0_DP |
530 |
|
|
|
531 |
|
|
d2y = a*yppa(j,atype) + b*yppa(j+1,atype) |
532 |
|
|
|
533 |
|
|
return |
534 |
|
|
end subroutine eam_splint |
535 |
|
|
|
536 |
|
|
subroutine eam_spline(atype, nx, xa, ya, yppa, yp1, ypn, boundary) |
537 |
|
|
|
538 |
|
|
! include 'headers/sizes.h' |
539 |
|
|
|
540 |
|
|
|
541 |
|
|
! yp1 and ypn are the first derivatives of y at the two endpoints |
542 |
|
|
! if boundary is 'L' the lower derivative is used |
543 |
|
|
! if boundary is 'U' the upper derivative is used |
544 |
|
|
! if boundary is 'B' then both derivatives are used |
545 |
|
|
! if boundary is anything else, then both derivatives are assumed to be 0 |
546 |
|
|
|
547 |
|
|
integer nx, i, k, atype, max_array_size |
548 |
|
|
|
549 |
|
|
real( kind = DP ), dimension(:,:) :: xa |
550 |
|
|
real( kind = DP ), dimension(:,:) :: ya |
551 |
|
|
real( kind = DP ), dimension(:,:) :: yppa |
552 |
|
|
real( kind = DP ), allocatable, dimension(:) :: u |
553 |
|
|
real( kind = DP ) :: yp1,ypn,un,qn,sig,p |
554 |
|
|
character boundary |
555 |
|
|
|
556 |
|
|
max_array_size = size(xa,1) |
557 |
|
|
allocate(u(max_array_size)) |
558 |
|
|
|
559 |
|
|
|
560 |
|
|
if ((boundary.eq.'l').or.(boundary.eq.'L').or. & |
561 |
|
|
(boundary.eq.'b').or.(boundary.eq.'B')) then |
562 |
|
|
yppa(1, atype) = -0.5E0_DP |
563 |
|
|
u(1) = (3.0E0_DP/(xa(2,atype)-xa(1,atype)))*((ya(2,atype)-& |
564 |
|
|
ya(1,atype))/(xa(2,atype)-xa(1,atype))-yp1) |
565 |
|
|
else |
566 |
|
|
yppa(1,atype) = 0.0E0_DP |
567 |
|
|
u(1) = 0.0E0_DP |
568 |
|
|
endif |
569 |
|
|
|
570 |
|
|
do i = 2, nx - 1 |
571 |
|
|
sig = (xa(i,atype) - xa(i-1,atype)) / (xa(i+1,atype) - xa(i-1,atype)) |
572 |
|
|
p = sig * yppa(i-1,atype) + 2.0E0_DP |
573 |
|
|
yppa(i,atype) = (sig - 1.0E0_DP) / p |
574 |
|
|
u(i) = (6.0E0_DP*((ya(i+1,atype)-ya(i,atype))/(xa(i+1,atype)-xa(i,atype)) - & |
575 |
|
|
(ya(i,atype)-ya(i-1,atype))/(xa(i,atype)-xa(i-1,atype)))/ & |
576 |
|
|
(xa(i+1,atype)-xa(i-1,atype)) - sig * u(i-1))/p |
577 |
|
|
enddo |
578 |
|
|
|
579 |
|
|
if ((boundary.eq.'u').or.(boundary.eq.'U').or. & |
580 |
|
|
(boundary.eq.'b').or.(boundary.eq.'B')) then |
581 |
|
|
qn = 0.5E0_DP |
582 |
|
|
un = (3.0E0_DP/(xa(nx,atype)-xa(nx-1,atype)))* & |
583 |
|
|
(ypn-(ya(nx,atype)-ya(nx-1,atype))/(xa(nx,atype)-xa(nx-1,atype))) |
584 |
|
|
else |
585 |
|
|
qn = 0.0E0_DP |
586 |
|
|
un = 0.0E0_DP |
587 |
|
|
endif |
588 |
|
|
|
589 |
|
|
yppa(nx,atype)=(un-qn*u(nx-1))/(qn*yppa(nx-1,atype)+1.0E0_DP) |
590 |
|
|
|
591 |
|
|
do k = nx-1, 1, -1 |
592 |
|
|
yppa(k,atype)=yppa(k,atype)*yppa(k+1,atype)+u(k) |
593 |
|
|
enddo |
594 |
|
|
|
595 |
|
|
deallocate(u) |
596 |
|
|
return |
597 |
|
|
end subroutine eam_spline |
598 |
|
|
|
599 |
|
|
|
600 |
|
|
|
601 |
|
|
|
602 |
|
|
end module calc_eam |