ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 2530
Committed: Fri Dec 30 00:18:28 2005 UTC (18 years, 6 months ago) by chuckv
File size: 19029 byte(s)
Log Message:
Sutton-Chen bug fixes. Almost there...

File Contents

# User Rev Content
1 chuckv 2401 !!
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     !! Impliments Sutton-Chen Metallic Potential
43     !! See A.P.SUTTON and J.CHEN,PHIL MAG LETT 61,139-146,1990
44    
45    
46     module suttonchen
47     use simulation
48     use force_globals
49     use status
50     use atype_module
51     use vector_class
52 chuckv 2497 use fForceOptions
53 chuckv 2401 #ifdef IS_MPI
54     use mpiSimulation
55     #endif
56     implicit none
57     PRIVATE
58     #define __FORTRAN90
59     #include "UseTheForce/DarkSide/fInteractionMap.h"
60    
61     INTEGER, PARAMETER :: DP = selected_real_kind(15)
62    
63     logical, save :: SC_FF_initialized = .false.
64     integer, save :: SC_Mixing_Policy
65     real(kind = dp), save :: SC_rcut
66     logical, save :: haveRcut = .false.
67 chuckv 2427 logical, save :: haveMixingMap = .false.
68     logical, save :: useGeometricDistanceMixing = .false.
69 chuckv 2401
70 chuckv 2427
71    
72    
73 chuckv 2401 character(len = statusMsgSize) :: errMesg
74 chuckv 2413 integer :: sc_err
75 chuckv 2401
76     character(len = 200) :: errMsg
77     character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE"
78     !! Logical that determines if eam arrays should be zeroed
79     logical :: cleanme = .true.
80     logical :: nmflag = .false.
81    
82    
83     type, private :: SCtype
84     integer :: atid
85 chuckv 2406 real(kind=dp) :: c
86 chuckv 2401 real(kind=dp) :: m
87     real(kind=dp) :: n
88     real(kind=dp) :: alpha
89     real(kind=dp) :: epsilon
90 chuckv 2427 real(kind=dp) :: sc_rcut
91 chuckv 2401 end type SCtype
92    
93    
94     !! Arrays for derivatives used in force calculation
95 chuckv 2412 real( kind = dp), dimension(:), allocatable :: rho
96 chuckv 2401 real( kind = dp), dimension(:), allocatable :: frho
97     real( kind = dp), dimension(:), allocatable :: dfrhodrho
98    
99    
100 chuckv 2427
101 chuckv 2401 !! Arrays for MPI storage
102     #ifdef IS_MPI
103     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
104     real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
105     real( kind = dp),save, dimension(:), allocatable :: frho_row
106     real( kind = dp),save, dimension(:), allocatable :: frho_col
107     real( kind = dp),save, dimension(:), allocatable :: rho_row
108     real( kind = dp),save, dimension(:), allocatable :: rho_col
109     real( kind = dp),save, dimension(:), allocatable :: rho_tmp
110     #endif
111    
112     type, private :: SCTypeList
113     integer :: nSCTypes = 0
114     integer :: currentSCtype = 0
115    
116     type (SCtype), pointer :: SCtypes(:) => null()
117     integer, pointer :: atidToSCtype(:) => null()
118     end type SCTypeList
119    
120     type (SCTypeList), save :: SCList
121    
122    
123    
124 chuckv 2406
125     type :: MixParameters
126     real(kind=DP) :: alpha
127     real(kind=DP) :: epsilon
128 chuckv 2427 real(kind=DP) :: m
129     real(Kind=DP) :: n
130     real(kind=DP) :: vpair_pot
131 chuckv 2406 real(kind=dp) :: rCut
132     logical :: rCutWasSet = .false.
133 chuckv 2427
134 chuckv 2406 end type MixParameters
135    
136     type(MixParameters), dimension(:,:), allocatable :: MixingMap
137    
138    
139    
140 chuckv 2401 public :: setCutoffSC
141     public :: do_SC_pair
142     public :: newSCtype
143     public :: calc_SC_prepair_rho
144 chuckv 2434 public :: calc_SC_preforce_Frho
145 chuckv 2401 public :: clean_SC
146     public :: destroySCtypes
147     public :: getSCCut
148 chuckv 2427 ! public :: setSCDefaultCutoff
149     ! public :: setSCUniformCutoff
150 chuckv 2497
151 chuckv 2401
152     contains
153    
154    
155 chuckv 2427 subroutine newSCtype(c_ident,c,m,n,alpha,epsilon,status)
156     real (kind = dp ) :: c ! Density Scaling
157     real (kind = dp ) :: m ! Density Exponent
158     real (kind = dp ) :: n ! Pair Potential Exponent
159     real (kind = dp ) :: alpha ! Length Scaling
160     real (kind = dp ) :: epsilon ! Energy Scaling
161    
162    
163 chuckv 2401 integer :: c_ident
164     integer :: status
165    
166 chuckv 2412 integer :: nAtypes,nSCTypes,myATID
167 chuckv 2401 integer :: maxVals
168     integer :: alloc_stat
169     integer :: current
170     integer,pointer :: Matchlist(:) => null()
171    
172     status = 0
173    
174    
175     !! Assume that atypes has already been set and get the total number of types in atypes
176    
177    
178     ! check to see if this is the first time into
179 chuckv 2427 if (.not.associated(SCList%SCTypes)) then
180 chuckv 2530 call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList)
181 chuckv 2427 SCList%nSCtypes = nSCtypes
182     allocate(SCList%SCTypes(nSCTypes))
183 chuckv 2401 nAtypes = getSize(atypes)
184 chuckv 2427 allocate(SCList%atidToSCType(nAtypes))
185 chuckv 2401 end if
186    
187 chuckv 2427 SCList%currentSCType = SCList%currentSCType + 1
188     current = SCList%currentSCType
189 chuckv 2401
190     myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
191 chuckv 2427 SCList%atidToSCType(myATID) = current
192 chuckv 2401
193    
194    
195 chuckv 2427 SCList%SCTypes(current)%atid = c_ident
196     SCList%SCTypes(current)%alpha = alpha
197     SCList%SCTypes(current)%c = c
198     SCList%SCTypes(current)%m = m
199     SCList%SCTypes(current)%n = n
200     SCList%SCTypes(current)%epsilon = epsilon
201 chuckv 2412 end subroutine newSCtype
202 chuckv 2401
203 chuckv 2412
204     subroutine destroySCTypes()
205 chuckv 2427 if (associated(SCList%SCtypes)) then
206     deallocate(SCList%SCTypes)
207     SCList%SCTypes=>null()
208     end if
209     if (associated(SCList%atidToSCtype)) then
210     deallocate(SCList%atidToSCtype)
211     SCList%atidToSCtype=>null()
212     end if
213 chuckv 2401
214    
215 chuckv 2412 end subroutine destroySCTypes
216 chuckv 2401
217    
218 chuckv 2412
219     function getSCCut(atomID) result(cutValue)
220 chuckv 2401 integer, intent(in) :: atomID
221 chuckv 2427 integer :: scID
222 chuckv 2401 real(kind=dp) :: cutValue
223    
224 chuckv 2427 scID = SCList%atidToSCType(atomID)
225 chuckv 2530 cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha
226 chuckv 2412 end function getSCCut
227 chuckv 2401
228    
229    
230 chuckv 2412
231     subroutine createMixingMap()
232     integer :: nSCtypes, i, j
233     real ( kind = dp ) :: e1, e2,m1,m2,alpha1,alpha2,n1,n2
234     real ( kind = dp ) :: rcut6, tp6, tp12
235     logical :: isSoftCore1, isSoftCore2, doShift
236    
237     if (SCList%currentSCtype == 0) then
238     call handleError("SuttonChen", "No members in SCMap")
239 chuckv 2401 return
240     end if
241    
242 chuckv 2412 nSCtypes = SCList%nSCtypes
243 chuckv 2401
244 chuckv 2412 if (.not. allocated(MixingMap)) then
245     allocate(MixingMap(nSCtypes, nSCtypes))
246     endif
247 chuckv 2497 useGeometricDistanceMixing = usesGeometricDistanceMixing()
248 chuckv 2412 do i = 1, nSCtypes
249 chuckv 2401
250 chuckv 2412 e1 = SCList%SCtypes(i)%epsilon
251     m1 = SCList%SCtypes(i)%m
252     n1 = SCList%SCtypes(i)%n
253     alpha1 = SCList%SCtypes(i)%alpha
254 chuckv 2401
255 chuckv 2412 do j = i, nSCtypes
256    
257 chuckv 2401
258 chuckv 2412 e2 = SCList%SCtypes(j)%epsilon
259     m2 = SCList%SCtypes(j)%m
260     n2 = SCList%SCtypes(j)%n
261     alpha2 = SCList%SCtypes(j)%alpha
262 chuckv 2401
263 chuckv 2427 if (useGeometricDistanceMixing) then
264     MixingMap(i,j)%alpha = sqrt(alpha1 * alpha2) !SC formulation
265     else
266     MixingMap(i,j)%alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation
267     endif
268    
269     MixingMap(i,j)%epsilon = sqrt(e1 * e2)
270 chuckv 2412 MixingMap(i,j)%m = 0.5_dp*(m1+m2)
271     MixingMap(i,j)%n = 0.5_dp*(n1+n2)
272     MixingMap(i,j)%alpha = 0.5_dp*(alpha1+alpha2)
273     MixingMap(i,j)%rcut = 2.0_dp *MixingMap(i,j)%alpha
274 chuckv 2427 MixingMap(i,j)%vpair_pot = MixingMap(i,j)%epsilon* &
275     (MixingMap(i,j)%alpha/MixingMap(i,j)%rcut)**MixingMap(i,j)%n
276     if (i.ne.j) then
277     MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
278     MixingMap(j,i)%m = MixingMap(i,j)%m
279     MixingMap(j,i)%n = MixingMap(i,j)%n
280     MixingMap(j,i)%alpha = MixingMap(i,j)%alpha
281     MixingMap(j,i)%rcut = MixingMap(i,j)%rcut
282     MixingMap(j,i)%vpair_pot = MixingMap(i,j)%vpair_pot
283     endif
284 chuckv 2412 enddo
285 chuckv 2401 enddo
286 chuckv 2412
287     haveMixingMap = .true.
288    
289     end subroutine createMixingMap
290    
291 chuckv 2401
292     !! routine checks to see if array is allocated, deallocates array if allocated
293     !! and then creates the array to the required size
294 chuckv 2412 subroutine allocateSC(status)
295 chuckv 2401 integer, intent(out) :: status
296    
297     #ifdef IS_MPI
298     integer :: nAtomsInRow
299     integer :: nAtomsInCol
300     #endif
301     integer :: alloc_stat
302    
303    
304     status = 0
305     #ifdef IS_MPI
306     nAtomsInRow = getNatomsInRow(plan_atom_row)
307     nAtomsInCol = getNatomsInCol(plan_atom_col)
308     #endif
309    
310 chuckv 2412
311    
312 chuckv 2401 if (allocated(frho)) deallocate(frho)
313     allocate(frho(nlocal),stat=alloc_stat)
314 chuckv 2412 if (alloc_stat /= 0) then
315 chuckv 2401 status = -1
316     return
317     end if
318 chuckv 2412
319 chuckv 2401 if (allocated(rho)) deallocate(rho)
320     allocate(rho(nlocal),stat=alloc_stat)
321     if (alloc_stat /= 0) then
322     status = -1
323     return
324     end if
325    
326     if (allocated(dfrhodrho)) deallocate(dfrhodrho)
327     allocate(dfrhodrho(nlocal),stat=alloc_stat)
328     if (alloc_stat /= 0) then
329     status = -1
330     return
331     end if
332    
333     #ifdef IS_MPI
334    
335     if (allocated(rho_tmp)) deallocate(rho_tmp)
336     allocate(rho_tmp(nlocal),stat=alloc_stat)
337     if (alloc_stat /= 0) then
338     status = -1
339     return
340     end if
341    
342    
343     if (allocated(frho_row)) deallocate(frho_row)
344     allocate(frho_row(nAtomsInRow),stat=alloc_stat)
345     if (alloc_stat /= 0) then
346     status = -1
347     return
348     end if
349     if (allocated(rho_row)) deallocate(rho_row)
350     allocate(rho_row(nAtomsInRow),stat=alloc_stat)
351     if (alloc_stat /= 0) then
352     status = -1
353     return
354     end if
355     if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
356     allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
357     if (alloc_stat /= 0) then
358     status = -1
359     return
360     end if
361    
362    
363     ! Now do column arrays
364    
365     if (allocated(frho_col)) deallocate(frho_col)
366     allocate(frho_col(nAtomsInCol),stat=alloc_stat)
367     if (alloc_stat /= 0) then
368     status = -1
369     return
370     end if
371     if (allocated(rho_col)) deallocate(rho_col)
372     allocate(rho_col(nAtomsInCol),stat=alloc_stat)
373     if (alloc_stat /= 0) then
374     status = -1
375     return
376     end if
377     if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
378     allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
379     if (alloc_stat /= 0) then
380     status = -1
381     return
382     end if
383    
384     #endif
385    
386 chuckv 2412 end subroutine allocateSC
387 chuckv 2401
388     !! C sets rcut to be the largest cutoff of any atype
389     !! present in this simulation. Doesn't include all atypes
390     !! sim knows about, just those in the simulation.
391 chuckv 2412 subroutine setCutoffSC(rcut, status)
392 chuckv 2401 real(kind=dp) :: rcut
393     integer :: status
394     status = 0
395    
396 chuckv 2412 SC_rcut = rcut
397 chuckv 2401
398 chuckv 2412 end subroutine setCutoffSC
399 chuckv 2401
400 chuckv 2427 subroutine clean_SC()
401    
402 chuckv 2401 ! clean non-IS_MPI first
403     frho = 0.0_dp
404     rho = 0.0_dp
405     dfrhodrho = 0.0_dp
406     ! clean MPI if needed
407     #ifdef IS_MPI
408     frho_row = 0.0_dp
409     frho_col = 0.0_dp
410     rho_row = 0.0_dp
411     rho_col = 0.0_dp
412     rho_tmp = 0.0_dp
413     dfrhodrho_row = 0.0_dp
414     dfrhodrho_col = 0.0_dp
415     #endif
416 chuckv 2427 end subroutine clean_SC
417 chuckv 2401
418    
419    
420     !! Calculates rho_r
421 gezelter 2461 subroutine calc_sc_prepair_rho(atom1,atom2,d,r,rijsq, rcut)
422 chuckv 2401 integer :: atom1,atom2
423     real(kind = dp), dimension(3) :: d
424     real(kind = dp), intent(inout) :: r
425 gezelter 2461 real(kind = dp), intent(inout) :: rijsq, rcut
426 chuckv 2401 ! value of electron density rho do to atom i at atom j
427     real(kind = dp) :: rho_i_at_j
428     ! value of electron density rho do to atom j at atom i
429     real(kind = dp) :: rho_j_at_i
430    
431     ! we don't use the derivatives, dummy variables
432 chuckv 2427 real( kind = dp) :: drho
433 chuckv 2413 integer :: sc_err
434 chuckv 2401
435     integer :: atid1,atid2 ! Global atid
436     integer :: myid_atom1 ! EAM atid
437     integer :: myid_atom2
438    
439    
440     ! check to see if we need to be cleaned at the start of a force loop
441    
442    
443    
444    
445     #ifdef IS_MPI
446     Atid1 = Atid_row(Atom1)
447     Atid2 = Atid_col(Atom2)
448     #else
449     Atid1 = Atid(Atom1)
450     Atid2 = Atid(Atom2)
451     #endif
452    
453 chuckv 2427 Myid_atom1 = SCList%atidtoSCtype(Atid1)
454     Myid_atom2 = SCList%atidtoSCtype(Atid2)
455 chuckv 2401
456    
457    
458 chuckv 2412 rho_i_at_j = (MixingMap(Myid_atom1,Myid_atom2)%alpha/r)&
459     **MixingMap(Myid_atom1,Myid_atom2)%m
460     rho_j_at_i = rho_i_at_j
461 chuckv 2401
462     #ifdef IS_MPI
463     rho_col(atom2) = rho_col(atom2) + rho_i_at_j
464 chuckv 2412 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
465 chuckv 2401 #else
466     rho(atom2) = rho(atom2) + rho_i_at_j
467     rho(atom1) = rho(atom1) + rho_j_at_i
468     #endif
469    
470 chuckv 2412 end subroutine calc_sc_prepair_rho
471 chuckv 2401
472    
473 chuckv 2413 !! Calculate the rho_a for all local atoms
474 chuckv 2427 subroutine calc_sc_preforce_Frho(nlocal,pot)
475 chuckv 2401 integer :: nlocal
476     real(kind=dp) :: pot
477     integer :: i,j
478     integer :: atom
479     real(kind=dp) :: U,U1,U2
480     integer :: atype1
481 chuckv 2427 integer :: atid1
482     integer :: myid
483 chuckv 2401
484    
485     cleanme = .true.
486     !! Scatter the electron density from pre-pair calculation back to local atoms
487     #ifdef IS_MPI
488 chuckv 2413 call scatter(rho_row,rho,plan_atom_row,sc_err)
489     if (sc_err /= 0 ) then
490 chuckv 2401 write(errMsg,*) " Error scattering rho_row into rho"
491     call handleError(RoutineName,errMesg)
492     endif
493 chuckv 2413 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
494     if (sc_err /= 0 ) then
495 chuckv 2401 write(errMsg,*) " Error scattering rho_col into rho"
496     call handleError(RoutineName,errMesg)
497 chuckv 2412
498 chuckv 2401 endif
499    
500     rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
501     #endif
502    
503    
504    
505 chuckv 2412 !! Calculate F(rho) and derivative
506 chuckv 2401 do atom = 1, nlocal
507 chuckv 2427 Myid = SCList%atidtoSctype(Atid(atom))
508     frho(atom) = -SCList%SCTypes(Myid)%c * &
509     SCList%SCTypes(Myid)%epsilon * sqrt(rho(i))
510    
511 chuckv 2412 dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom)
512 chuckv 2401 pot = pot + u
513     enddo
514    
515 chuckv 2430 #ifdef IS_MPI
516 chuckv 2401 !! communicate f(rho) and derivatives back into row and column arrays
517 chuckv 2413 call gather(frho,frho_row,plan_atom_row, sc_err)
518     if (sc_err /= 0) then
519 chuckv 2401 call handleError("cal_eam_forces()","MPI gather frho_row failure")
520     endif
521 chuckv 2413 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
522     if (sc_err /= 0) then
523 chuckv 2401 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
524     endif
525 chuckv 2413 call gather(frho,frho_col,plan_atom_col, sc_err)
526     if (sc_err /= 0) then
527 chuckv 2401 call handleError("cal_eam_forces()","MPI gather frho_col failure")
528     endif
529 chuckv 2413 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
530     if (sc_err /= 0) then
531 chuckv 2401 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
532     endif
533     #endif
534 chuckv 2412
535    
536 chuckv 2427 end subroutine calc_sc_preforce_Frho
537 chuckv 2401
538    
539 chuckv 2412 !! Does Sutton-Chen pairwise Force calculation.
540 gezelter 2461 subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
541 chuckv 2401 pot, f, do_pot)
542     !Arguments
543     integer, intent(in) :: atom1, atom2
544 gezelter 2461 real( kind = dp ), intent(in) :: rij, r2, rcut
545 chuckv 2401 real( kind = dp ) :: pot, sw, vpair
546     real( kind = dp ), dimension(3,nLocal) :: f
547     real( kind = dp ), intent(in), dimension(3) :: d
548     real( kind = dp ), intent(inout), dimension(3) :: fpair
549    
550     logical, intent(in) :: do_pot
551    
552     real( kind = dp ) :: drdx,drdy,drdz
553 chuckv 2427 real( kind = dp ) :: dvpdr
554     real( kind = dp ) :: drhodr
555 chuckv 2401 real( kind = dp ) :: dudr
556 chuckv 2413 real( kind = dp ) :: rcij
557 chuckv 2401 real( kind = dp ) :: drhoidr,drhojdr
558     real( kind = dp ) :: Fx,Fy,Fz
559     real( kind = dp ) :: r,d2pha,phb,d2phb
560 chuckv 2427 real( kind = dp ) :: pot_temp,vptmp
561     real( kind = dp ) :: epsilonij,aij,nij,mij,vcij
562 chuckv 2401 integer :: id1,id2
563     integer :: mytype_atom1
564     integer :: mytype_atom2
565     integer :: atid1,atid2
566     !Local Variables
567    
568     ! write(*,*) "Frho: ", Frho(atom1)
569    
570 chuckv 2427
571 chuckv 2401 dvpdr = 0.0E0_DP
572    
573    
574     #ifdef IS_MPI
575     atid1 = atid_row(atom1)
576     atid2 = atid_col(atom2)
577     #else
578     atid1 = atid(atom1)
579     atid2 = atid(atom2)
580     #endif
581    
582 chuckv 2427 mytype_atom1 = SCList%atidToSCType(atid1)
583     mytype_atom2 = SCList%atidTOSCType(atid2)
584 chuckv 2401
585     drdx = d(1)/rij
586     drdy = d(2)/rij
587     drdz = d(3)/rij
588    
589 chuckv 2413
590     epsilonij = MixingMap(mytype_atom1,mytype_atom2)%epsilon
591     aij = MixingMap(mytype_atom1,mytype_atom2)%alpha
592     nij = MixingMap(mytype_atom1,mytype_atom2)%n
593     mij = MixingMap(mytype_atom1,mytype_atom2)%m
594     vcij = MixingMap(mytype_atom1,mytype_atom2)%vpair_pot
595 chuckv 2401
596 chuckv 2427 vptmp = epsilonij*((aij/r)**nij)
597 chuckv 2401
598    
599 chuckv 2413 dvpdr = -nij*vptmp/r
600     drhodr = -mij*((aij/r)**mij)/r
601 chuckv 2427
602 chuckv 2413
603 chuckv 2427 dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) &
604 chuckv 2413 + dvpdr
605 chuckv 2427
606     pot_temp = vptmp + vcij
607 chuckv 2413
608 chuckv 2401
609     #ifdef IS_MPI
610 chuckv 2413 dudr = drhodr*(dfrhodrho_row(atom1)+dfrhodrho_col(atom2)) &
611 chuckv 2401 + dvpdr
612    
613     #else
614 chuckv 2413 dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) &
615 chuckv 2401 + dvpdr
616     #endif
617    
618 chuckv 2413
619 chuckv 2401 fx = dudr * drdx
620     fy = dudr * drdy
621     fz = dudr * drdz
622    
623    
624     #ifdef IS_MPI
625     if (do_pot) then
626 chuckv 2427 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5
627     pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5
628 chuckv 2401 end if
629    
630     f_Row(1,atom1) = f_Row(1,atom1) + fx
631     f_Row(2,atom1) = f_Row(2,atom1) + fy
632     f_Row(3,atom1) = f_Row(3,atom1) + fz
633    
634     f_Col(1,atom2) = f_Col(1,atom2) - fx
635     f_Col(2,atom2) = f_Col(2,atom2) - fy
636     f_Col(3,atom2) = f_Col(3,atom2) - fz
637     #else
638    
639     if(do_pot) then
640 chuckv 2427 pot = pot + pot_temp
641 chuckv 2401 end if
642    
643     f(1,atom1) = f(1,atom1) + fx
644     f(2,atom1) = f(2,atom1) + fy
645     f(3,atom1) = f(3,atom1) + fz
646    
647     f(1,atom2) = f(1,atom2) - fx
648     f(2,atom2) = f(2,atom2) - fy
649     f(3,atom2) = f(3,atom2) - fz
650     #endif
651    
652 chuckv 2427
653 chuckv 2401 #ifdef IS_MPI
654     id1 = AtomRowToGlobal(atom1)
655     id2 = AtomColToGlobal(atom2)
656     #else
657     id1 = atom1
658     id2 = atom2
659     #endif
660    
661     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
662    
663     fpair(1) = fpair(1) + fx
664     fpair(2) = fpair(2) + fy
665     fpair(3) = fpair(3) + fz
666    
667     endif
668    
669    
670 chuckv 2427 end subroutine do_sc_pair
671 chuckv 2401
672    
673    
674 chuckv 2412 end module suttonchen