1 |
!! |
2 |
!! Copyright (c) 2005, 2009 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. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
|
42 |
!! Implements 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 |
implicit none |
56 |
PRIVATE |
57 |
#define __FORTRAN90 |
58 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
59 |
|
60 |
!! number of points for the spline approximations |
61 |
INTEGER, PARAMETER :: np = 3000 |
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 |
logical, save :: haveMixingMap = .false. |
68 |
logical, save :: useGeometricDistanceMixing = .false. |
69 |
logical, save :: cleanArrays = .true. |
70 |
logical, save :: arraysAllocated = .false. |
71 |
|
72 |
|
73 |
character(len = statusMsgSize) :: errMesg |
74 |
integer :: sc_err |
75 |
|
76 |
character(len = 200) :: errMsg |
77 |
character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE" |
78 |
|
79 |
type, private :: SCtype |
80 |
integer :: atid |
81 |
real(kind=dp) :: c |
82 |
real(kind=dp) :: m |
83 |
real(kind=dp) :: n |
84 |
real(kind=dp) :: alpha |
85 |
real(kind=dp) :: epsilon |
86 |
real(kind=dp) :: sc_rcut |
87 |
end type SCtype |
88 |
|
89 |
|
90 |
type, private :: SCTypeList |
91 |
integer :: nSCTypes = 0 |
92 |
integer :: currentSCtype = 0 |
93 |
type (SCtype), pointer :: SCtypes(:) => null() |
94 |
integer, pointer :: atidToSCtype(:) => null() |
95 |
end type SCTypeList |
96 |
|
97 |
type (SCTypeList), save :: SCList |
98 |
|
99 |
type:: MixParameters |
100 |
real(kind=DP) :: alpha |
101 |
real(kind=DP) :: epsilon |
102 |
real(kind=DP) :: m |
103 |
real(Kind=DP) :: n |
104 |
real(kind=dp) :: rCut |
105 |
real(kind=dp) :: vCut |
106 |
logical :: rCutWasSet = .false. |
107 |
type(cubicSpline) :: V |
108 |
type(cubicSpline) :: phi |
109 |
end type MixParameters |
110 |
|
111 |
type(MixParameters), dimension(:,:), allocatable :: MixingMap |
112 |
|
113 |
public :: setCutoffSC |
114 |
public :: do_SC_pair |
115 |
public :: newSCtype |
116 |
public :: calc_SC_prepair_rho |
117 |
public :: calc_SC_preforce_Frho |
118 |
public :: destroySCtypes |
119 |
public :: getSCCut |
120 |
! public :: setSCDefaultCutoff |
121 |
! public :: setSCUniformCutoff |
122 |
|
123 |
|
124 |
contains |
125 |
|
126 |
|
127 |
subroutine newSCtype(c_ident,c,m,n,alpha,epsilon,status) |
128 |
real (kind = dp ) :: c ! Density Scaling |
129 |
real (kind = dp ) :: m ! Density Exponent |
130 |
real (kind = dp ) :: n ! Pair Potential Exponent |
131 |
real (kind = dp ) :: alpha ! Length Scaling |
132 |
real (kind = dp ) :: epsilon ! Energy Scaling |
133 |
integer :: c_ident |
134 |
integer :: status |
135 |
integer :: nAtypes,nSCTypes,myATID |
136 |
integer :: maxVals |
137 |
integer :: alloc_stat |
138 |
integer :: current |
139 |
integer,pointer :: Matchlist(:) => null() |
140 |
|
141 |
status = 0 |
142 |
|
143 |
|
144 |
!! Assume that atypes has already been set and get the total number of types in atypes |
145 |
|
146 |
! check to see if this is the first time into |
147 |
if (.not.associated(SCList%SCTypes)) then |
148 |
call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList) |
149 |
SCList%nSCtypes = nSCtypes |
150 |
allocate(SCList%SCTypes(nSCTypes)) |
151 |
nAtypes = getSize(atypes) |
152 |
allocate(SCList%atidToSCType(nAtypes)) |
153 |
SCList%atidToSCType = -1 |
154 |
end if |
155 |
|
156 |
SCList%currentSCType = SCList%currentSCType + 1 |
157 |
current = SCList%currentSCType |
158 |
|
159 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
160 |
SCList%atidToSCType(myATID) = current |
161 |
|
162 |
|
163 |
SCList%SCTypes(current)%atid = c_ident |
164 |
SCList%SCTypes(current)%alpha = alpha |
165 |
SCList%SCTypes(current)%c = c |
166 |
SCList%SCTypes(current)%m = m |
167 |
SCList%SCTypes(current)%n = n |
168 |
SCList%SCTypes(current)%epsilon = epsilon |
169 |
end subroutine newSCtype |
170 |
|
171 |
|
172 |
subroutine destroySCTypes() |
173 |
if (associated(SCList%SCtypes)) then |
174 |
deallocate(SCList%SCTypes) |
175 |
SCList%SCTypes=>null() |
176 |
end if |
177 |
if (associated(SCList%atidToSCtype)) then |
178 |
deallocate(SCList%atidToSCtype) |
179 |
SCList%atidToSCtype=>null() |
180 |
end if |
181 |
! Reset Capacity |
182 |
SCList%nSCTypes = 0 |
183 |
SCList%currentSCtype=0 |
184 |
|
185 |
end subroutine destroySCTypes |
186 |
|
187 |
function getSCCut(atomID) result(cutValue) |
188 |
integer, intent(in) :: atomID |
189 |
integer :: scID |
190 |
real(kind=dp) :: cutValue |
191 |
|
192 |
scID = SCList%atidToSCType(atomID) |
193 |
cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha |
194 |
end function getSCCut |
195 |
|
196 |
subroutine createMixingMap() |
197 |
integer :: nSCtypes, i, j, k |
198 |
real ( kind = dp ) :: e1, e2, m1, m2, alpha1, alpha2, n1, n2 |
199 |
real ( kind = dp ) :: epsilon, m, n, alpha, rCut, vCut, dr, r |
200 |
real ( kind = dp ), dimension(np) :: rvals, vvals, phivals |
201 |
|
202 |
if (SCList%currentSCtype == 0) then |
203 |
call handleError("SuttonChen", "No members in SCMap") |
204 |
return |
205 |
end if |
206 |
|
207 |
nSCtypes = SCList%nSCtypes |
208 |
|
209 |
if (.not. allocated(MixingMap)) then |
210 |
allocate(MixingMap(nSCtypes, nSCtypes)) |
211 |
endif |
212 |
useGeometricDistanceMixing = usesGeometricDistanceMixing() |
213 |
do i = 1, nSCtypes |
214 |
|
215 |
e1 = SCList%SCtypes(i)%epsilon |
216 |
m1 = SCList%SCtypes(i)%m |
217 |
n1 = SCList%SCtypes(i)%n |
218 |
alpha1 = SCList%SCtypes(i)%alpha |
219 |
|
220 |
do j = 1, nSCtypes |
221 |
|
222 |
e2 = SCList%SCtypes(j)%epsilon |
223 |
m2 = SCList%SCtypes(j)%m |
224 |
n2 = SCList%SCtypes(j)%n |
225 |
alpha2 = SCList%SCtypes(j)%alpha |
226 |
|
227 |
if (useGeometricDistanceMixing) then |
228 |
alpha = sqrt(alpha1 * alpha2) ! SC formulation |
229 |
else |
230 |
alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation |
231 |
endif |
232 |
rcut = 2.0_dp * alpha |
233 |
epsilon = sqrt(e1 * e2) |
234 |
m = 0.5_dp*(m1+m2) |
235 |
n = 0.5_dp*(n1+n2) |
236 |
|
237 |
dr = (rCut) / dble(np-1) |
238 |
rvals(1) = 0.0_dp |
239 |
vvals(1) = 0.0_dp |
240 |
phivals(1) = 0.0_dp |
241 |
|
242 |
do k = 2, np |
243 |
r = dble(k-1)*dr |
244 |
rvals(k) = r |
245 |
vvals(k) = epsilon*((alpha/r)**n) |
246 |
phivals(k) = (alpha/r)**m |
247 |
enddo |
248 |
|
249 |
vCut = epsilon*((alpha/rCut)**n) |
250 |
|
251 |
call newSpline(MixingMap(i,j)%V, rvals, vvals, .true.) |
252 |
call newSpline(MixingMap(i,j)%phi, rvals, phivals, .true.) |
253 |
|
254 |
MixingMap(i,j)%epsilon = epsilon |
255 |
MixingMap(i,j)%m = m |
256 |
MixingMap(i,j)%n = n |
257 |
MixingMap(i,j)%alpha = alpha |
258 |
MixingMap(i,j)%rCut = rcut |
259 |
MixingMap(i,j)%vCut = vCut |
260 |
enddo |
261 |
enddo |
262 |
|
263 |
haveMixingMap = .true. |
264 |
|
265 |
end subroutine createMixingMap |
266 |
|
267 |
|
268 |
|
269 |
subroutine setCutoffSC(rcut) |
270 |
real(kind=dp) :: rcut |
271 |
SC_rcut = rcut |
272 |
end subroutine setCutoffSC |
273 |
|
274 |
|
275 |
!! Calculates rho_r |
276 |
subroutine calc_sc_prepair_rho(atid1, atid2, d, r, rijsq, rho_i_at_j, rho_j_at_i) |
277 |
integer :: atid1, atid2 |
278 |
real(kind = dp), dimension(3) :: d |
279 |
real(kind = dp), intent(inout) :: r |
280 |
real(kind = dp), intent(inout) :: rijsq |
281 |
! value of electron density rho do to atom i at atom j |
282 |
real(kind = dp), intent(inout) :: rho_i_at_j |
283 |
! value of electron density rho do to atom j at atom i |
284 |
real(kind = dp), intent(inout) :: rho_j_at_i |
285 |
! we don't use the derivatives, dummy variables |
286 |
real( kind = dp) :: drho |
287 |
integer :: sc_err |
288 |
|
289 |
integer :: myid_atom1 ! SC atid |
290 |
integer :: myid_atom2 |
291 |
|
292 |
! check to see if we need to be cleaned at the start of a force loop |
293 |
|
294 |
if (.not.haveMixingMap) call createMixingMap() |
295 |
haveMixingMap = .true. |
296 |
|
297 |
Myid_atom1 = SCList%atidtoSCtype(Atid1) |
298 |
Myid_atom2 = SCList%atidtoSCtype(Atid2) |
299 |
|
300 |
call lookupUniformSpline(MixingMap(myid_atom1,myid_atom2)%phi, r, & |
301 |
rho_i_at_j) |
302 |
rho_j_at_i = rho_i_at_j |
303 |
|
304 |
end subroutine calc_sc_prepair_rho |
305 |
|
306 |
|
307 |
!! Calculate the rho_a for all local atoms |
308 |
subroutine calc_sc_preforce_Frho(nlocal, pot, particle_pot) |
309 |
integer :: nlocal |
310 |
real(kind=dp) :: pot |
311 |
real(kind=dp), dimension(nlocal) :: particle_pot |
312 |
integer :: i,j |
313 |
integer :: atom |
314 |
integer :: atype1 |
315 |
integer :: atid1 |
316 |
integer :: myid |
317 |
|
318 |
!! Calculate F(rho) and derivative |
319 |
do atom = 1, nlocal |
320 |
Myid = SCList%atidtoSctype(Atid(atom)) |
321 |
! Myid is set to -1 for non SC atoms. |
322 |
! Punt if we are a non-SC atom type. |
323 |
if (Myid == -1) then |
324 |
frho(atom) = 0.0_dp |
325 |
dfrhodrho(atom) = 0.0_dp |
326 |
else |
327 |
frho(atom) = - SCList%SCTypes(Myid)%c * & |
328 |
SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom)) |
329 |
|
330 |
dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom) |
331 |
end if |
332 |
pot = pot + frho(atom) |
333 |
particle_pot(atom) = particle_pot(atom) + frho(atom) |
334 |
enddo |
335 |
|
336 |
end subroutine calc_sc_preforce_Frho |
337 |
|
338 |
!! Does Sutton-Chen pairwise Force calculation. |
339 |
subroutine do_sc_pair(atid1, atid2, d, rij, r2, sw, vpair, & |
340 |
pot, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, & |
341 |
fshift_i, fshift_j) |
342 |
!Arguments |
343 |
integer, intent(in) :: atid1, atid2 |
344 |
real( kind = dp ), intent(in) :: rij, r2 |
345 |
real( kind = dp ) :: pot, sw, vpair |
346 |
real( kind = dp ), dimension(3) :: f1 |
347 |
real( kind = dp ), intent(in), dimension(3) :: d |
348 |
real( kind = dp ), intent(inout) :: dfrhodrho_i, dfrhodrho_j |
349 |
real( kind = dp ), intent(inout) :: rho_i, rho_j |
350 |
real( kind = dp ), intent(inout):: fshift_i, fshift_j |
351 |
|
352 |
real( kind = dp ) :: drdx, drdy, drdz |
353 |
real( kind = dp ) :: dvpdr |
354 |
real( kind = dp ) :: rhtmp, drhodr |
355 |
real( kind = dp ) :: dudr |
356 |
real( kind = dp ) :: Fx,Fy,Fz |
357 |
real( kind = dp ) :: pot_temp, vptmp |
358 |
real( kind = dp ) :: rcij, vcij |
359 |
|
360 |
integer :: id1, id2 |
361 |
integer :: mytype_atom1 |
362 |
integer :: mytype_atom2 |
363 |
|
364 |
mytype_atom1 = SCList%atidToSCType(atid1) |
365 |
mytype_atom2 = SCList%atidTOSCType(atid2) |
366 |
|
367 |
drdx = d(1)/rij |
368 |
drdy = d(2)/rij |
369 |
drdz = d(3)/rij |
370 |
|
371 |
rcij = MixingMap(mytype_atom1,mytype_atom2)%rCut |
372 |
vcij = MixingMap(mytype_atom1,mytype_atom2)%vCut |
373 |
|
374 |
call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%phi, & |
375 |
rij, rhtmp, drhodr) |
376 |
|
377 |
call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%V, & |
378 |
rij, vptmp, dvpdr) |
379 |
|
380 |
dudr = drhodr*(dfrhodrho_i + dfrhodrho_j) + dvpdr |
381 |
|
382 |
pot_temp = vptmp - vcij |
383 |
|
384 |
vpair = vpair + pot_temp |
385 |
|
386 |
fx = dudr * drdx |
387 |
fy = dudr * drdy |
388 |
fz = dudr * drdz |
389 |
|
390 |
! particle_pot is the difference between the full potential |
391 |
! and the full potential without the presence of a particular |
392 |
! particle (atom1). |
393 |
! |
394 |
! This reduces the density at other particle locations, so |
395 |
! we need to recompute the density at atom2 assuming atom1 |
396 |
! didn't contribute. This then requires recomputing the |
397 |
! density functional for atom2 as well. |
398 |
! |
399 |
! Most of the particle_pot heavy lifting comes from the |
400 |
! pair interaction, and will be handled by vpair. |
401 |
|
402 |
fshift_i = -SCList%SCTypes(mytype_atom1)%c * & |
403 |
SCList%SCTypes(mytype_atom1)%epsilon * sqrt(rho_i-rhtmp) |
404 |
fshift_j = -SCList%SCTypes(mytype_atom2)%c * & |
405 |
SCList%SCTypes(mytype_atom2)%epsilon * sqrt(rho_j-rhtmp) |
406 |
|
407 |
pot = pot + pot_temp |
408 |
|
409 |
f1(1) = f1(1) + fx |
410 |
f1(2) = f1(2) + fy |
411 |
f1(3) = f1(3) + fz |
412 |
|
413 |
end subroutine do_sc_pair |
414 |
end module suttonchen |