ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/eam.F90
Revision: 2717
Committed: Mon Apr 17 21:49:12 2006 UTC (18 years, 5 months ago) by gezelter
File size: 18933 byte(s)
Log Message:
Many performance improvements

File Contents

# User Rev Content
1 gezelter 1930 !!
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 gezelter 1608 module eam
43     use simulation
44     use force_globals
45     use status
46     use atype_module
47 chrisfen 2277 use vector_class
48 gezelter 2717 use interpolation
49 gezelter 1608 #ifdef IS_MPI
50     use mpiSimulation
51     #endif
52     implicit none
53     PRIVATE
54 chuckv 2355 #define __FORTRAN90
55     #include "UseTheForce/DarkSide/fInteractionMap.h"
56 gezelter 1608
57 chuckv 1624 INTEGER, PARAMETER :: DP = selected_real_kind(15)
58 gezelter 1608
59     logical, save :: EAM_FF_initialized = .false.
60     integer, save :: EAM_Mixing_Policy
61     real(kind = dp), save :: EAM_rcut
62     logical, save :: haveRcut = .false.
63    
64     character(len = statusMsgSize) :: errMesg
65     integer :: eam_err
66    
67     character(len = 200) :: errMsg
68     character(len=*), parameter :: RoutineName = "EAM MODULE"
69 gezelter 2204 !! Logical that determines if eam arrays should be zeroed
70 gezelter 1608 logical :: cleanme = .true.
71    
72     type, private :: EAMtype
73     integer :: eam_atype
74     real( kind = DP ) :: eam_lattice
75     real( kind = DP ) :: eam_rcut
76     integer :: eam_atype_map
77 gezelter 2717 type(cubicSpline) :: rho
78     type(cubicSpline) :: Z
79     type(cubicSpline) :: F
80     type(cubicSpline) :: phi
81 gezelter 1608 end type EAMtype
82    
83     !! Arrays for derivatives used in force calculation
84     real( kind = dp), dimension(:), allocatable :: frho
85     real( kind = dp), dimension(:), allocatable :: rho
86     real( kind = dp), dimension(:), allocatable :: dfrhodrho
87    
88 gezelter 2204 !! Arrays for MPI storage
89 gezelter 1608 #ifdef IS_MPI
90     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
91     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
92     real( kind = dp),save, dimension(:), allocatable :: frho_row
93     real( kind = dp),save, dimension(:), allocatable :: frho_col
94     real( kind = dp),save, dimension(:), allocatable :: rho_row
95     real( kind = dp),save, dimension(:), allocatable :: rho_col
96     real( kind = dp),save, dimension(:), allocatable :: rho_tmp
97     #endif
98    
99     type, private :: EAMTypeList
100     integer :: n_eam_types = 0
101     integer :: currentAddition = 0
102 gezelter 2204
103 gezelter 1608 type (EAMtype), pointer :: EAMParams(:) => null()
104 chuckv 2276 integer, pointer :: atidToEAMType(:) => null()
105 gezelter 1608 end type EAMTypeList
106    
107     type (eamTypeList), save :: EAMList
108    
109     !! standard eam stuff
110    
111    
112     public :: init_EAM_FF
113     public :: setCutoffEAM
114     public :: do_eam_pair
115     public :: newEAMtype
116     public :: calc_eam_prepair_rho
117     public :: calc_eam_preforce_Frho
118     public :: clean_EAM
119 chuckv 2169 public :: destroyEAMTypes
120 chrisfen 2277 public :: getEAMCut
121 gezelter 1608
122     contains
123    
124     subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
125 gezelter 2717 eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho, c_ident, status)
126 gezelter 1608 real (kind = dp ) :: lattice_constant
127     integer :: eam_nrho
128     real (kind = dp ) :: eam_drho
129     integer :: eam_nr
130     real (kind = dp ) :: eam_dr
131     real (kind = dp ) :: rcut
132 gezelter 2717 real (kind = dp ), dimension(eam_nr) :: eam_Z_r, rvals
133     real (kind = dp ), dimension(eam_nr) :: eam_rho_r, eam_phi_r
134     real (kind = dp ), dimension(eam_nrho) :: eam_F_rho, rhovals
135 chrisfen 2278 integer :: c_ident
136 gezelter 1608 integer :: status
137    
138 chuckv 2276 integer :: nAtypes,nEAMTypes,myATID
139 gezelter 1608 integer :: maxVals
140     integer :: alloc_stat
141 gezelter 2717 integer :: current, j
142 gezelter 1608 integer,pointer :: Matchlist(:) => null()
143    
144     status = 0
145    
146     !! Assume that atypes has already been set and get the total number of types in atypes
147     !! Also assume that every member of atypes is a EAM model.
148    
149     ! check to see if this is the first time into
150     if (.not.associated(EAMList%EAMParams)) then
151 chuckv 2276 call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList)
152     EAMList%n_eam_types = nEAMtypes
153     allocate(EAMList%EAMParams(nEAMTypes))
154     nAtypes = getSize(atypes)
155     allocate(EAMList%atidToEAMType(nAtypes))
156 gezelter 1608 end if
157    
158     EAMList%currentAddition = EAMList%currentAddition + 1
159     current = EAMList%currentAddition
160    
161 chrisfen 2278 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
162 chuckv 2276 EAMList%atidToEAMType(myATID) = current
163 gezelter 2204
164 chrisfen 2278 EAMList%EAMParams(current)%eam_atype = c_ident
165 gezelter 1608 EAMList%EAMParams(current)%eam_lattice = lattice_constant
166     EAMList%EAMParams(current)%eam_rcut = rcut
167    
168 gezelter 2717 ! Build array of r values
169     do j = 1, eam_nr
170     rvals(j) = real(j-1,kind=dp) * eam_dr
171     end do
172    
173     ! Build array of rho values
174     do j = 1, eam_nrho
175     rhovals(j) = real(j-1,kind=dp) * eam_drho
176     end do
177    
178     ! convert from eV to kcal / mol:
179     do j = 1, eam_nrho
180     eam_F_rho(j) = eam_F_rho(j) * 23.06054E0_DP
181     end do
182    
183     ! precompute the pair potential and get it into kcal / mol:
184     eam_phi_r(1) = 0.0E0_DP
185     do j = 2, eam_nr
186     eam_phi_r(j) = 331.999296E0_DP * (eam_Z_r(j)**2) / rvals(j)
187     enddo
188    
189     call newSpline(EAMList%EAMParams(current)%rho, rvals, rhovals, &
190     0.0d0, 0.0d0, .true.)
191     call newSpline(EAMList%EAMParams(current)%Z, rvals, eam_Z_r, &
192     0.0d0, 0.0d0, .true.)
193     call newSpline(EAMList%EAMParams(current)%F, rhovals, eam_F_rho, &
194     0.0d0, 0.0d0, .true.)
195     call newSpline(EAMList%EAMParams(current)%phi, rvals, eam_phi_r, &
196     0.0d0, 0.0d0, .true.)
197 gezelter 1608 end subroutine newEAMtype
198    
199    
200 chuckv 2169 ! kills all eam types entered and sets simulation to uninitalized
201     subroutine destroyEAMtypes()
202     integer :: i
203     type(EAMType), pointer :: tempEAMType=>null()
204 gezelter 1608
205 chuckv 2169 do i = 1, EAMList%n_eam_types
206     tempEAMType => eamList%EAMParams(i)
207     call deallocate_EAMType(tempEAMType)
208     end do
209     if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
210     eamList%EAMParams => null()
211 gezelter 2204
212 chuckv 2169 eamList%n_eam_types = 0
213     eamList%currentAddition = 0
214     end subroutine destroyEAMtypes
215    
216 chrisfen 2277 function getEAMCut(atomID) result(cutValue)
217     integer, intent(in) :: atomID
218 chrisfen 2278 integer :: eamID
219 chrisfen 2277 real(kind=dp) :: cutValue
220    
221 chrisfen 2278 eamID = EAMList%atidToEAMType(atomID)
222     cutValue = EAMList%EAMParams(eamID)%eam_rcut
223 chrisfen 2277 end function getEAMCut
224    
225 gezelter 1608 subroutine init_EAM_FF(status)
226     integer :: status
227     integer :: i,j
228     real(kind=dp) :: current_rcut_max
229 gezelter 2717 #ifdef IS_MPI
230     integer :: nAtomsInRow
231     integer :: nAtomsInCol
232     #endif
233 gezelter 1608 integer :: alloc_stat
234     integer :: number_r, number_rho
235    
236     status = 0
237     if (EAMList%currentAddition == 0) then
238     call handleError("init_EAM_FF","No members in EAMList")
239     status = -1
240     return
241     end if
242 gezelter 2717
243 gezelter 2204 !! Allocate arrays for force calculation
244 gezelter 2717
245 gezelter 1608 #ifdef IS_MPI
246     nAtomsInRow = getNatomsInRow(plan_atom_row)
247     nAtomsInCol = getNatomsInCol(plan_atom_col)
248     #endif
249    
250     if (allocated(frho)) deallocate(frho)
251     allocate(frho(nlocal),stat=alloc_stat)
252     if (alloc_stat /= 0) then
253     status = -1
254     return
255     end if
256 gezelter 2717
257 gezelter 1608 if (allocated(rho)) deallocate(rho)
258     allocate(rho(nlocal),stat=alloc_stat)
259     if (alloc_stat /= 0) then
260     status = -1
261     return
262     end if
263    
264     if (allocated(dfrhodrho)) deallocate(dfrhodrho)
265     allocate(dfrhodrho(nlocal),stat=alloc_stat)
266     if (alloc_stat /= 0) then
267     status = -1
268     return
269     end if
270    
271     #ifdef IS_MPI
272    
273     if (allocated(rho_tmp)) deallocate(rho_tmp)
274     allocate(rho_tmp(nlocal),stat=alloc_stat)
275     if (alloc_stat /= 0) then
276     status = -1
277     return
278     end if
279    
280     if (allocated(frho_row)) deallocate(frho_row)
281     allocate(frho_row(nAtomsInRow),stat=alloc_stat)
282     if (alloc_stat /= 0) then
283     status = -1
284     return
285     end if
286     if (allocated(rho_row)) deallocate(rho_row)
287     allocate(rho_row(nAtomsInRow),stat=alloc_stat)
288     if (alloc_stat /= 0) then
289     status = -1
290     return
291     end if
292     if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
293     allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
294     if (alloc_stat /= 0) then
295     status = -1
296     return
297     end if
298    
299 gezelter 2204 ! Now do column arrays
300 gezelter 1608
301     if (allocated(frho_col)) deallocate(frho_col)
302     allocate(frho_col(nAtomsInCol),stat=alloc_stat)
303     if (alloc_stat /= 0) then
304     status = -1
305     return
306     end if
307     if (allocated(rho_col)) deallocate(rho_col)
308     allocate(rho_col(nAtomsInCol),stat=alloc_stat)
309     if (alloc_stat /= 0) then
310     status = -1
311     return
312     end if
313     if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
314     allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
315     if (alloc_stat /= 0) then
316     status = -1
317     return
318     end if
319 gezelter 2204
320 gezelter 1608 #endif
321    
322 gezelter 2717 end subroutine init_EAM_FF
323 gezelter 1608
324 gezelter 2717 subroutine setCutoffEAM(rcut)
325 gezelter 1608 real(kind=dp) :: rcut
326     EAM_rcut = rcut
327     end subroutine setCutoffEAM
328    
329     subroutine clean_EAM()
330 gezelter 2204
331     ! clean non-IS_MPI first
332 gezelter 1608 frho = 0.0_dp
333     rho = 0.0_dp
334     dfrhodrho = 0.0_dp
335 gezelter 2204 ! clean MPI if needed
336 gezelter 1608 #ifdef IS_MPI
337     frho_row = 0.0_dp
338     frho_col = 0.0_dp
339     rho_row = 0.0_dp
340     rho_col = 0.0_dp
341     rho_tmp = 0.0_dp
342     dfrhodrho_row = 0.0_dp
343     dfrhodrho_col = 0.0_dp
344     #endif
345     end subroutine clean_EAM
346    
347     subroutine deallocate_EAMType(thisEAMType)
348     type (EAMtype), pointer :: thisEAMType
349    
350 gezelter 2717 call deleteSpline(thisEAMType%F)
351     call deleteSpline(thisEAMType%rho)
352     call deleteSpline(thisEAMType%phi)
353     call deleteSpline(thisEAMType%Z)
354 gezelter 2204
355 gezelter 1608 end subroutine deallocate_EAMType
356    
357 gezelter 2204 !! Calculates rho_r
358 gezelter 1608 subroutine calc_eam_prepair_rho(atom1,atom2,d,r,rijsq)
359 gezelter 2717 integer :: atom1, atom2
360 gezelter 1608 real(kind = dp), dimension(3) :: d
361     real(kind = dp), intent(inout) :: r
362     real(kind = dp), intent(inout) :: rijsq
363     ! value of electron density rho do to atom i at atom j
364     real(kind = dp) :: rho_i_at_j
365     ! value of electron density rho do to atom j at atom i
366     real(kind = dp) :: rho_j_at_i
367     integer :: eam_err
368 chuckv 2291
369 gezelter 2717 integer :: atid1, atid2 ! Global atid
370 chuckv 2291 integer :: myid_atom1 ! EAM atid
371     integer :: myid_atom2
372 gezelter 2204
373     ! check to see if we need to be cleaned at the start of a force loop
374 gezelter 1608
375     #ifdef IS_MPI
376 chuckv 2291 Atid1 = Atid_row(Atom1)
377     Atid2 = Atid_col(Atom2)
378 gezelter 1608 #else
379 chuckv 2291 Atid1 = Atid(Atom1)
380     Atid2 = Atid(Atom2)
381 gezelter 1608 #endif
382    
383 chuckv 2291 Myid_atom1 = Eamlist%atidtoeamtype(Atid1)
384     Myid_atom2 = Eamlist%atidtoeamtype(Atid2)
385    
386 gezelter 1608 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
387    
388 gezelter 2717 call lookupUniformSpline(EAMList%EAMParams(myid_atom1)%rho, r, &
389     rho_i_at_j)
390 gezelter 1608
391     #ifdef IS_MPI
392     rho_col(atom2) = rho_col(atom2) + rho_i_at_j
393     #else
394     rho(atom2) = rho(atom2) + rho_i_at_j
395     #endif
396 gezelter 2204 endif
397 gezelter 1608
398 gezelter 2204 if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then
399 gezelter 1608
400 gezelter 2717 call lookupUniformSpline(EAMList%EAMParams(myid_atom2)%rho, r, &
401     rho_j_at_i)
402 gezelter 2204
403 gezelter 1608 #ifdef IS_MPI
404 gezelter 2204 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
405 gezelter 1608 #else
406 gezelter 2204 rho(atom1) = rho(atom1) + rho_j_at_i
407 gezelter 1608 #endif
408 gezelter 2204 endif
409 gezelter 1608
410     end subroutine calc_eam_prepair_rho
411    
412    
413     !! Calculate the functional F(rho) for all local atoms
414 gezelter 2717 subroutine calc_eam_preforce_Frho(nlocal, pot)
415 gezelter 1608 integer :: nlocal
416     real(kind=dp) :: pot
417 gezelter 2717 integer :: i, j
418 gezelter 1608 integer :: atom
419 gezelter 2717 real(kind=dp) :: U,U1
420 gezelter 1608 integer :: atype1
421 gezelter 2717 integer :: me, atid1
422 gezelter 1608
423     cleanme = .true.
424 gezelter 2717 !! Scatter the electron density from pre-pair calculation back to
425     !! local atoms
426 gezelter 1608 #ifdef IS_MPI
427     call scatter(rho_row,rho,plan_atom_row,eam_err)
428     if (eam_err /= 0 ) then
429 gezelter 2204 write(errMsg,*) " Error scattering rho_row into rho"
430     call handleError(RoutineName,errMesg)
431     endif
432 gezelter 1608 call scatter(rho_col,rho_tmp,plan_atom_col,eam_err)
433     if (eam_err /= 0 ) then
434 gezelter 2204 write(errMsg,*) " Error scattering rho_col into rho"
435     call handleError(RoutineName,errMesg)
436     endif
437 gezelter 1608
438 gezelter 2204 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
439 gezelter 1608 #endif
440    
441 gezelter 2204 !! Calculate F(rho) and derivative
442 gezelter 1608 do atom = 1, nlocal
443 chuckv 2291 atid1 = atid(atom)
444     me = eamList%atidToEAMtype(atid1)
445 gezelter 1608
446 gezelter 2717 call lookupUniformSpline1d(EAMList%EAMParams(me)%F, rho(atom), &
447     u, u1)
448    
449 gezelter 1608 frho(atom) = u
450     dfrhodrho(atom) = u1
451     pot = pot + u
452     enddo
453    
454     #ifdef IS_MPI
455     !! communicate f(rho) and derivatives back into row and column arrays
456     call gather(frho,frho_row,plan_atom_row, eam_err)
457     if (eam_err /= 0) then
458     call handleError("cal_eam_forces()","MPI gather frho_row failure")
459     endif
460     call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, eam_err)
461     if (eam_err /= 0) then
462     call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
463     endif
464     call gather(frho,frho_col,plan_atom_col, eam_err)
465     if (eam_err /= 0) then
466     call handleError("cal_eam_forces()","MPI gather frho_col failure")
467     endif
468     call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, eam_err)
469     if (eam_err /= 0) then
470     call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
471     endif
472     #endif
473    
474 gezelter 2204
475 gezelter 1608 end subroutine calc_eam_preforce_Frho
476 gezelter 2717
477 gezelter 1608 !! Does EAM pairwise Force calculation.
478     subroutine do_eam_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
479     pot, f, do_pot)
480     !Arguments
481     integer, intent(in) :: atom1, atom2
482     real( kind = dp ), intent(in) :: rij, r2
483     real( kind = dp ) :: pot, sw, vpair
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 gezelter 2204
490 gezelter 2717 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 gezelter 1608 real( kind = dp ) :: dudr
495 gezelter 2717 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 gezelter 1608
500 gezelter 2717 integer :: id1, id2
501 gezelter 1608 integer :: mytype_atom1
502     integer :: mytype_atom2
503 gezelter 2717 integer :: atid1, atid2
504 gezelter 1608
505     phab = 0.0E0_DP
506     dvpdr = 0.0E0_DP
507    
508     if (rij .lt. EAM_rcut) then
509    
510     #ifdef IS_MPI
511 chuckv 2291 atid1 = atid_row(atom1)
512     atid2 = atid_col(atom2)
513 gezelter 1608 #else
514 chuckv 2291 atid1 = atid(atom1)
515     atid2 = atid(atom2)
516 gezelter 1608 #endif
517 chuckv 2291
518     mytype_atom1 = EAMList%atidToEAMType(atid1)
519     mytype_atom2 = EAMList%atidTOEAMType(atid2)
520    
521    
522 gezelter 1608 ! get cutoff for atom 1
523     rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
524     ! get type specific cutoff for atom 2
525     rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
526 gezelter 2204
527 gezelter 1608 drdx = d(1)/rij
528     drdy = d(2)/rij
529     drdz = d(3)/rij
530 gezelter 2204
531 gezelter 1608 if (rij.lt.rci) then
532 gezelter 2717
533     ! Calculate rho and drho for atom1
534    
535     call lookupUniformSpline1d(EAMList%EAMParams(mytype_atom1)%rho, &
536     rij, rha, drha)
537    
538     ! Calculate Phi(r) for atom1.
539    
540     call lookupUniformSpline1d(EAMList%EAMParams(mytype_atom1)%phi, &
541     rij, pha, dpha)
542    
543 gezelter 1608 endif
544    
545     if (rij.lt.rcj) then
546 gezelter 2204
547 gezelter 2717 ! Calculate rho and drho for atom2
548    
549     call lookupUniformSpline1d(EAMList%EAMParams(mytype_atom2)%rho, &
550     rij, rhb, drhb)
551    
552     ! Calculate Phi(r) for atom2.
553    
554     call lookupUniformSpline1d(EAMList%EAMParams(mytype_atom2)%phi, &
555     rij, phb, dphb)
556    
557 gezelter 1608 endif
558    
559     if (rij.lt.rci) then
560     phab = phab + 0.5E0_DP*(rhb/rha)*pha
561     dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
562     pha*((drhb/rha) - (rhb*drha/rha/rha)))
563     endif
564 gezelter 2204
565 gezelter 1608 if (rij.lt.rcj) then
566     phab = phab + 0.5E0_DP*(rha/rhb)*phb
567     dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
568     phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
569     endif
570 gezelter 2204
571 gezelter 1608 drhoidr = drha
572     drhojdr = drhb
573    
574     #ifdef IS_MPI
575     dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) &
576     + dvpdr
577    
578     #else
579     dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) &
580     + dvpdr
581     #endif
582    
583     fx = dudr * drdx
584     fy = dudr * drdy
585     fz = dudr * drdz
586    
587    
588     #ifdef IS_MPI
589     if (do_pot) then
590 gezelter 2361 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + phab*0.5
591     pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + phab*0.5
592 gezelter 1608 end if
593    
594     f_Row(1,atom1) = f_Row(1,atom1) + fx
595     f_Row(2,atom1) = f_Row(2,atom1) + fy
596     f_Row(3,atom1) = f_Row(3,atom1) + fz
597 gezelter 2204
598 gezelter 1608 f_Col(1,atom2) = f_Col(1,atom2) - fx
599     f_Col(2,atom2) = f_Col(2,atom2) - fy
600     f_Col(3,atom2) = f_Col(3,atom2) - fz
601     #else
602    
603     if(do_pot) then
604     pot = pot + phab
605     end if
606    
607     f(1,atom1) = f(1,atom1) + fx
608     f(2,atom1) = f(2,atom1) + fy
609     f(3,atom1) = f(3,atom1) + fz
610 gezelter 2204
611 gezelter 1608 f(1,atom2) = f(1,atom2) - fx
612     f(2,atom2) = f(2,atom2) - fy
613     f(3,atom2) = f(3,atom2) - fz
614     #endif
615 gezelter 2204
616 gezelter 1608 vpair = vpair + phab
617 gezelter 2717
618 gezelter 1608 #ifdef IS_MPI
619     id1 = AtomRowToGlobal(atom1)
620     id2 = AtomColToGlobal(atom2)
621     #else
622     id1 = atom1
623     id2 = atom2
624     #endif
625 gezelter 2204
626 gezelter 1608 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
627 gezelter 2204
628 gezelter 1608 fpair(1) = fpair(1) + fx
629     fpair(2) = fpair(2) + fy
630     fpair(3) = fpair(3) + fz
631 gezelter 2204
632 gezelter 1608 endif
633 gezelter 2204 endif
634 gezelter 1608 end subroutine do_eam_pair
635    
636     end module eam