forked from aromanro/PythonCompphys
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathhartree-fock.py
529 lines (324 loc) · 15 KB
/
hartree-fock.py
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
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
#!/usr/bin/env python
# coding: utf-8
# In[1]:
import math as m
import numpy as np
from scipy.linalg import eigh
import matplotlib.pyplot as plt
get_ipython().run_line_magic('matplotlib', 'inline')
# ### References
#
# The referenced chapters, formulae and problems are from the book [Computational Physics, by Jos Thijssen](https://www.cambridge.org/core/books/computational-physics/BEE73B0139D4A9993193B57CDC62096E#fndtn-information).
#
# This is only introductory material for more advanced topics/code presented on the [Computational Physics Blog](https://compphys.go.ro).
#
# Related with Hartree-Fock I would recommend visiting:
# [How to solve a quantum many body problem](https://compphys.go.ro/how-to-solve-a-quantum-many-body-problem/) - on Born-Oppenheimer approximation, diagonalization and variational primciple (on this one you may find several other topics on the blog).
#
# [The Hartree-Fock program](https://compphys.go.ro/the-hartree-fock-program/) - describes a C++ project I have on GitHub: [HartreeFock](https://github.com/aromanro/HartreeFock) which goes way beyond the following simple issues.
#
#
# ### From here, the variational Hydrogen computation, chapter 3
#
# We're starting using Gaussians as basis functions. Here only $\Phi(r)=e^{-\alpha r^2}$. For the more complex ones and the details on the many centers cases, please check out the blog and the C++ Hartree-Fock project.
# In[2]:
# 3.26, 4.16
def Gaussian(alpha, r):
return m.exp(-alpha*r*r)
# In[3]:
alpha=(13.00773, 1.962079, 0.444529, 0.1219492) # from 3.27
#alpha=(14, 2, 0.5, 0.1) # try something different here to see how it affects the results
# Let's see how they look like:
# In[4]:
f = np.zeros((100, 4))
r = np.linspace(0, 4, 100)
for i in range(100):
for j in range(4):
f[i, j] = Gaussian(alpha[j], r[i])
plt.plot(r, f)
plt.show()
# Integrals are very simple, the gaussians are centered in the same point (the nucleus), see eqn 3.29:
# In[5]:
def Overlap(alpha, ind1, ind2):
return m.pow(m.pi/(alpha[ind1] + alpha[ind2]), 3./2.)
# In[6]:
def Kinetic(alpha, ind1, ind2):
return 3*alpha[ind1]*alpha[ind2]*m.pow(m.pi, 3./2.)/m.pow(alpha[ind1] + alpha[ind2],5./2.)
# In[7]:
def Coulomb(alpha, ind1, ind2):
return -2 * m.pi / (alpha[ind1] + alpha[ind2])
# The basis functions are not orthogonal, so we have to solve the generalized eigenvalue problem, using the overlap:
# In[8]:
basisSize = 4
# In[9]:
H = np.zeros((basisSize, basisSize))
O = np.zeros((basisSize, basisSize))
for i in range(basisSize):
H[i, i] = Kinetic(alpha, i, i) + Coulomb(alpha, i, i)
O[i, i] = Overlap(alpha, i, i)
for j in range(i):
H[i, j] = Kinetic(alpha, i, j) + Coulomb(alpha, i, j)
O[i, j] = Overlap(alpha, i, j)
H[j, i] = H[i, j]
O[j, i] = O[i, j]
# In[10]:
eigvals, eigvecs = eigh(H, O, eigvals_only=False)
# In[11]:
print(eigvals[0])
# The exact result is -0.5 (here we use atomic units, the energy is in Hartrees).
#
# The result here is quite good (you can often do that if you know the exact result in advance). DFTAtom gets -0.445 (that's the expected result for LDA, check out [NIST](https://www.nist.gov/pml/atomic-reference-data-electronic-structure-calculations/atomic-reference-data-electronic-7-0)) and HartreeFock (unrestricted) with STO6G gets -0.471. With the Variational Quantum Monte Carlo project with STO6G I also got -0.471. As the problem is analytically solvable and one knows the result, it's no surprise that the result is so good. One just has to ensure that the space spanned by the basis vectors includes as much as possible of the real solution.
#
# Now let's look at the computed wavefunction together with the exact solution:
# In[12]:
def Wavefunction(alpha, eigenvecs, r, ind = 0):
v = 0.
for i in range(basisSize):
v += eigenvecs[i, ind] * Gaussian(alpha[i], r)
return v
# In[13]:
r = np.linspace(0, 2, 100)
w = np.zeros((100, 2))
for i in range(100):
w[i, 0] = np.abs(Wavefunction(alpha, eigvecs, r[i]))
w[i, 1] = 1. / np.sqrt(np.pi) * np.exp(-r[i])
# In[14]:
plt.plot(r, w)
plt.show()
# ### Chap 4, Helium computation
# In[15]:
alpha = (0.298073, 1.242567, 5.782948, 38.47497)
# In[16]:
# 4.11
def TwoElectronSingleCenter(alpha, p, r, q, s):
return 2. * m.pow(m.pi, 5. / 2.) / ((alpha[p]+alpha[q])*(alpha[r]+alpha[s])*m.sqrt(alpha[p]+alpha[q]+alpha[r]+alpha[s]))
# In[17]:
H = np.zeros((basisSize, basisSize))
O = np.zeros((basisSize, basisSize))
Q = np.zeros((basisSize, basisSize, basisSize, basisSize))
for i in range(basisSize):
for j in range(basisSize):
H[i, j] = Kinetic(alpha, i, j) + 2. * Coulomb(alpha, i, j) #the 2 is due of Z=2 for Helium
O[i, j] = Overlap(alpha, i, j)
for k in range(basisSize):
for l in range(basisSize):
Q[i, j, k, l]=TwoElectronSingleCenter(alpha, i, j, k, l)
# In[18]:
v = 1. / m.sqrt(O.sum())
C = np.array([v, v, v, v]) # a choice for C to start with. Check the commented one instead
#C = np.array([1, 1, 1, 1])
# In[19]:
F = np.zeros((basisSize, basisSize))
oldE = 100
for cycle in range(100):
# 4.18
#for i in range(basisSize):
# for j in range(basisSize):
# F[i, j] = H[i, j]
# for k in range(basisSize):
# for l in range(basisSize):
# F[i, j] += Q[i, k, j, l] * C[k] * C[l]
F = H + np.einsum('ikjl,k,l', Q, C, C)
# 4.20
eigvals, eigvecs = eigh(F, O, eigvals_only=False)
C = eigvecs[:,0]
# 4.21
#Eg = 0
#for i in range(basisSize):
# for j in range(basisSize):
# Eg += 2 * C[i] * C[j] * H[i, j]
# for k in range(basisSize):
# for l in range(basisSize):
# Eg += Q[i, k, j, l] * C[i] * C[j] * C[k] * C[l]
Eg = 2 * np.einsum('ij,i,j', H, C, C) + np.einsum('ikjl,i,j,k,l', Q, C, C, C, C)
if abs(oldE-Eg) < 1E-10:
break
oldE = Eg
# In[20]:
print(Eg)
# HartreeFock with STO6G gets here -2.846 (but with a basis set with more Gaussians, you certainly can get better results), DFTAtom -2.8348 (again, that's normal for LDA: [NIST](https://www.nist.gov/pml/atomic-reference-data-electronic-structure-calculations/atomic-reference-data-electronic-7-1)) and the Variational Quantum Monte Carlo, -2.8759. The later 'beats' the result obtained here.
# ### The first problem, H2+ (problem 4.8)
#
# Now we have two centers so intergrals get more complex. You may simplify some integrals computation by using the already computed overlap, such optimizations are left out from here, they exist in the C++ project. For the H2+, as there is a single electron, we don't need a self-consistency loop.
# In[21]:
# see 4.114 and 4.116
def F0(t):
if t==0:
return 1.
p = m.sqrt(t)
a = 1. / p
return a * m.sqrt(m.pi) / 2. * m.erf(p)
# In[22]:
# 4.98
def Rp(alpha, beta, Ra, Rb):
return (alpha*Ra + beta*Rb) / (alpha + beta)
# In[23]:
# 4.100
def OverlapTwoCenters(alpha, beta, Ra, Rb):
difR = Ra - Rb
len2 = difR.dot(difR)
aplusb = alpha + beta
ab = alpha * beta / aplusb
return m.pow(m.pi / aplusb, 3./2.) * m.exp(-ab * len2)
# In[24]:
# 4.103
def KineticTwoCenters(alpha, beta, Ra, Rb):
difR = Ra - Rb
len2 = difR.dot(difR)
aplusb = alpha + beta
ab = alpha * beta / aplusb
O = m.pow(m.pi/aplusb, 3./2.) * m.exp(-ab * len2) # it's actually the overlap, check the OverlapTwoCenters
return ab * (3. - 2. * ab * len2) * O #this can be optimized with already computed overlap, see above
# In[25]:
# 4.115
def Nuclear(alpha, beta, Ra, Rb, Rc, Z = 1.):
aplusb = alpha + beta
ab = alpha * beta / aplusb
difR = Ra - Rb
len2 = difR.dot(difR)
difRc = Rp(alpha, beta, Ra, Rb) - Rc
len2c = difRc.dot(difRc)
K = m.exp(-ab*len2)
return -2. * m.pi * Z / aplusb * K * F0(aplusb*len2c)
# In[26]:
# 4.123
def TwoElectronTwoCenter(alpha, beta, gamma, delta, Ra, Rb, Rc, Rd):
RP = Rp(alpha, gamma, Ra, Rc)
RQ = Rp(beta, delta, Rb, Rd)
alphaplusgamma = alpha + gamma
betaplusdelta = beta + delta
Rac = Ra - Rc
Rbd = Rb - Rd
Rpq = RP - RQ
Racl2 = Rac.dot(Rac)
Rbdl2 = Rbd.dot(Rbd)
Rpql2 = Rpq.dot(Rpq)
return 2. * m.pow(m.pi, 5./2.) / (alphaplusgamma * betaplusdelta * m.sqrt(alphaplusgamma+betaplusdelta)) * m.exp(-alpha*gamma/alphaplusgamma*Racl2 - beta*delta/betaplusdelta*Rbdl2) * F0(alphaplusgamma*betaplusdelta / (alphaplusgamma+betaplusdelta) * Rpql2)
# In[27]:
alpha=(13.00773, 1.962079, 0.444529, 0.1219492)
# In[28]:
basisSize = 4 # for each atom
# In[29]:
H = np.zeros((basisSize * 2, basisSize * 2))
O = np.zeros((basisSize * 2, basisSize * 2))
R0 = np.array([0, 0, 0])
R1 = np.array([1, 0, 0])
for i in range(basisSize):
a = alpha[i]
basisSizei = basisSize + i
for j in range(basisSize):
b = alpha[j]
basisSizej = basisSize + j
O[i, j] = OverlapTwoCenters(a, b, R0, R0)
O[basisSizei, j] = OverlapTwoCenters(a, b, R1, R0)
O[i, basisSizej] = OverlapTwoCenters(a, b, R0, R1)
O[basisSizei, basisSizej] = OverlapTwoCenters(a, b, R1, R1)
H[i, j] = KineticTwoCenters(a, b, R0, R0) + Nuclear(a, b, R0, R0, R0) + Nuclear(a, b, R0, R0, R1)
H[basisSizei, j] = KineticTwoCenters(a, b, R1, R0) + Nuclear(a, b, R1, R0, R0) + Nuclear(a, b, R1, R0, R1)
H[i, basisSizej] = KineticTwoCenters(a, b, R0, R1) + Nuclear(a, b, R0, R1, R0) + Nuclear(a, b, R0, R1, R1)
H[basisSizei, basisSizej] = KineticTwoCenters(a, b, R1, R1) + Nuclear(a, b, R1, R1, R0) + Nuclear(a, b, R1, R1, R1)
# In[30]:
eigvals, eigvecs = eigh(H, O, eigvals_only=False)
# In[31]:
print(eigvals)
# ### Now the next problem, for H2 (problem 4.9) - here only the Hartree equation is solved
#
# Filling the Q tensor looks ugly but I think it's better to leave it like that. If you want, you may make it prettier by observing that you can count from 0 to 0xF and use the proper bit to select to pass R0 or R1 and to select the proper sector of the tensor.
# In[32]:
Q = np.zeros((basisSize*2, basisSize*2, basisSize*2, basisSize*2))
for i in range(basisSize):
a = alpha[i]
basisSizei = basisSize + i
for j in range(basisSize):
b = alpha[j]
basisSizej = basisSize + j
for k in range(basisSize):
c = alpha[k]
basisSizek = basisSize + k
for l in range(basisSize):
basisSizel = basisSize + l
d = alpha[l]
Q[i, j, k, l]=TwoElectronTwoCenter(a, b, c, d, R0, R0, R0, R0)
Q[i, j, k, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R0, R0, R0, R1)
Q[i, j, basisSizek, l]=TwoElectronTwoCenter(a, b, c, d, R0, R0, R1, R0)
Q[i, j, basisSizek, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R0, R0, R1, R1)
Q[i, basisSizej, k, l]=TwoElectronTwoCenter(a, b, c, d, R0, R1, R0, R0)
Q[i, basisSizej, k, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R0, R1, R0, R1)
Q[i, basisSizej, basisSizek, l]=TwoElectronTwoCenter(a, b, c, d, R0, R1, R1, R0)
Q[i, basisSizej, basisSizek, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R0, R1, R1, R1)
Q[basisSizei, j, k, l]=TwoElectronTwoCenter(a, b, c, d, R1, R0, R0, R0)
Q[basisSizei, j, k, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R1, R0, R0, R1)
Q[basisSizei, j, basisSizek, l]=TwoElectronTwoCenter(a, b, c, d, R1, R0, R1, R0)
Q[basisSizei, j, basisSizek, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R1, R0, R1, R1)
Q[basisSizei, basisSizej, k, l]=TwoElectronTwoCenter(a, b, c, d, R1, R1, R0, R0)
Q[basisSizei, basisSizej, k, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R1, R1, R0, R1)
Q[basisSizei, basisSizej, basisSizek, l]=TwoElectronTwoCenter(a, b, c, d, R1, R1, R1, R0)
Q[basisSizei, basisSizej, basisSizek, basisSizel]=TwoElectronTwoCenter(a, b, c, d, R1, R1, R1, R1)
# In[33]:
v = 1. / m.sqrt(O.sum())
C = np.array([v, v, v, v, v, v, v, v])
# In[34]:
F = np.zeros((2*basisSize, 2*basisSize))
oldE = 100
for cycle in range(100):
#for i in range(2*basisSize):
# for j in range(2*basisSize):
# F[i, j] = H[i, j]
# for k in range(2*basisSize):
# for l in range(2*basisSize):
# F[i, j] += Q[i, k, j, l] * C[k] * C[l]
F = H + np.einsum('ikjl,k,l', Q, C, C)
eigvals, eigvecs = eigh(F, O, eigvals_only=False)
C = eigvecs[:,0]
#Eg = 0
#for i in range(2*basisSize):
# for j in range(2*basisSize):
# Eg += 2 * C[i] * C[j] * H[i, j]
# for k in range(2*basisSize):
# for l in range(2*basisSize):
# Eg += Q[i, k, j, l] * C[i] * C[j] * C[k] * C[l]
Eg = C.dot(H + F).dot(C)
if abs(oldE-Eg) < 1E-10:
break
oldE = Eg
# In[35]:
print(Eg + 1) # +1 for nuclear repulsion energy
# In[36]:
print(C)
# ### Now the next problem, for H2 (problem 4.12) - but now with Hartree-Fock
#
# Here there are a lot of optimizations suggested in the book - for the electron-electron integrals, one can reduce roughly to 1/8 the number of computed integrals using symmetry. Instead of solving the general eigenvalue problem each selfconsistency step you can go with solving only the regular eigenvalue problem, for details please check the book and the C++ HartreeFock project. I didn't want to go with those here.
# In[37]:
Qt = np.zeros((basisSize*2, basisSize*2, basisSize*2, basisSize*2))
for p in range(2*basisSize):
for q in range(2*basisSize):
for r in range(2*basisSize):
for s in range(2*basisSize):
Qt[p, q, r, s] = 2. * Q[p, q, r, s] - Q[p, r, s, q]
# In[38]:
C = np.array([v, v, v, v, v, v, v, v]) #reinitialize C
# In[39]:
# as above, but Hartree-Fock with Qt instead of Q
F = np.zeros((2*basisSize, 2*basisSize))
oldE = 100
for cycle in range(100):
#for i in range(2*basisSize):
# for j in range(2*basisSize):
# F[i, j] = H[i, j]
# for k in range(2*basisSize):
# for l in range(2*basisSize):
# F[i, j] += Qt[i, k, j, l] * C[k] * C[l]
F = H + np.einsum('ikjl,k,l', Q, C, C)
eigvals, eigvecs = eigh(F, O, eigvals_only=False)
C = eigvecs[:,0]
Eg = C.dot(H + F).dot(C)
if abs(oldE-Eg) < 1E-10:
break
oldE = Eg
# In[40]:
print(Eg + 1) # +1 for nuclear repulsion energy
# We obtain the same result as above, despite using Hartree-Fock, because we have only two electrons, one with spin up, one with down, so no spin-exchange.
#
# With the C++ Hartree-Fock project, restricted method with STO6G, the result is -1.0758. Again, with a better basis set the result could be much better. With the Variational Quantum Monte Carlo method I got -1.092, again, better than the result above.
# In[41]:
print(C)