ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 939
Committed: Thu Apr 20 18:24:24 2006 UTC (19 years, 1 month ago) by gezelter
File size: 18824 byte(s)
Log Message:
Complete rewrite of spline code and everything that uses it.

File Contents

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