1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
!! publication of scientific results based in part on use of the |
11 |
!! program. An acceptable form of acknowledgement is citation of |
12 |
!! the article in which the program was described (Matthew |
13 |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
!! |
18 |
!! 2. Redistributions of source code must retain the above copyright |
19 |
!! notice, this list of conditions and the following disclaimer. |
20 |
!! |
21 |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
!! notice, this list of conditions and the following disclaimer in the |
23 |
!! documentation and/or other materials provided with the |
24 |
!! distribution. |
25 |
!! |
26 |
!! This software is provided "AS IS," without a warranty of any |
27 |
!! kind. All express or implied conditions, representations and |
28 |
!! warranties, including any implied warranty of merchantability, |
29 |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
!! be liable for any damages suffered by licensee as a result of |
32 |
!! using, modifying or distributing the software or its |
33 |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
!! damages, however caused and regardless of the theory of liability, |
37 |
!! arising out of the use of or inability to use software, even if the |
38 |
!! University of Notre Dame has been advised of the possibility of |
39 |
!! such damages. |
40 |
!! |
41 |
|
42 |
!! 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 definitions |
48 |
use simulation |
49 |
use force_globals |
50 |
use status |
51 |
use atype_module |
52 |
use vector_class |
53 |
use fForceOptions |
54 |
use interpolation |
55 |
#ifdef IS_MPI |
56 |
use mpiSimulation |
57 |
#endif |
58 |
implicit none |
59 |
PRIVATE |
60 |
#define __FORTRAN90 |
61 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
62 |
|
63 |
!! number of points for the spline approximations |
64 |
INTEGER, PARAMETER :: np = 3000 |
65 |
|
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 |
logical, save :: haveMixingMap = .false. |
71 |
logical, save :: useGeometricDistanceMixing = .false. |
72 |
logical, save :: cleanArrays = .true. |
73 |
logical, save :: arraysAllocated = .false. |
74 |
|
75 |
|
76 |
character(len = statusMsgSize) :: errMesg |
77 |
integer :: sc_err |
78 |
|
79 |
character(len = 200) :: errMsg |
80 |
character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE" |
81 |
|
82 |
type, private :: SCtype |
83 |
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 |
end type SCtype |
91 |
|
92 |
|
93 |
!! Arrays for derivatives used in force calculation |
94 |
real( kind = dp), dimension(:), allocatable :: rho |
95 |
real( kind = dp), dimension(:), allocatable :: frho |
96 |
real( kind = dp), dimension(:), allocatable :: dfrhodrho |
97 |
|
98 |
!! 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 |
integer :: currentSCtype = 0 |
112 |
type (SCtype), pointer :: SCtypes(:) => null() |
113 |
integer, pointer :: atidToSCtype(:) => null() |
114 |
end type SCTypeList |
115 |
|
116 |
type (SCTypeList), save :: SCList |
117 |
|
118 |
type:: MixParameters |
119 |
real(kind=DP) :: alpha |
120 |
real(kind=DP) :: epsilon |
121 |
real(kind=DP) :: m |
122 |
real(Kind=DP) :: n |
123 |
real(kind=dp) :: rCut |
124 |
real(kind=dp) :: vCut |
125 |
logical :: rCutWasSet = .false. |
126 |
type(cubicSpline) :: V |
127 |
type(cubicSpline) :: phi |
128 |
end type MixParameters |
129 |
|
130 |
type(MixParameters), dimension(:,:), allocatable :: MixingMap |
131 |
|
132 |
public :: setCutoffSC |
133 |
public :: do_SC_pair |
134 |
public :: newSCtype |
135 |
public :: calc_SC_prepair_rho |
136 |
public :: calc_SC_preforce_Frho |
137 |
public :: clean_SC |
138 |
public :: destroySCtypes |
139 |
public :: getSCCut |
140 |
! public :: setSCDefaultCutoff |
141 |
! public :: setSCUniformCutoff |
142 |
|
143 |
|
144 |
contains |
145 |
|
146 |
|
147 |
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 |
integer :: c_ident |
154 |
integer :: status |
155 |
integer :: nAtypes,nSCTypes,myATID |
156 |
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 |
if (.not.associated(SCList%SCTypes)) then |
169 |
call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList) |
170 |
SCList%nSCtypes = nSCtypes |
171 |
allocate(SCList%SCTypes(nSCTypes)) |
172 |
nAtypes = getSize(atypes) |
173 |
allocate(SCList%atidToSCType(nAtypes)) |
174 |
SCList%atidToSCType = -1 |
175 |
end if |
176 |
|
177 |
SCList%currentSCType = SCList%currentSCType + 1 |
178 |
current = SCList%currentSCType |
179 |
|
180 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
181 |
SCList%atidToSCType(myATID) = current |
182 |
|
183 |
|
184 |
SCList%SCTypes(current)%atid = c_ident |
185 |
SCList%SCTypes(current)%alpha = alpha |
186 |
SCList%SCTypes(current)%c = c |
187 |
SCList%SCTypes(current)%m = m |
188 |
SCList%SCTypes(current)%n = n |
189 |
SCList%SCTypes(current)%epsilon = epsilon |
190 |
end subroutine newSCtype |
191 |
|
192 |
|
193 |
subroutine destroySCTypes() |
194 |
if (associated(SCList%SCtypes)) then |
195 |
deallocate(SCList%SCTypes) |
196 |
SCList%SCTypes=>null() |
197 |
end if |
198 |
if (associated(SCList%atidToSCtype)) then |
199 |
deallocate(SCList%atidToSCtype) |
200 |
SCList%atidToSCtype=>null() |
201 |
end if |
202 |
! Reset Capacity |
203 |
SCList%nSCTypes = 0 |
204 |
SCList%currentSCtype=0 |
205 |
|
206 |
end subroutine destroySCTypes |
207 |
|
208 |
function getSCCut(atomID) result(cutValue) |
209 |
integer, intent(in) :: atomID |
210 |
integer :: scID |
211 |
real(kind=dp) :: cutValue |
212 |
|
213 |
scID = SCList%atidToSCType(atomID) |
214 |
cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha |
215 |
end function getSCCut |
216 |
|
217 |
subroutine createMixingMap() |
218 |
integer :: nSCtypes, i, j, k |
219 |
real ( kind = dp ) :: e1, e2, m1, m2, alpha1, alpha2, n1, n2 |
220 |
real ( kind = dp ) :: epsilon, m, n, alpha, rCut, vCut, dr, r |
221 |
real ( kind = dp ), dimension(np) :: rvals, vvals, phivals |
222 |
|
223 |
if (SCList%currentSCtype == 0) then |
224 |
call handleError("SuttonChen", "No members in SCMap") |
225 |
return |
226 |
end if |
227 |
|
228 |
nSCtypes = SCList%nSCtypes |
229 |
|
230 |
if (.not. allocated(MixingMap)) then |
231 |
allocate(MixingMap(nSCtypes, nSCtypes)) |
232 |
endif |
233 |
useGeometricDistanceMixing = usesGeometricDistanceMixing() |
234 |
do i = 1, nSCtypes |
235 |
|
236 |
e1 = SCList%SCtypes(i)%epsilon |
237 |
m1 = SCList%SCtypes(i)%m |
238 |
n1 = SCList%SCtypes(i)%n |
239 |
alpha1 = SCList%SCtypes(i)%alpha |
240 |
|
241 |
do j = 1, nSCtypes |
242 |
|
243 |
e2 = SCList%SCtypes(j)%epsilon |
244 |
m2 = SCList%SCtypes(j)%m |
245 |
n2 = SCList%SCtypes(j)%n |
246 |
alpha2 = SCList%SCtypes(j)%alpha |
247 |
|
248 |
if (useGeometricDistanceMixing) then |
249 |
alpha = sqrt(alpha1 * alpha2) !SC formulation |
250 |
else |
251 |
alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation |
252 |
endif |
253 |
rcut = 2.0_dp * alpha |
254 |
epsilon = sqrt(e1 * e2) |
255 |
m = 0.5_dp*(m1+m2) |
256 |
n = 0.5_dp*(n1+n2) |
257 |
|
258 |
dr = (rCut) / dble(np-1) |
259 |
rvals(1) = 0.0_dp |
260 |
vvals(1) = 0.0_dp |
261 |
phivals(1) = 0.0_dp |
262 |
|
263 |
do k = 2, np |
264 |
r = dble(k-1)*dr |
265 |
rvals(k) = r |
266 |
vvals(k) = epsilon*((alpha/r)**n) |
267 |
phivals(k) = (alpha/r)**m |
268 |
enddo |
269 |
|
270 |
vCut = epsilon*((alpha/rCut)**n) |
271 |
|
272 |
call newSpline(MixingMap(i,j)%V, rvals, vvals, .true.) |
273 |
call newSpline(MixingMap(i,j)%phi, rvals, phivals, .true.) |
274 |
|
275 |
MixingMap(i,j)%epsilon = epsilon |
276 |
MixingMap(i,j)%m = m |
277 |
MixingMap(i,j)%n = n |
278 |
MixingMap(i,j)%alpha = alpha |
279 |
MixingMap(i,j)%rCut = rcut |
280 |
MixingMap(i,j)%vCut = vCut |
281 |
enddo |
282 |
enddo |
283 |
|
284 |
haveMixingMap = .true. |
285 |
|
286 |
end subroutine createMixingMap |
287 |
|
288 |
|
289 |
!! routine checks to see if array is allocated, deallocates array if allocated |
290 |
!! and then creates the array to the required size |
291 |
subroutine allocateSC() |
292 |
integer :: status |
293 |
|
294 |
#ifdef IS_MPI |
295 |
integer :: nAtomsInRow |
296 |
integer :: nAtomsInCol |
297 |
#endif |
298 |
integer :: alloc_stat |
299 |
|
300 |
|
301 |
status = 0 |
302 |
#ifdef IS_MPI |
303 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
304 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
305 |
#endif |
306 |
|
307 |
if (allocated(frho)) deallocate(frho) |
308 |
allocate(frho(nlocal),stat=alloc_stat) |
309 |
if (alloc_stat /= 0) then |
310 |
status = -1 |
311 |
end if |
312 |
|
313 |
if (allocated(rho)) deallocate(rho) |
314 |
allocate(rho(nlocal),stat=alloc_stat) |
315 |
if (alloc_stat /= 0) then |
316 |
status = -1 |
317 |
end if |
318 |
|
319 |
if (allocated(dfrhodrho)) deallocate(dfrhodrho) |
320 |
allocate(dfrhodrho(nlocal),stat=alloc_stat) |
321 |
if (alloc_stat /= 0) then |
322 |
status = -1 |
323 |
end if |
324 |
|
325 |
#ifdef IS_MPI |
326 |
|
327 |
if (allocated(rho_tmp)) deallocate(rho_tmp) |
328 |
allocate(rho_tmp(nlocal),stat=alloc_stat) |
329 |
if (alloc_stat /= 0) then |
330 |
status = -1 |
331 |
end if |
332 |
|
333 |
|
334 |
if (allocated(frho_row)) deallocate(frho_row) |
335 |
allocate(frho_row(nAtomsInRow),stat=alloc_stat) |
336 |
if (alloc_stat /= 0) then |
337 |
status = -1 |
338 |
end if |
339 |
if (allocated(rho_row)) deallocate(rho_row) |
340 |
allocate(rho_row(nAtomsInRow),stat=alloc_stat) |
341 |
if (alloc_stat /= 0) then |
342 |
status = -1 |
343 |
end if |
344 |
if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row) |
345 |
allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat) |
346 |
if (alloc_stat /= 0) then |
347 |
status = -1 |
348 |
end if |
349 |
|
350 |
|
351 |
! Now do column arrays |
352 |
|
353 |
if (allocated(frho_col)) deallocate(frho_col) |
354 |
allocate(frho_col(nAtomsInCol),stat=alloc_stat) |
355 |
if (alloc_stat /= 0) then |
356 |
status = -1 |
357 |
end if |
358 |
if (allocated(rho_col)) deallocate(rho_col) |
359 |
allocate(rho_col(nAtomsInCol),stat=alloc_stat) |
360 |
if (alloc_stat /= 0) then |
361 |
status = -1 |
362 |
end if |
363 |
if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col) |
364 |
allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat) |
365 |
if (alloc_stat /= 0) then |
366 |
status = -1 |
367 |
end if |
368 |
|
369 |
#endif |
370 |
if (status == -1) then |
371 |
call handleError("SuttonChen:allocateSC","Error in allocating SC arrays") |
372 |
end if |
373 |
arraysAllocated = .true. |
374 |
end subroutine allocateSC |
375 |
|
376 |
subroutine setCutoffSC(rcut) |
377 |
real(kind=dp) :: rcut |
378 |
SC_rcut = rcut |
379 |
end subroutine setCutoffSC |
380 |
|
381 |
!! This array allocates module arrays if needed and builds mixing map. |
382 |
subroutine clean_SC() |
383 |
if (.not.arraysAllocated) call allocateSC() |
384 |
if (.not.haveMixingMap) call createMixingMap() |
385 |
! clean non-IS_MPI first |
386 |
frho = 0.0_dp |
387 |
rho = 0.0_dp |
388 |
dfrhodrho = 0.0_dp |
389 |
! clean MPI if needed |
390 |
#ifdef IS_MPI |
391 |
frho_row = 0.0_dp |
392 |
frho_col = 0.0_dp |
393 |
rho_row = 0.0_dp |
394 |
rho_col = 0.0_dp |
395 |
rho_tmp = 0.0_dp |
396 |
dfrhodrho_row = 0.0_dp |
397 |
dfrhodrho_col = 0.0_dp |
398 |
#endif |
399 |
end subroutine clean_SC |
400 |
|
401 |
!! Calculates rho_r |
402 |
subroutine calc_sc_prepair_rho(atom1, atom2, d, r, rijsq, rcut) |
403 |
integer :: atom1,atom2 |
404 |
real(kind = dp), dimension(3) :: d |
405 |
real(kind = dp), intent(inout) :: r |
406 |
real(kind = dp), intent(inout) :: rijsq, rcut |
407 |
! value of electron density rho do to atom i at atom j |
408 |
real(kind = dp) :: rho_i_at_j |
409 |
! value of electron density rho do to atom j at atom i |
410 |
real(kind = dp) :: rho_j_at_i |
411 |
|
412 |
! we don't use the derivatives, dummy variables |
413 |
real( kind = dp) :: drho |
414 |
integer :: sc_err |
415 |
|
416 |
integer :: atid1,atid2 ! Global atid |
417 |
integer :: myid_atom1 ! SC atid |
418 |
integer :: myid_atom2 |
419 |
|
420 |
|
421 |
! check to see if we need to be cleaned at the start of a force loop |
422 |
|
423 |
if (cleanArrays) call clean_SC() |
424 |
cleanArrays = .false. |
425 |
|
426 |
|
427 |
|
428 |
#ifdef IS_MPI |
429 |
Atid1 = Atid_row(Atom1) |
430 |
Atid2 = Atid_col(Atom2) |
431 |
#else |
432 |
Atid1 = Atid(Atom1) |
433 |
Atid2 = Atid(Atom2) |
434 |
#endif |
435 |
|
436 |
Myid_atom1 = SCList%atidtoSCtype(Atid1) |
437 |
Myid_atom2 = SCList%atidtoSCtype(Atid2) |
438 |
|
439 |
call lookupUniformSpline(MixingMap(myid_atom1,myid_atom2)%phi, r, & |
440 |
rho_i_at_j) |
441 |
rho_j_at_i = rho_i_at_j |
442 |
|
443 |
#ifdef IS_MPI |
444 |
rho_col(atom2) = rho_col(atom2) + rho_i_at_j |
445 |
rho_row(atom1) = rho_row(atom1) + rho_j_at_i |
446 |
#else |
447 |
rho(atom2) = rho(atom2) + rho_i_at_j |
448 |
rho(atom1) = rho(atom1) + rho_j_at_i |
449 |
#endif |
450 |
|
451 |
end subroutine calc_sc_prepair_rho |
452 |
|
453 |
|
454 |
!! Calculate the rho_a for all local atoms |
455 |
subroutine calc_sc_preforce_Frho(nlocal, pot, particle_pot) |
456 |
integer :: nlocal |
457 |
real(kind=dp) :: pot |
458 |
real(kind=dp), dimension(nlocal) :: particle_pot |
459 |
integer :: i,j |
460 |
integer :: atom |
461 |
integer :: atype1 |
462 |
integer :: atid1 |
463 |
integer :: myid |
464 |
|
465 |
!! Scatter the electron density from pre-pair calculation back to |
466 |
!! local atoms |
467 |
|
468 |
|
469 |
if (cleanArrays) call clean_SC() |
470 |
cleanArrays = .false. |
471 |
|
472 |
|
473 |
#ifdef IS_MPI |
474 |
call scatter(rho_row,rho,plan_atom_row,sc_err) |
475 |
if (sc_err /= 0 ) then |
476 |
write(errMsg,*) " Error scattering rho_row into rho" |
477 |
call handleError(RoutineName,errMesg) |
478 |
endif |
479 |
call scatter(rho_col,rho_tmp,plan_atom_col,sc_err) |
480 |
if (sc_err /= 0 ) then |
481 |
write(errMsg,*) " Error scattering rho_col into rho" |
482 |
call handleError(RoutineName,errMesg) |
483 |
|
484 |
endif |
485 |
|
486 |
rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) |
487 |
#endif |
488 |
|
489 |
!! Calculate F(rho) and derivative |
490 |
do atom = 1, nlocal |
491 |
Myid = SCList%atidtoSctype(Atid(atom)) |
492 |
! Myid is set to -1 for non SC atoms. |
493 |
! Punt if we are a non-SC atom type. |
494 |
if (Myid == -1) then |
495 |
frho(atom) = 0.0_dp |
496 |
dfrhodrho(atom) = 0.0_dp |
497 |
else |
498 |
frho(atom) = - SCList%SCTypes(Myid)%c * & |
499 |
SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom)) |
500 |
|
501 |
dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom) |
502 |
end if |
503 |
|
504 |
pot = pot + frho(atom) |
505 |
particle_pot(atom) = particle_pot(atom) + frho(atom) |
506 |
enddo |
507 |
|
508 |
#ifdef IS_MPI |
509 |
!! communicate f(rho) and derivatives back into row and column arrays |
510 |
call gather(frho,frho_row,plan_atom_row, sc_err) |
511 |
if (sc_err /= 0) then |
512 |
call handleError("calc_sc_forces()","MPI gather frho_row failure") |
513 |
endif |
514 |
call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err) |
515 |
if (sc_err /= 0) then |
516 |
call handleError("calc_sc_forces()","MPI gather dfrhodrho_row failure") |
517 |
endif |
518 |
call gather(frho,frho_col,plan_atom_col, sc_err) |
519 |
if (sc_err /= 0) then |
520 |
call handleError("calc_sc_forces()","MPI gather frho_col failure") |
521 |
endif |
522 |
call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err) |
523 |
if (sc_err /= 0) then |
524 |
call handleError("calc_sc_forces()","MPI gather dfrhodrho_col failure") |
525 |
endif |
526 |
#endif |
527 |
|
528 |
|
529 |
end subroutine calc_sc_preforce_Frho |
530 |
|
531 |
!! Does Sutton-Chen pairwise Force calculation. |
532 |
subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, & |
533 |
pot, f, do_pot) |
534 |
!Arguments |
535 |
integer, intent(in) :: atom1, atom2 |
536 |
real( kind = dp ), intent(in) :: rij, r2, rcut |
537 |
real( kind = dp ) :: pot, sw, vpair |
538 |
real( kind = dp ), dimension(3,nLocal) :: f |
539 |
real( kind = dp ), intent(in), dimension(3) :: d |
540 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
541 |
|
542 |
logical, intent(in) :: do_pot |
543 |
|
544 |
real( kind = dp ) :: drdx, drdy, drdz |
545 |
real( kind = dp ) :: dvpdr |
546 |
real( kind = dp ) :: rhtmp, drhodr |
547 |
real( kind = dp ) :: dudr |
548 |
real( kind = dp ) :: Fx,Fy,Fz |
549 |
real( kind = dp ) :: pot_temp, vptmp |
550 |
real( kind = dp ) :: rcij, vcij |
551 |
integer :: id1, id2 |
552 |
integer :: mytype_atom1 |
553 |
integer :: mytype_atom2 |
554 |
integer :: atid1, atid2 |
555 |
!Local Variables |
556 |
|
557 |
cleanArrays = .true. |
558 |
|
559 |
#ifdef IS_MPI |
560 |
atid1 = atid_row(atom1) |
561 |
atid2 = atid_col(atom2) |
562 |
#else |
563 |
atid1 = atid(atom1) |
564 |
atid2 = atid(atom2) |
565 |
#endif |
566 |
|
567 |
mytype_atom1 = SCList%atidToSCType(atid1) |
568 |
mytype_atom2 = SCList%atidTOSCType(atid2) |
569 |
|
570 |
drdx = d(1)/rij |
571 |
drdy = d(2)/rij |
572 |
drdz = d(3)/rij |
573 |
|
574 |
rcij = MixingMap(mytype_atom1,mytype_atom2)%rCut |
575 |
vcij = MixingMap(mytype_atom1,mytype_atom2)%vCut |
576 |
|
577 |
call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%phi, & |
578 |
rij, rhtmp, drhodr) |
579 |
|
580 |
call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%V, & |
581 |
rij, vptmp, dvpdr) |
582 |
|
583 |
#ifdef IS_MPI |
584 |
dudr = drhodr*(dfrhodrho_row(atom1) + dfrhodrho_col(atom2)) + dvpdr |
585 |
#else |
586 |
dudr = drhodr*(dfrhodrho(atom1)+ dfrhodrho(atom2)) + dvpdr |
587 |
#endif |
588 |
|
589 |
pot_temp = vptmp - vcij |
590 |
|
591 |
vpair = vpair + pot_temp |
592 |
|
593 |
fx = dudr * drdx |
594 |
fy = dudr * drdy |
595 |
fz = dudr * drdz |
596 |
|
597 |
#ifdef IS_MPI |
598 |
if (do_pot) then |
599 |
pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5 |
600 |
pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5 |
601 |
end if |
602 |
|
603 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
604 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
605 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
606 |
|
607 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
608 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
609 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
610 |
#else |
611 |
|
612 |
if(do_pot) then |
613 |
pot = pot + pot_temp |
614 |
end if |
615 |
|
616 |
f(1,atom1) = f(1,atom1) + fx |
617 |
f(2,atom1) = f(2,atom1) + fy |
618 |
f(3,atom1) = f(3,atom1) + fz |
619 |
|
620 |
f(1,atom2) = f(1,atom2) - fx |
621 |
f(2,atom2) = f(2,atom2) - fy |
622 |
f(3,atom2) = f(3,atom2) - fz |
623 |
#endif |
624 |
|
625 |
|
626 |
#ifdef IS_MPI |
627 |
id1 = AtomRowToGlobal(atom1) |
628 |
id2 = AtomColToGlobal(atom2) |
629 |
#else |
630 |
id1 = atom1 |
631 |
id2 = atom2 |
632 |
#endif |
633 |
|
634 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
635 |
|
636 |
fpair(1) = fpair(1) + fx |
637 |
fpair(2) = fpair(2) + fy |
638 |
fpair(3) = fpair(3) + fz |
639 |
|
640 |
endif |
641 |
end subroutine do_sc_pair |
642 |
end module suttonchen |