-
Notifications
You must be signed in to change notification settings - Fork 4
/
real_space_electrostatic_sum.f90
361 lines (287 loc) · 12.2 KB
/
real_space_electrostatic_sum.f90
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
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
! MIT License
!
! Copyright (c) 2019-2020 William C. Witt
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
module real_space_electrostatic_sum
implicit none
integer, parameter :: dp = 8
real(dp), parameter :: pi = 3.14159265358979323846264338327950288419_dp
real(dp), parameter :: sqrt_pi = sqrt(pi)
real(dp), parameter :: one_third = 1.0_dp/3.0_dp
contains
subroutine energy(a1, a2, a3, n, rx, ry, rz, z, rc, rd, e)
!______________________________________________________________________________
!
implicit none
real(dp), intent(in) :: a1(3), a2(3), a3(3)
integer, intent(in) :: n
real(dp), intent(in) :: rx(n), ry(n), rz(n)
real(dp), intent(in) :: z(n)
real(dp), intent(in) :: rc
real(dp), intent(in) :: rd
real(dp), intent(out) :: e
real(dp) :: vol, rho, ei, qi, rij, a(3,3), bt(3,3), &
d_100, d_010, d_001, origin_j(3), xyz_ij(3), ra
integer :: i, j, shift1, shift2, shift3, &
shift1max, shift2max, shift3max
!______________________________________________________________________________
!
! compute cell volume and average density
vol = a1(1) * (a2(2)*a3(3) - a2(3)*a3(2)) + &
a1(2) * (a2(3)*a3(1) - a2(1)*a3(3)) + &
a1(3) * (a2(1)*a3(2) - a2(2)*a3(1))
vol = abs(vol) ! for left-handed coordinate systems
rho = sum(z) / vol
! compute reciprocal lattice vectors
a(:,1) = a1; a(:,2) = a2; a(:,3) = a3
call invert_3x3(a, bt) ! note: bt still missing factor of 2*pi
! compute distances between planes (accounts for missing 2*pi in bt)
d_100 = 1.0_dp / sqrt(sum(bt(1,:) * bt(1,:)))
d_010 = 1.0_dp / sqrt(sum(bt(2,:) * bt(2,:)))
d_001 = 1.0_dp / sqrt(sum(bt(3,:) * bt(3,:)))
! compute the number of cells to include along each direction
shift1max = ceiling(rc / d_100)
shift2max = ceiling(rc / d_010)
shift3max = ceiling(rc / d_001)
! prepare for loop over ions in cell
e = 0.0_dp
! loop over ions in cell
do i = 1, n
! prepare for loop over neighboring ions
ei = 0.0_dp
qi = z(i) ! b/c the i==j part of the sum is skipped below
! loop over cells
do shift3 = -shift3max, shift3max
do shift2 = -shift2max, shift2max
do shift1 = -shift1max, shift1max
! get origin of shifted cell and loop over ions in that cell
origin_j = shift1*a1 + shift2*a2 + shift3*a3
do j = 1, n
! exclude i==j
if (i==j .and. shift1==0 .and. shift2==0 .and. shift3==0) cycle
! compute distance between points
xyz_ij = (/rx(i), ry(i), rz(i)/) &
- (origin_j + (/rx(j), ry(j), rz(j)/))
rij = sqrt(sum(xyz_ij * xyz_ij))
! only proceed if rij < rc
if (rij > rc) cycle
! update energy and charge
ei = ei + z(j) * erfc(rij/rd) / rij
qi = qi + z(j)
end do ! j
end do ! shift1
end do ! shift2
end do ! shift3
! apply 1/2 z(i) factor to energy
ei = 0.5_dp * z(i) * ei
! add correction terms
ra = (3.0_dp * qi / (4.0_dp * pi * rho))**one_third ! adaptive cutoff
ei = ei - pi * z(i) * rho * ra * ra &
+ pi * z(i) * rho * (ra*ra - rd*rd/2.0_dp) * erf(ra/rd) &
+ sqrt_pi * z(i) * rho * ra * rd * exp(-ra*ra/(rd*rd)) &
- 1.0/(sqrt_pi * rd) * z(i) * z(i)
! increment the total energy
e = e + ei
end do ! i
end subroutine
subroutine force(a1, a2, a3, n, rx, ry, rz, z, rc, rd, fx, fy, fz)
!______________________________________________________________________________
!
implicit none
real(dp), intent(in) :: a1(3), a2(3), a3(3)
integer, intent(in) :: n
real(dp), intent(in) :: rx(n), ry(n), rz(n)
real(dp), intent(in) :: z(n)
real(dp), intent(in) :: rc
real(dp), intent(in) :: rd
real(dp), intent(out) :: fx(n), fy(n), fz(n)
real(dp) :: rij, rijrij, rij_rd, a(3,3), bt(3,3), &
d_100, d_010, d_001, origin_j(3), xyz_ij(3), t
integer :: i, j, shift1, shift2, shift3, &
shift1max, shift2max, shift3max
!______________________________________________________________________________
!
! compute reciprocal lattice vectors
a(:,1) = a1; a(:,2) = a2; a(:,3) = a3
call invert_3x3(a, bt) ! note: bt still missing factor of 2*pi
! compute distances between planes (accounts for missing 2*pi in bt)
d_100 = 1.0_dp / sqrt(sum(bt(1,:) * bt(1,:)))
d_010 = 1.0_dp / sqrt(sum(bt(2,:) * bt(2,:)))
d_001 = 1.0_dp / sqrt(sum(bt(3,:) * bt(3,:)))
! compute the number of cells to include along each direction
shift1max = ceiling(rc / d_100)
shift2max = ceiling(rc / d_010)
shift3max = ceiling(rc / d_001)
! prepare for loop over ions in cell
fx = 0.0_dp
fy = 0.0_dp
fz = 0.0_dp
! loop over ions in cell
do i = 1, n
! loop over cells
do shift3 = -shift3max, shift3max
do shift2 = -shift2max, shift2max
do shift1 = -shift1max, shift1max
! get origin of shifted cell and loop over ions in that cell
origin_j = shift1*a1 + shift2*a2 + shift3*a3
do j = 1, n
! exclude i==j
if (i==j .and. shift1==0 .and. shift2==0 .and. shift3==0) cycle
! compute distance between points
xyz_ij = (/rx(i), ry(i), rz(i)/) &
- (origin_j + (/rx(j), ry(j), rz(j)/))
rijrij = sum(xyz_ij * xyz_ij)
rij = sqrt(rijrij)
rij_rd = rij / rd
! only proceed if rij < rc
if (rij > rc) cycle
! compute common term
t = z(j) * &
(2.0_dp / sqrt_pi * rij_rd * exp(-rij_rd * rij_rd) &
+ erfc(rij_rd)) / (rij * rijrij)
! add contributions to forces
fx(i) = fx(i) + t * xyz_ij(1)
fy(i) = fy(i) + t * xyz_ij(2)
fz(i) = fz(i) + t * xyz_ij(3)
end do ! j
end do ! shift1
end do ! shift2
end do ! shift3
! apply z(i) factor
fx(i) = z(i) * fx(i)
fy(i) = z(i) * fy(i)
fz(i) = z(i) * fz(i)
end do ! i
end subroutine
subroutine stress(a1, a2, a3, n, rx, ry, rz, z, rc, rd, s)
!______________________________________________________________________________
!
implicit none
real(dp), intent(in) :: a1(3), a2(3), a3(3)
integer, intent(in) :: n
real(dp), intent(in) :: rx(n), ry(n), rz(n)
real(dp), intent(in) :: z(n)
real(dp), intent(in) :: rc
real(dp), intent(in) :: rd
real(dp), intent(out) :: s(6)
real(dp) :: vol, rho, qi, rij, rijrij, rij_rd, a(3,3), bt(3,3), &
d_100, d_010, d_001, origin_j(3), xyz_ij(3), t, si(6), &
ra, ra_rd
integer :: i, j, shift1, shift2, shift3, &
shift1max, shift2max, shift3max
!______________________________________________________________________________
!
! compute cell volume and average density
vol = a1(1) * (a2(2)*a3(3) - a2(3)*a3(2)) + &
a1(2) * (a2(3)*a3(1) - a2(1)*a3(3)) + &
a1(3) * (a2(1)*a3(2) - a2(2)*a3(1))
vol = abs(vol) ! for left-handed coordinate systems
rho = sum(z) / vol
! compute reciprocal lattice vectors
a(:,1) = a1; a(:,2) = a2; a(:,3) = a3
call invert_3x3(a, bt) ! note: bt still missing factor of 2*pi
! compute distances between planes (accounts for missing 2*pi in bt)
d_100 = 1.0_dp / sqrt(sum(bt(1,:) * bt(1,:)))
d_010 = 1.0_dp / sqrt(sum(bt(2,:) * bt(2,:)))
d_001 = 1.0_dp / sqrt(sum(bt(3,:) * bt(3,:)))
! compute the number of cells to include along each direction
shift1max = ceiling(rc / d_100)
shift2max = ceiling(rc / d_010)
shift3max = ceiling(rc / d_001)
! prepare for loop over ions in cell
s = 0.0_dp
! loop over ions in cell
do i = 1, n
! prepare for loop over neighboring ions
si = 0.0_dp
qi = z(i) ! b/c the i==j part of the sum is skipped below
! loop over cells
do shift3 = -shift3max, shift3max
do shift2 = -shift2max, shift2max
do shift1 = -shift1max, shift1max
! get origin of shifted cell and loop over ions in that cell
origin_j = shift1*a1 + shift2*a2 + shift3*a3
do j = 1, n
! exclude i==j
if (i==j .and. shift1==0 .and. shift2==0 .and. shift3==0) cycle
! compute distance between points
xyz_ij = (/rx(i), ry(i), rz(i)/) &
- (origin_j + (/rx(j), ry(j), rz(j)/))
rijrij = sum(xyz_ij * xyz_ij)
rij = sqrt(rijrij)
rij_rd = rij / rd
! only proceed if rij < rc
if (rij > rc) cycle
! compute common term
t = z(j) * &
(2.0_dp / sqrt_pi * rij_rd * exp(-rij_rd * rij_rd) &
+ erfc(rij_rd)) / (rij * rijrij)
! add contributions to stresses and charge
si(1:3) = si(1:3) + t * xyz_ij(1:3) * xyz_ij(1:3)
si(4) = si(4) + t * xyz_ij(2) * xyz_ij(3)
si(5) = si(5) + t * xyz_ij(1) * xyz_ij(3)
si(6) = si(6) + t * xyz_ij(1) * xyz_ij(2)
qi = qi + z(j)
end do ! j
end do ! shift1
end do ! shift2
end do ! shift3
! apply factor of -z(i) / (2 * volume)
si = -z(i) / (2.0_dp * vol) * si
! add remaining terms
ra = (3.0_dp * qi / (4.0_dp * pi * rho))**one_third ! adaptive cutoff
ra_rd = ra / rd
t = pi / (2.0_dp * vol) * rd * rd * rho * z(i) &
* (1.0_dp - 2.0/sqrt_pi * ra_rd * exp(-ra_rd * ra_rd) &
+ (2.0_dp / 3.0_dp * ra_rd * ra_rd - 1.0_dp) * erfc(ra_rd))
si(1:3) = si(1:3) + t
! increment the total stress
s = s + si
end do ! i
end subroutine
subroutine invert_3x3(a, b)
!______________________________________________________________________________
!
! inverts a 3x3 matrix. not hyper-optimized.
!______________________________________________________________________________
!
implicit none
real(dp), intent(in) :: a(3,3)
real(dp), intent(out) :: b(3,3)
real(dp) :: det
!______________________________________________________________________________
!
! compute cofactor matrix
b(1,1) = a(2,2)*a(3,3) - a(3,2)*a(2,3)
b(2,1) = a(2,3)*a(3,1) - a(3,3)*a(2,1)
b(3,1) = a(2,1)*a(3,2) - a(3,1)*a(2,2)
b(1,2) = a(1,3)*a(3,2) - a(3,3)*a(1,2)
b(2,2) = a(1,1)*a(3,3) - a(3,1)*a(1,3)
b(3,2) = a(1,2)*a(3,1) - a(3,2)*a(1,1)
b(1,3) = a(1,2)*a(2,3) - a(2,2)*a(1,3)
b(2,3) = a(1,3)*a(2,1) - a(2,3)*a(1,1)
b(3,3) = a(1,1)*a(2,2) - a(2,1)*a(1,2)
! compute determinant and divide
det = a(1,1)*(a(2,2)*a(3,3) - a(3,2)*a(2,3))
det = det - a(1,2)*(a(2,1)*a(3,3) - a(3,1)*a(2,3))
det = det + a(1,3)*(a(2,1)*a(3,2) - a(3,1)*a(2,2))
b = b / det
end subroutine
end module