-
Notifications
You must be signed in to change notification settings - Fork 1
/
pair_potential_util.F
325 lines (287 loc) · 16.1 KB
/
pair_potential_util.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \par History
!> September 2005 - Introduced the Born-Mayer-Huggins-Fumi-Tosi Potential (BMHTF)
!> 2006 - Major rewriting of the routines.. Linear scaling setup of splines
!> \author CJM
! **************************************************************************************************
MODULE pair_potential_util
USE fparser, ONLY: EvalErrType,&
evalf
USE kinds, ONLY: dp
USE mathconstants, ONLY: ifac
USE pair_potential_types, ONLY: &
b4_type, bm_type, ea_type, ft_type, ftd_type, gp_type, gw_type, ip_type, lj_charmm_type, &
lj_type, not_initialized, pair_potential_single_type, tab_type, wl_type
USE physcon, ONLY: bohr,&
evolt
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'pair_potential_util'
REAL(KIND=dp), PARAMETER, PRIVATE :: MIN_HICUT_VALUE = 1.0E-15_dp, &
DEFAULT_HICUT_VALUE = 1.0E3_dp
PUBLIC :: ener_pot, ener_zbl, zbl_matching_polinomial
CONTAINS
! **************************************************************************************************
!> \brief Evaluates the nonbond potential energy for the implemented FF kinds
!> \param pot ...
!> \param r ...
!> \param energy_cutoff ...
!> \return ...
! **************************************************************************************************
FUNCTION ener_pot(pot, r, energy_cutoff) RESULT(value)
TYPE(pair_potential_single_type), POINTER :: pot
REAL(KIND=dp), INTENT(IN) :: r, energy_cutoff
REAL(KIND=dp) :: value
INTEGER :: i, index, index1, index2, j, n
REAL(KIND=dp) :: bd6, bd8, dampsum, f6, f8, lvalue, pp, &
qq, scale, xf
value = 0.0_dp
DO j = 1, SIZE(pot%type)
! A lower boundary for the potential definition was defined
IF ((pot%set(j)%rmin /= not_initialized) .AND. (r < pot%set(j)%rmin)) CYCLE
! An upper boundary for the potential definition was defined
IF ((pot%set(j)%rmax /= not_initialized) .AND. (r >= pot%set(j)%rmax)) CYCLE
! If within limits let's compute the potential...
IF (pot%type(j) == lj_charmm_type) THEN
lvalue = &
4.0_dp*pot%set(j)%lj%epsilon*(pot%set(j)%lj%sigma12*r**(-12) - pot%set(j)%lj% &
sigma6*r**(-6))
ELSE IF (pot%type(j) == lj_type) THEN
lvalue = pot%set(j)%lj%epsilon* &
(pot%set(j)%lj%sigma12*r**(-12) - pot%set(j)%lj%sigma6*r**(-6))
ELSE IF (pot%type(j) == ip_type) THEN
lvalue = 0._dp
IF (r > pot%set(j)%ipbv%rcore) THEN
DO i = 2, 15
lvalue = lvalue + pot%set(j)%ipbv%a(i)/(r**(i - 1)*REAL(i - 1, dp))
END DO
ELSE
! use a linear potential
lvalue = pot%set(j)%ipbv%m*r + pot%set(j)%ipbv%b
END IF
lvalue = lvalue
ELSE IF (pot%type(j) == wl_type) THEN
lvalue = pot%set(j)%willis%a*EXP(-pot%set(j)%willis%b*r) - pot%set(j)%willis%c/r**6
ELSE IF (pot%type(j) == gw_type) THEN
scale = EXP(pot%set(j)%goodwin%m*(-(r/pot%set(j)%goodwin%dc)**pot%set(j)%goodwin%mc + &
(pot%set(j)%goodwin%d/pot%set(j)%goodwin%dc)**pot%set(j)%goodwin%mc))
lvalue = scale*pot%set(j)%goodwin%vr0*(pot%set(j)%goodwin%d/r)**pot%set(j)%goodwin%m
ELSE IF (pot%type(j) == ft_type) THEN
lvalue = pot%set(j)%ft%a*EXP(-pot%set(j)%ft%b*r) - pot%set(j)%ft%c/r**6 - pot%set(j)%ft%d/r**8
ELSE IF (pot%type(j) == ftd_type) THEN
! Compute 6th order dispersion correction term
bd6 = pot%set(j)%ftd%bd(1)
dampsum = 1.0_dp
xf = 1.0_dp
DO i = 1, 6
xf = xf*bd6*r
dampsum = dampsum + xf*ifac(i)
END DO
f6 = 1.0_dp - EXP(-bd6*r)*dampsum
! Compute 8th order dispersion correction term
bd8 = pot%set(j)%ftd%bd(2)
dampsum = 1.0_dp
xf = 1.0_dp
DO i = 1, 8
xf = xf*bd8*r
dampsum = dampsum + xf*ifac(i)
END DO
f8 = 1.0_dp - EXP(-bd8*r)*dampsum
lvalue = pot%set(j)%ftd%a*EXP(-pot%set(j)%ftd%b*r) - f6*pot%set(j)%ftd%c/r**6 - f8*pot%set(j)%ftd%d/r**8
ELSE IF (pot%type(j) == ea_type) THEN
index = INT(r/pot%set(j)%eam%drar) + 1
IF (index > pot%set(j)%eam%npoints) THEN
index = pot%set(j)%eam%npoints
ELSEIF (index < 1) THEN
index = 1
END IF
qq = r - pot%set(j)%eam%rval(index)
pp = pot%set(j)%eam%phi(index) + &
qq*pot%set(j)%eam%phip(index)
lvalue = pp
ELSE IF (pot%type(j) == b4_type) THEN
IF (r <= pot%set(j)%buck4r%r1) THEN
pp = pot%set(j)%buck4r%a*EXP(-pot%set(j)%buck4r%b*r)
ELSEIF (r > pot%set(j)%buck4r%r1 .AND. r <= pot%set(j)%buck4r%r2) THEN
pp = 0.0_dp
DO n = 0, pot%set(j)%buck4r%npoly1
pp = pp + pot%set(j)%buck4r%poly1(n)*r**n
END DO
ELSEIF (r > pot%set(j)%buck4r%r2 .AND. r <= pot%set(j)%buck4r%r3) THEN
pp = 0.0_dp
DO n = 0, pot%set(j)%buck4r%npoly2
pp = pp + pot%set(j)%buck4r%poly2(n)*r**n
END DO
ELSEIF (r > pot%set(j)%buck4r%r3) THEN
pp = -pot%set(j)%buck4r%c/r**6
END IF
lvalue = pp
ELSE IF (pot%type(j) == tab_type) THEN
index1 = FLOOR((r - pot%set(j)%tab%r(1))/pot%set(j)%tab%dr) + 1
index2 = index1 + 1
IF (index2 > pot%set(j)%tab%npoints) THEN
index2 = pot%set(j)%tab%npoints
index1 = index2 - 1
ELSEIF (index1 < 1) THEN
index1 = 1
index2 = 2
END IF
pp = pot%set(j)%tab%e(index1) + (r - pot%set(j)%tab%r(index1))* &
(pot%set(j)%tab%e(index2) - pot%set(j)%tab%e(index1))/ &
(pot%set(j)%tab%r(index2) - pot%set(j)%tab%r(index1))
lvalue = pp
ELSE IF (pot%type(j) == bm_type) THEN
lvalue = pot%set(j)%buckmo%f0*(pot%set(j)%buckmo%b1 + pot%set(j)%buckmo%b2)* &
EXP((pot%set(j)%buckmo%a1 + pot%set(j)%buckmo%a2 - r)/(pot%set(j)%buckmo%b1 + pot%set(j)%buckmo%b2)) &
- pot%set(j)%buckmo%c/r**6 &
+ pot%set(j)%buckmo%d*(EXP(-2._dp*pot%set(j)%buckmo%beta*(r - pot%set(j)%buckmo%r0)) - &
2.0_dp*EXP(-pot%set(j)%buckmo%beta*(r - pot%set(j)%buckmo%r0)))
ELSE IF (pot%type(j) == gp_type) THEN
pot%set(j)%gp%values(1) = r
lvalue = evalf(pot%set(j)%gp%myid, pot%set(j)%gp%values)
IF (EvalErrType > 0) &
CPABORT("Error evaluating generic potential energy function")
ELSE
lvalue = 0.0_dp
END IF
value = value + lvalue
END DO
value = value - energy_cutoff
END FUNCTION ener_pot
! **************************************************************************************************
!> \brief Evaluates the ZBL scattering potential, very short range
!> Only shell-model for interactions among pairs without repulsive term
!> \param pot ...
!> \param r ...
!> \return ...
!> \author i soliti ignoti
! **************************************************************************************************
FUNCTION ener_zbl(pot, r)
TYPE(pair_potential_single_type), POINTER :: pot
REAL(KIND=dp), INTENT(IN) :: r
REAL(KIND=dp) :: ener_zbl
REAL(KIND=dp) :: au, fac, x
ener_zbl = 0.0_dp
IF (r <= pot%zbl_rcut(1)) THEN
au = 0.88534_dp*bohr/(pot%z1**0.23_dp + pot%z2**0.23_dp)
x = r/au
fac = pot%z1*pot%z2/evolt
ener_zbl = fac/r*(0.1818_dp*EXP(-3.2_dp*x) + 0.5099_dp*EXP(-0.9423_dp*x) + &
0.2802_dp*EXP(-0.4029_dp*x) + 0.02817_dp*EXP(-0.2016_dp*x))
ELSEIF (r > pot%zbl_rcut(1) .AND. r <= pot%zbl_rcut(2)) THEN
ener_zbl = pot%zbl_poly(0) + pot%zbl_poly(1)*r + pot%zbl_poly(2)*r*r + pot%zbl_poly(3)*r*r*r + &
pot%zbl_poly(4)*r*r*r*r + pot%zbl_poly(5)*r*r*r*r*r
ELSE
ener_zbl = 0.0_dp
END IF
END FUNCTION ener_zbl
! **************************************************************************************************
!> \brief Determine the polinomial coefficients used to set to zero the zbl potential
!> at the cutoff radius, with continuity in function, first and second derivative
!> Only shell-model for interactions among pairs without repulsive term
!> \param pot ...
!> \param rcov1 ...
!> \param rcov2 ...
!> \param z1 ...
!> \param z2 ...
!> \author i soliti ignoti
! **************************************************************************************************
SUBROUTINE zbl_matching_polinomial(pot, rcov1, rcov2, z1, z2)
TYPE(pair_potential_single_type), POINTER :: pot
REAL(KIND=dp), INTENT(IN) :: rcov1, rcov2, z1, z2
REAL(KIND=dp) :: au, d1, d2, dd1, dd2, fac, v1, v2, x, &
x1, x2
pot%zbl_rcut(1) = (rcov1 + rcov2)*(1.0_dp - 0.2_dp)*bohr
pot%zbl_rcut(2) = (rcov1 + rcov2)*bohr
x1 = pot%zbl_rcut(1)
x2 = pot%zbl_rcut(2)
pot%z1 = z1
pot%z2 = z2
au = 0.88534_dp*bohr/(z1**0.23_dp + z2**0.23_dp)
x = x1/au
fac = z1*z2/evolt
v1 = fac/x1*(0.1818_dp*EXP(-3.2_dp*x) + 0.5099_dp*EXP(-0.9423_dp*x) + &
0.2802_dp*EXP(-0.4029_dp*x) + 0.02817_dp*EXP(-0.2016_dp*x))
d1 = fac/x1/au*(-3.2_dp*0.1818_dp*EXP(-3.2_dp*x) - 0.9423_dp*0.5099_dp*EXP(-0.9423_dp*x) &
- 0.4029_dp*0.2802_dp*EXP(-0.4029_dp*x) - 0.2016_dp*0.02817_dp*EXP(-0.2016_dp*x)) &
- fac/x1/x1*(0.1818_dp*EXP(-3.2_dp*x) + 0.5099_dp*EXP(-0.9423_dp*x) + &
0.2802_dp*EXP(-0.4029_dp*x) + 0.02817_dp*EXP(-0.2016_dp*x))
dd1 = 2.0_dp*fac/x1**3*(0.1818_dp*EXP(-0.32E1_dp*x) &
+ 0.5099_dp*EXP(-0.9423_dp*x) + 0.2802_dp*EXP(-0.4029_dp*x) &
+ 0.2817E-1_dp*EXP(-0.2016_dp*x)) &
- 0.2E1_dp*fac/x1**2/au*(-0.58176_dp*EXP(-0.32E1_dp*x) - 0.48047877_dp*EXP(-0.9423_dp*x) &
- 0.11289258_dp*EXP(-0.4029_dp*x) - 0.5679072E-2_dp*EXP(-0.2016_dp*x)) &
+ fac/x1/au**2*(0.1861632E1_dp*EXP(-0.32E1_dp*x) + &
0.4527551450_dp*EXP(-0.9423_dp*x) + 0.4548442048E-1_dp*EXP(-0.4029_dp*x) + &
0.1144900915E-2_dp*EXP(-0.2016_dp*x))
v2 = 0.0_dp
d2 = 0.0_dp
dd2 = 0.0_dp
CALL compute_polinomial_5th(x1, v1, d1, dd1, x2, v2, d2, dd2, pot%zbl_poly)
END SUBROUTINE zbl_matching_polinomial
! **************************************************************************************************
!> \brief ...
!> \param r1 ...
!> \param v1 ...
!> \param d1 ...
!> \param dd1 ...
!> \param r2 ...
!> \param v2 ...
!> \param d2 ...
!> \param dd2 ...
!> \param poly ...
! **************************************************************************************************
SUBROUTINE compute_polinomial_5th(r1, v1, d1, dd1, r2, v2, d2, dd2, poly)
REAL(KIND=dp) :: r1, v1, d1, dd1, r2, v2, d2, dd2, &
poly(0:5)
REAL(KIND=dp) :: a0, a1, a2, a3, a4, a5
! 5th order
a0 = .5_dp*(2._dp*r1**5*v2 - 2._dp*v1*r2**5 + 10._dp*v1*r2**4*r1 - 20._dp*v1*r1**2*r2**3 - r1**2*dd1*r2**5 - &
r1**4*r2**3*dd1 + 20._dp*r1**3*r2**2*v2 + 2._dp*r1**3*r2**4*dd1 + r1**3*r2**4*dd2 - 8._dp*r1**3*r2**3*d2 - &
2._dp*r1**4*r2**3*dd2 + 10._dp*r1**4*r2**2*d2 - 10._dp*r1**4*r2*v2 + 2._dp*r1*d1*r2**5 - 10._dp*r1**2*d1*r2**4 + &
8._dp*r1**3*d1*r2**3 - 2._dp*r1**5*r2*d2 + r1**5*r2**2*dd2)/ &
(10.*r2**2*r1**3 - 5._dp*r2*r1**4 - 10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
a1 = -.5_dp*(-4._dp*r2**3*r1**3*dd2 + 24._dp*r2**2*r1**3*d1 + 4._dp*r2**3*r1**3*dd1 + 3._dp*r2**4*r1**2*dd2 + &
r2**4*r1**2*dd1 - 2._dp*r2**5*r1*dd1 - 10._dp*r2**4*r1*d1 + 10._dp*r1**4*r2*d2 - &
r1**4*r2**2*dd2 - 3._dp*r1**4*r2**2*dd1 + &
2._dp*r1**5*r2*dd2 - 24._dp*r2**3*r1**2*d2 - 16._dp*r2**3*r1**2*d1 + &
16._dp*r2**2*r1**3*d2 - 2._dp*r1**5*d2 + 2._dp*r2**5*d1 - &
60._dp*r1**2*r2**2*v1 + 60._dp*r1**2*r2**2*v2)/ &
(10._dp*r2**2*r1**3 - 5._dp*r2*r1**4 - 10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
a2 = .5_dp*(60._dp*r1**2*r2*v2 - 60._dp*v1*r1*r2**2 - 12._dp*r1**2*r2**2*d2 - 36._dp*r1*d1*r2**3 + 3._dp*r2**4*r1*dd2 - &
24._dp*r2**3*r1*d2 - 4._dp*r2**4*r1*dd1 + 12._dp*r1**2*r2**2*d1 - 8._dp*r1**3*r2**2*dd2 + 24._dp*r1**3*r2*d1 + &
4._dp*r1**4*r2*dd2 + 36._dp*r1**3*r2*d2 - 3._dp*r1**4*r2*dd1 + 8._dp*r2**3*r1**2*dd1 + 60._dp*r2**2*r1*v2 - &
60._dp*r1**2*v1*r2 + r1**5*dd2 - r2**5*dd1)/ &
(10._dp*r2**2*r1**3 - 5._dp*r2*r1**4 - 10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
a3 = -.5_dp*(3._dp*r1**4*dd2 - r1**4*dd1 + 8.*r1**3*d1 - 4.*r1**3*r2*dd1 + &
12._dp*r1**3*d2 + 32._dp*r1**2*r2*d1 - 8._dp*r1**2*r2**2*dd2 - &
20._dp*r1**2*v1 + 8._dp*r1**2*r2**2*dd1 + 28._dp*r1**2*r2*d2 + &
20._dp*r1**2*v2 + 80._dp*r1*r2*v2 - 28._dp*r2**2*r1*d1 - 80._dp*r1*v1*r2 - &
32._dp*r2**2*r1*d2 + 4._dp*r1*r2**3*dd2 - 8._dp*r2**3*d2 - 12._dp*r2**3*d1 + &
r2**4*dd2 - 3._dp*r2**4*dd1 + 20._dp*r2**2*v2 - 20._dp*r2**2*v1)/ &
(10._dp*r2**2*r1**3 - 5._dp*r2*r1**4 - 10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
a4 = .5_dp*(3._dp*r1**3*dd2 - 2._dp*r1**3*dd1 + r1**2*r2*dd1 + 14.*r1**2*d1 - 4._dp*r1**2*r2*dd2 + &
16._dp*r1**2*d2 - 2._dp*r1*r2*d2 - r1*r2**2*dd2 - &
30._dp*r1*v1 + 30.*r1*v2 + 2._dp*r1*r2*d1 + 4._dp*r1*r2**2*dd1 - 16._dp*r2**2*d1 + &
2._dp*r2**3*dd2 - 14._dp*r2**2*d2 + 30._dp*r2*v2 - 30._dp*v1*r2 - &
3._dp*r2**3*dd1)/(10._dp*r2**2*r1**3 - 5._dp*r2*r1**4 - &
10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
a5 = -.5_dp*(6._dp*r1*d1 + 2._dp*r2*r1*dd1 + 6._dp*r1*d2 - 2.*r2*r1*dd2 - &
r2**2*dd1 - r1**2*dd1 - 12.*v1 + 12._dp*v2 + r1**2*dd2 - &
6._dp*r2*d1 + r2**2*dd2 - 6._dp*r2*d2)/ &
(10._dp*r2**2*r1**3 - 5._dp*r2*r1**4 - 10._dp*r1**2*r2**3 + 5._dp*r2**4*r1 - r2**5 + r1**5)
poly(0) = a0
poly(1) = a1
poly(2) = a2
poly(3) = a3
poly(4) = a4
poly(5) = a5
END SUBROUTINE compute_polinomial_5th
END MODULE pair_potential_util