ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/eam.F90
Revision: 1309
Committed: Tue Oct 21 18:23:31 2008 UTC (16 years, 9 months ago) by gezelter
File size: 22278 byte(s)
Log Message:
many changes to keep track of particle potentials for both bonded and non-bonded interactions

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
40 !!
41
42 module eam
43 use definitions
44 use simulation
45 use force_globals
46 use status
47 use atype_module
48 use vector_class
49 use interpolation
50 #ifdef IS_MPI
51 use mpiSimulation
52 #endif
53 implicit none
54 PRIVATE
55 #define __FORTRAN90
56 #include "UseTheForce/DarkSide/fInteractionMap.h"
57
58 logical, save :: EAM_FF_initialized = .false.
59 integer, save :: EAM_Mixing_Policy
60 real(kind = dp), save :: EAM_rcut
61 logical, save :: haveRcut = .false.
62
63 character(len = statusMsgSize) :: errMesg
64 integer :: eam_err
65
66 character(len = 200) :: errMsg
67 character(len=*), parameter :: RoutineName = "EAM MODULE"
68 !! Logical that determines if eam arrays should be zeroed
69 logical :: cleanme = .true.
70
71 type, private :: EAMtype
72 integer :: eam_atype
73 real( kind = DP ) :: eam_lattice
74 real( kind = DP ) :: eam_rcut
75 integer :: eam_atype_map
76 type(cubicSpline) :: rho
77 type(cubicSpline) :: Z
78 type(cubicSpline) :: F
79 type(cubicSpline) :: phi
80 end type EAMtype
81
82 !! Arrays for derivatives used in force calculation
83 real( kind = dp), dimension(:), allocatable :: frho
84 real( kind = dp), dimension(:), allocatable :: rho
85 real( kind = dp), dimension(:), allocatable :: dfrhodrho
86
87 !! Arrays for MPI storage
88 #ifdef IS_MPI
89 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
90 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
91 real( kind = dp),save, dimension(:), allocatable :: frho_row
92 real( kind = dp),save, dimension(:), allocatable :: frho_col
93 real( kind = dp),save, dimension(:), allocatable :: rho_row
94 real( kind = dp),save, dimension(:), allocatable :: rho_col
95 real( kind = dp),save, dimension(:), allocatable :: rho_tmp
96 #endif
97
98 type, private :: EAMTypeList
99 integer :: n_eam_types = 0
100 integer :: currentAddition = 0
101
102 type (EAMtype), pointer :: EAMParams(:) => null()
103 integer, pointer :: atidToEAMType(:) => null()
104 end type EAMTypeList
105
106 type (eamTypeList), save :: EAMList
107
108 !! standard eam stuff
109
110
111 public :: init_EAM_FF
112 public :: setCutoffEAM
113 public :: do_eam_pair
114 public :: newEAMtype
115 public :: calc_eam_prepair_rho
116 public :: calc_eam_preforce_Frho
117 public :: clean_EAM
118 public :: destroyEAMTypes
119 public :: getEAMCut
120 public :: lookupEAMSpline
121 public :: lookupEAMSpline1d
122
123 contains
124
125 subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
126 eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho, c_ident, status)
127 real (kind = dp ) :: lattice_constant
128 integer :: eam_nrho
129 real (kind = dp ) :: eam_drho
130 integer :: eam_nr
131 real (kind = dp ) :: eam_dr
132 real (kind = dp ) :: rcut
133 real (kind = dp ), dimension(eam_nr) :: eam_Z_r, rvals
134 real (kind = dp ), dimension(eam_nr) :: eam_rho_r, eam_phi_r
135 real (kind = dp ), dimension(eam_nrho) :: eam_F_rho, rhovals
136 integer :: c_ident
137 integer :: status
138
139 integer :: nAtypes,nEAMTypes,myATID
140 integer :: maxVals
141 integer :: alloc_stat
142 integer :: current, j
143 integer,pointer :: Matchlist(:) => null()
144
145 status = 0
146
147 !! Assume that atypes has already been set and get the total number of types in atypes
148 !! Also assume that every member of atypes is a EAM model.
149
150 ! check to see if this is the first time into
151 if (.not.associated(EAMList%EAMParams)) then
152 call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList)
153 EAMList%n_eam_types = nEAMtypes
154 allocate(EAMList%EAMParams(nEAMTypes))
155 nAtypes = getSize(atypes)
156 allocate(EAMList%atidToEAMType(nAtypes))
157 end if
158
159 EAMList%currentAddition = EAMList%currentAddition + 1
160 current = EAMList%currentAddition
161
162 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
163 EAMList%atidToEAMType(myATID) = current
164
165 EAMList%EAMParams(current)%eam_atype = c_ident
166 EAMList%EAMParams(current)%eam_lattice = lattice_constant
167 EAMList%EAMParams(current)%eam_rcut = rcut
168
169 ! Build array of r values
170 do j = 1, eam_nr
171 rvals(j) = real(j-1,kind=dp) * eam_dr
172 end do
173
174 ! Build array of rho values
175 do j = 1, eam_nrho
176 rhovals(j) = real(j-1,kind=dp) * eam_drho
177 end do
178
179 ! convert from eV to kcal / mol:
180 do j = 1, eam_nrho
181 eam_F_rho(j) = eam_F_rho(j) * 23.06054E0_DP
182 end do
183
184 ! precompute the pair potential and get it into kcal / mol:
185 eam_phi_r(1) = 0.0E0_DP
186 do j = 2, eam_nr
187 eam_phi_r(j) = 331.999296E0_DP * (eam_Z_r(j)**2) / rvals(j)
188 enddo
189
190 call newSpline(EAMList%EAMParams(current)%rho, rvals, eam_rho_r, .true.)
191 call newSpline(EAMList%EAMParams(current)%Z, rvals, eam_Z_r, .true.)
192 call newSpline(EAMList%EAMParams(current)%F, rhovals, eam_F_rho, .true.)
193 call newSpline(EAMList%EAMParams(current)%phi, rvals, eam_phi_r, .true.)
194 end subroutine newEAMtype
195
196
197 ! kills all eam types entered and sets simulation to uninitalized
198 subroutine destroyEAMtypes()
199 integer :: i
200 type(EAMType), pointer :: tempEAMType=>null()
201
202 do i = 1, EAMList%n_eam_types
203 tempEAMType => eamList%EAMParams(i)
204 call deallocate_EAMType(tempEAMType)
205 end do
206 if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
207 eamList%EAMParams => null()
208
209 eamList%n_eam_types = 0
210 eamList%currentAddition = 0
211 end subroutine destroyEAMtypes
212
213 function getEAMCut(atomID) result(cutValue)
214 integer, intent(in) :: atomID
215 integer :: eamID
216 real(kind=dp) :: cutValue
217
218 eamID = EAMList%atidToEAMType(atomID)
219 cutValue = EAMList%EAMParams(eamID)%eam_rcut
220 end function getEAMCut
221
222 subroutine init_EAM_FF(status)
223 integer :: status
224 integer :: i,j
225 real(kind=dp) :: current_rcut_max
226 #ifdef IS_MPI
227 integer :: nAtomsInRow
228 integer :: nAtomsInCol
229 #endif
230 integer :: alloc_stat
231 integer :: number_r, number_rho
232
233 status = 0
234 if (EAMList%currentAddition == 0) then
235 call handleError("init_EAM_FF","No members in EAMList")
236 status = -1
237 return
238 end if
239
240 !! Allocate arrays for force calculation
241
242 #ifdef IS_MPI
243 nAtomsInRow = getNatomsInRow(plan_atom_row)
244 nAtomsInCol = getNatomsInCol(plan_atom_col)
245 #endif
246
247 if (allocated(frho)) deallocate(frho)
248 allocate(frho(nlocal),stat=alloc_stat)
249 if (alloc_stat /= 0) then
250 status = -1
251 return
252 end if
253
254 if (allocated(rho)) deallocate(rho)
255 allocate(rho(nlocal),stat=alloc_stat)
256 if (alloc_stat /= 0) then
257 status = -1
258 return
259 end if
260
261 if (allocated(dfrhodrho)) deallocate(dfrhodrho)
262 allocate(dfrhodrho(nlocal),stat=alloc_stat)
263 if (alloc_stat /= 0) then
264 status = -1
265 return
266 end if
267
268 #ifdef IS_MPI
269
270 if (allocated(rho_tmp)) deallocate(rho_tmp)
271 allocate(rho_tmp(nlocal),stat=alloc_stat)
272 if (alloc_stat /= 0) then
273 status = -1
274 return
275 end if
276
277 if (allocated(frho_row)) deallocate(frho_row)
278 allocate(frho_row(nAtomsInRow),stat=alloc_stat)
279 if (alloc_stat /= 0) then
280 status = -1
281 return
282 end if
283 if (allocated(rho_row)) deallocate(rho_row)
284 allocate(rho_row(nAtomsInRow),stat=alloc_stat)
285 if (alloc_stat /= 0) then
286 status = -1
287 return
288 end if
289 if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
290 allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
291 if (alloc_stat /= 0) then
292 status = -1
293 return
294 end if
295
296 ! Now do column arrays
297
298 if (allocated(frho_col)) deallocate(frho_col)
299 allocate(frho_col(nAtomsInCol),stat=alloc_stat)
300 if (alloc_stat /= 0) then
301 status = -1
302 return
303 end if
304 if (allocated(rho_col)) deallocate(rho_col)
305 allocate(rho_col(nAtomsInCol),stat=alloc_stat)
306 if (alloc_stat /= 0) then
307 status = -1
308 return
309 end if
310 if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
311 allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
312 if (alloc_stat /= 0) then
313 status = -1
314 return
315 end if
316
317 #endif
318
319 end subroutine init_EAM_FF
320
321 subroutine setCutoffEAM(rcut)
322 real(kind=dp) :: rcut
323 EAM_rcut = rcut
324 end subroutine setCutoffEAM
325
326 subroutine clean_EAM()
327
328 ! clean non-IS_MPI first
329 frho = 0.0_dp
330 rho = 0.0_dp
331 dfrhodrho = 0.0_dp
332 ! clean MPI if needed
333 #ifdef IS_MPI
334 frho_row = 0.0_dp
335 frho_col = 0.0_dp
336 rho_row = 0.0_dp
337 rho_col = 0.0_dp
338 rho_tmp = 0.0_dp
339 dfrhodrho_row = 0.0_dp
340 dfrhodrho_col = 0.0_dp
341 #endif
342 end subroutine clean_EAM
343
344 subroutine deallocate_EAMType(thisEAMType)
345 type (EAMtype), pointer :: thisEAMType
346
347 call deleteSpline(thisEAMType%F)
348 call deleteSpline(thisEAMType%rho)
349 call deleteSpline(thisEAMType%phi)
350 call deleteSpline(thisEAMType%Z)
351
352 end subroutine deallocate_EAMType
353
354 !! Calculates rho_r
355 subroutine calc_eam_prepair_rho(atom1,atom2,d,r,rijsq)
356 integer :: atom1, atom2
357 real(kind = dp), dimension(3) :: d
358 real(kind = dp), intent(inout) :: r
359 real(kind = dp), intent(inout) :: rijsq
360 ! value of electron density rho do to atom i at atom j
361 real(kind = dp) :: rho_i_at_j
362 ! value of electron density rho do to atom j at atom i
363 real(kind = dp) :: rho_j_at_i
364 integer :: eam_err
365
366 integer :: atid1, atid2 ! Global atid
367 integer :: myid_atom1 ! EAM atid
368 integer :: myid_atom2
369
370 ! check to see if we need to be cleaned at the start of a force loop
371
372 #ifdef IS_MPI
373 Atid1 = Atid_row(Atom1)
374 Atid2 = Atid_col(Atom2)
375 #else
376 Atid1 = Atid(Atom1)
377 Atid2 = Atid(Atom2)
378 #endif
379
380 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
381 Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
382
383 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
384
385 call lookupEAMSpline(EAMList%EAMParams(myid_atom1)%rho, r, &
386 rho_i_at_j)
387
388 #ifdef IS_MPI
389 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
390 #else
391 rho(atom2) = rho(atom2) + rho_i_at_j
392 #endif
393 endif
394
395 if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then
396
397 call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
398 rho_j_at_i)
399
400 #ifdef IS_MPI
401 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
402 #else
403 rho(atom1) = rho(atom1) + rho_j_at_i
404 #endif
405 endif
406 end subroutine calc_eam_prepair_rho
407
408
409 !! Calculate the functional F(rho) for all local atoms
410 subroutine calc_eam_preforce_Frho(nlocal, pot, particle_pot)
411 integer :: nlocal
412 real(kind=dp) :: pot
413 integer :: i, j
414 integer :: atom
415 real(kind=dp) :: U,U1
416 real( kind = dp ), dimension(nlocal) :: particle_pot
417 integer :: atype1
418 integer :: me, atid1
419
420 cleanme = .true.
421 !! Scatter the electron density from pre-pair calculation back to
422 !! local atoms
423 #ifdef IS_MPI
424 call scatter(rho_row,rho,plan_atom_row,eam_err)
425 if (eam_err /= 0 ) then
426 write(errMsg,*) " Error scattering rho_row into rho"
427 call handleError(RoutineName,errMesg)
428 endif
429 call scatter(rho_col,rho_tmp,plan_atom_col,eam_err)
430 if (eam_err /= 0 ) then
431 write(errMsg,*) " Error scattering rho_col into rho"
432 call handleError(RoutineName,errMesg)
433 endif
434
435 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
436 #endif
437
438 !! Calculate F(rho) and derivative
439 do atom = 1, nlocal
440 atid1 = atid(atom)
441 me = eamList%atidToEAMtype(atid1)
442
443 call lookupEAMSpline1d(EAMList%EAMParams(me)%F, rho(atom), &
444 u, u1)
445
446 frho(atom) = u
447 dfrhodrho(atom) = u1
448 pot = pot + u
449 particle_pot(atom) = particle_pot(atom) + u
450
451 enddo
452
453 #ifdef IS_MPI
454 !! communicate f(rho) and derivatives back into row and column arrays
455 call gather(frho,frho_row,plan_atom_row, eam_err)
456 if (eam_err /= 0) then
457 call handleError("cal_eam_forces()","MPI gather frho_row failure")
458 endif
459 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, eam_err)
460 if (eam_err /= 0) then
461 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
462 endif
463 call gather(frho,frho_col,plan_atom_col, eam_err)
464 if (eam_err /= 0) then
465 call handleError("cal_eam_forces()","MPI gather frho_col failure")
466 endif
467 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, eam_err)
468 if (eam_err /= 0) then
469 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
470 endif
471 #endif
472
473
474 end subroutine calc_eam_preforce_Frho
475
476 !! Does EAM pairwise Force calculation.
477 subroutine do_eam_pair(atom1, atom2, d, rij, r2, sw, vpair, particle_pot, &
478 fpair, pot, f, do_pot)
479 !Arguments
480 integer, intent(in) :: atom1, atom2
481 real( kind = dp ), intent(in) :: rij, r2
482 real( kind = dp ) :: pot, sw, vpair
483 real( kind = dp ), dimension(nLocal) :: particle_pot
484 real( kind = dp ), dimension(3,nLocal) :: f
485 real( kind = dp ), intent(in), dimension(3) :: d
486 real( kind = dp ), intent(inout), dimension(3) :: fpair
487
488 logical, intent(in) :: do_pot
489
490 real( kind = dp ) :: drdx, drdy, drdz
491 real( kind = dp ) :: phab, pha, dvpdr
492 real( kind = dp ) :: rha, drha, dpha
493 real( kind = dp ) :: rhb, drhb, dphb
494 real( kind = dp ) :: dudr
495 real( kind = dp ) :: rci, rcj
496 real( kind = dp ) :: drhoidr, drhojdr
497 real( kind = dp ) :: Fx, Fy, Fz
498 real( kind = dp ) :: r, phb
499 real( kind = dp ) :: fshift1, fshift2, u1, u2
500
501 integer :: id1, id2
502 integer :: mytype_atom1
503 integer :: mytype_atom2
504 integer :: atid1, atid2
505
506 phab = 0.0E0_DP
507 dvpdr = 0.0E0_DP
508
509 if (rij .lt. EAM_rcut) then
510
511 #ifdef IS_MPI
512 atid1 = atid_row(atom1)
513 atid2 = atid_col(atom2)
514 #else
515 atid1 = atid(atom1)
516 atid2 = atid(atom2)
517 #endif
518
519 mytype_atom1 = EAMList%atidToEAMType(atid1)
520 mytype_atom2 = EAMList%atidTOEAMType(atid2)
521
522
523 ! get cutoff for atom 1
524 rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
525 ! get type specific cutoff for atom 2
526 rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
527
528 drdx = d(1)/rij
529 drdy = d(2)/rij
530 drdz = d(3)/rij
531
532 if (rij.lt.rci) then
533
534 ! Calculate rho and drho for atom1
535
536 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
537 rij, rha, drha)
538
539 ! Calculate Phi(r) for atom1.
540
541 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
542 rij, pha, dpha)
543
544 endif
545
546 if (rij.lt.rcj) then
547
548 ! Calculate rho and drho for atom2
549
550 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
551 rij, rhb, drhb)
552
553 ! Calculate Phi(r) for atom2.
554
555 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
556 rij, phb, dphb)
557
558 endif
559
560 if (rij.lt.rci) then
561 phab = phab + 0.5E0_DP*(rhb/rha)*pha
562 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
563 pha*((drhb/rha) - (rhb*drha/rha/rha)))
564 endif
565
566 if (rij.lt.rcj) then
567 phab = phab + 0.5E0_DP*(rha/rhb)*phb
568 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
569 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
570 endif
571
572 drhoidr = drha
573 drhojdr = drhb
574
575 #ifdef IS_MPI
576 dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) &
577 + dvpdr
578
579 #else
580 dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) &
581 + dvpdr
582 #endif
583
584 fx = dudr * drdx
585 fy = dudr * drdy
586 fz = dudr * drdz
587
588
589 #ifdef IS_MPI
590 if (do_pot) then
591 ! particle_pot is the difference between the full potential
592 ! and the full potential without the presence of a particular
593 ! particle (atom1).
594 !
595 ! This reduces the density at other particle locations, so
596 ! we need to recompute the density at atom2 assuming atom1
597 ! didn't contribute. This then requires recomputing the
598 ! density functional for atom2 as well.
599 !
600 ! Most of the particle_pot heavy lifting comes from the
601 ! pair interaction, and will be handled by vpair.
602
603 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, &
604 rho_row(atom1)-rhb, &
605 fshift1, u1)
606 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, &
607 rho_col(atom2)-rha, &
608 fshift2, u2)
609
610 ppot_Row(atom1) = ppot_Row(atom1) - frho_row(atom2) + fshift2
611 ppot_Col(atom2) = ppot_Col(atom2) - frho_col(atom1) + fshift1
612
613
614 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + phab*0.5
615 pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + phab*0.5
616 end if
617
618 f_Row(1,atom1) = f_Row(1,atom1) + fx
619 f_Row(2,atom1) = f_Row(2,atom1) + fy
620 f_Row(3,atom1) = f_Row(3,atom1) + fz
621
622 f_Col(1,atom2) = f_Col(1,atom2) - fx
623 f_Col(2,atom2) = f_Col(2,atom2) - fy
624 f_Col(3,atom2) = f_Col(3,atom2) - fz
625 #else
626
627 if(do_pot) then
628 ! particle_pot is the difference between the full potential
629 ! and the full potential without the presence of a particular
630 ! particle (atom1).
631 !
632 ! This reduces the density at other particle locations, so
633 ! we need to recompute the density at atom2 assuming atom1
634 ! didn't contribute. This then requires recomputing the
635 ! density functional for atom2 as well.
636 !
637 ! Most of the particle_pot heavy lifting comes from the
638 ! pair interaction, and will be handled by vpair.
639
640 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%F, &
641 rho(atom1)-rhb, &
642 fshift1, u1)
643 call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%F, &
644 rho(atom2)-rha, &
645 fshift2, u2)
646
647 particle_pot(atom1) = particle_pot(atom1) - frho(atom2) + fshift2
648 particle_pot(atom2) = particle_pot(atom2) - frho(atom1) + fshift1
649
650 pot = pot + phab
651 end if
652
653 f(1,atom1) = f(1,atom1) + fx
654 f(2,atom1) = f(2,atom1) + fy
655 f(3,atom1) = f(3,atom1) + fz
656
657 f(1,atom2) = f(1,atom2) - fx
658 f(2,atom2) = f(2,atom2) - fy
659 f(3,atom2) = f(3,atom2) - fz
660 #endif
661
662 vpair = vpair + phab
663
664 #ifdef IS_MPI
665 id1 = AtomRowToGlobal(atom1)
666 id2 = AtomColToGlobal(atom2)
667 #else
668 id1 = atom1
669 id2 = atom2
670 #endif
671
672 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
673
674 fpair(1) = fpair(1) + fx
675 fpair(2) = fpair(2) + fy
676 fpair(3) = fpair(3) + fz
677
678 endif
679 endif
680 end subroutine do_eam_pair
681
682 subroutine lookupEAMSpline(cs, xval, yval)
683
684 implicit none
685
686 type (cubicSpline), intent(in) :: cs
687 real( kind = DP ), intent(in) :: xval
688 real( kind = DP ), intent(out) :: yval
689 real( kind = DP ) :: dx
690 integer :: i, j
691 !
692 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
693 ! or is nearest to xval.
694
695 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
696
697 dx = xval - cs%x(j)
698 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
699
700 return
701 end subroutine lookupEAMSpline
702
703 subroutine lookupEAMSpline1d(cs, xval, yval, dydx)
704
705 implicit none
706
707 type (cubicSpline), intent(in) :: cs
708 real( kind = DP ), intent(in) :: xval
709 real( kind = DP ), intent(out) :: yval, dydx
710 real( kind = DP ) :: dx
711 integer :: i, j
712
713 ! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains
714 ! or is nearest to xval.
715
716
717 j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1))
718
719 dx = xval - cs%x(j)
720 yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j)))
721
722 dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j))
723
724 return
725 end subroutine lookupEAMSpline1d
726
727 end module eam