ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 762
Committed: Mon Nov 21 22:59:02 2005 UTC (19 years, 6 months ago) by gezelter
File size: 19138 byte(s)
Log Message:
Cutoff mixing fixes have been made.

File Contents

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