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

# Content
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), save :: 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
168 status = 0
169 if (EAMList%currentAddition == 0) then
170 call handleError("init_EAM_FF","No members in EAMList")
171 status = -1
172 return
173 end if
174
175
176 do i = 1, EAMList%currentAddition
177
178 ! Build array of r values
179
180 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 ! 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 EAMList%EAMParams(i)%eam_phi_r(1) = 0.0E0_DP
196 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
202
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
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 integer :: nlocal
254 #ifdef IS_MPI
255 integer :: nrow
256 integer :: ncol
257 #endif
258 integer :: alloc_stat
259
260
261 nlocal = getNlocal()
262 status = 0
263 #ifdef IS_MPI
264 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 end if
274 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
281 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
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
295 #ifdef IS_MPI
296
297 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 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
322
323 ! 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 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 #endif
351
352 end subroutine allocateEAM
353
354 !! 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 subroutine setCutoffEAM(rcut, status)
358 real(kind=dp) :: rcut
359 integer :: status
360 status = 0
361
362 EAM_rcut = rcut
363
364 end subroutine setCutoffEAM
365
366
367
368 subroutine clean_EAM()
369
370 ! clean non-IS_MPI first
371 frho = 0.0_dp
372 rho = 0.0_dp
373 dfrhodrho = 0.0_dp
374 ! clean MPI if needed
375 #ifdef IS_MPI
376 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 #endif
383 end subroutine clean_EAM
384
385
386
387 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 type (EAMType) :: thisEAMType
391 integer, optional :: stat
392 integer :: alloc_stat
393
394
395
396 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
450 end subroutine allocate_EAMType
451
452
453 subroutine deallocate_EAMType(thisEAMType)
454 type (EAMtype), pointer :: thisEAMType
455
456 ! 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
470 !! 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
481 ! we don't use the derivatives, dummy variables
482 real( kind = dp) :: drho,d2rho
483 integer :: eam_err
484
485 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
490
491
492
493 #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
501 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
502
503
504
505 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 #else
516 rho(atom2) = rho(atom2) + rho_i_at_j
517 #endif
518 ! write(*,*) atom1,atom2,r,rho_i_at_j
519 endif
520
521 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
528
529
530
531 #ifdef IS_MPI
532 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
533 #else
534 rho(atom1) = rho(atom1) + rho_j_at_i
535 #endif
536 endif
537
538
539 end subroutine calc_eam_prepair_rho
540
541
542
543
544 !! Calculate the functional F(rho) for all local atoms
545 subroutine calc_eam_preforce_Frho(nlocal,pot)
546 integer :: nlocal
547 real(kind=dp) :: pot
548 integer :: i,j
549 integer :: atom
550 real(kind=dp) :: U,U1,U2
551 integer :: atype1
552 integer :: me
553 integer :: n_rho_points
554
555
556 cleanme = .true.
557 !! Scatter the electron density from pre-pair calculation back to local atoms
558 #ifdef IS_MPI
559 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
571
572 !! 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 frho(atom) = u
596 dfrhodrho(atom) = u1
597 d2frhodrhodrho(atom) = u2
598 pot = pot + u
599
600 enddo
601
602
603
604 #ifdef IS_MPI
605 !! communicate f(rho) and derivatives back into row and column arrays
606 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
623
624
625
626
627 if (nmflag) then
628 call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_row)
629 call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_col)
630 endif
631 #endif
632
633
634 end subroutine calc_eam_preforce_Frho
635
636
637
638
639 !! Does EAM pairwise Force calculation.
640 subroutine do_eam_pair(atom1,atom2,d,rij,r2,pot,f,do_pot,do_stress)
641 !Arguments
642 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
649 real( kind = dp ) :: drdx,drdy,drdz
650 real( kind = dp ) :: d2
651 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
662 integer :: id1,id2
663 integer :: mytype_atom1
664 integer :: mytype_atom2
665
666 !Local Variables
667
668 ! write(*,*) "Frho: ", Frho(atom1)
669
670 phab = 0.0E0_DP
671 dvpdr = 0.0E0_DP
672 d2vpdrdr = 0.0E0_DP
673
674 if (rij .lt. EAM_rcut) then
675 #ifdef IS_MPI
676 !!!!! FIX ME
677 mytype_atom1 = atid_row(atom1)
678 #else
679 mytype_atom1 = atid(atom1)
680 #endif
681
682 drdx = d(1)/rij
683 drdy = d(2)/rij
684 drdz = d(3)/rij
685
686
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 #ifdef IS_MPI
704 mytype_atom2 = atid_col(atom2)
705 #else
706 mytype_atom2 = atid(atom2)
707 #endif
708
709 ! 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 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
739
740 if (rij.lt.rcj) then
741 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
750 drhoidr = drha
751 drhojdr = drhb
752
753 d2rhoidrdr = d2rha
754 d2rhojdrdr = d2rhb
755
756
757 #ifdef IS_MPI
758 dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) &
759 + dvpdr
760
761 #else
762 dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) &
763 + dvpdr
764 ! write(*,*) "Atom1,Atom2, dfrhodrho(atom1) dfrhodrho(atom2): ", atom1,atom2,dfrhodrho(atom1),dfrhodrho(atom2)
765 #endif
766
767 fx = dudr * drdx
768 fy = dudr * drdy
769 fz = dudr * drdz
770
771
772 #ifdef IS_MPI
773 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
778 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
782 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
787 if(do_pot) then
788 pot = pot + phab
789 end if
790
791 f(1,atom1) = f(1,atom1) + fx
792 f(2,atom1) = f(2,atom1) + fy
793 f(3,atom1) = f(3,atom1) + fz
794
795 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
800 if (nmflag) then
801
802 drhoidr = drha
803 drhojdr = drhb
804 d2rhoidrdr = d2rha
805 d2rhojdrdr = d2rhb
806
807 #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
816 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
824
825
826
827 if (do_stress) then
828
829 #ifdef IS_MPI
830 id1 = tagRow(atom1)
831 id2 = tagColumn(atom2)
832 #else
833 id1 = atom1
834 id2 = atom2
835 #endif
836
837 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
838
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 virial_Temp = virial_Temp + &
850 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
851
852 endif
853 endif
854 endif
855
856
857 end subroutine do_eam_pair
858
859
860 subroutine eam_splint(nx, xa, ya, yppa, x, y, dy, d2y)
861
862 integer :: atype, nx, j
863 real( kind = DP ), dimension(:) :: xa
864 real( kind = DP ), dimension(:) :: ya
865 real( kind = DP ), dimension(:) :: yppa
866 real( kind = DP ) :: x, y
867 real( kind = DP ) :: dy, d2y
868 real( kind = DP ) :: del, h, a, b, c, d
869 integer :: pp_arraySize
870
871
872 ! 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 j = FLOOR(real((nx-1),kind=dp) * (x - xa(1)) / (xa(nx) - xa(1))) + 1
878
879 ! check to make sure we're inside the spline range:
880 if ((j.gt.nx).or.(j.lt.1)) then
881 write(errMSG,*) "EAM_splint: x is outside bounds of spline"
882 call handleError(routineName,errMSG)
883 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 write(errMSG,*) "EAM_splint:",x," x is outside bounding range"
888 call handleError(routineName,errMSG)
889 endif
890 endif
891
892 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
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
910 end subroutine eam_splint
911
912
913 subroutine eam_spline(nx, xa, ya, yppa, yp1, ypn, boundary)
914
915
916 ! 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 character(len=*) :: boundary
930
931 ! 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 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
965 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
971 end subroutine eam_spline
972
973
974
975
976 end module eam