-
Notifications
You must be signed in to change notification settings - Fork 7
/
AtomicOrbitals.py
1325 lines (1084 loc) · 45.3 KB
/
AtomicOrbitals.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
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
import numpy as np
from math import pi, log, factorial
from literature_data.parser import parse
import os
try:
import matplotlib.pyplot as plt
from matplotlib.colors import LinearSegmentedColormap
MAKE_MPL_COLOR_MAP = True
except ImportError:
MAKE_MPL_COLOR_MAP = False
try:
import pylibxc
HAVE_LIBXC = True
except ImportError:
HAVE_LIBXC = False
try:
import xcfun
HAVE_XCFUN = True
except ImportError:
HAVE_XCFUN = False
try:
import scipy
HAVE_SCIPY = True
except ImportError:
HAVE_SCIPY = False
"""Code to calculate the electron density, its gradient and laplacian,
as well as the local kinetic energy density from analytic Slater-type
orbital expansions of atomic Hartree-Fock wave functions published in
the literature.
Although the best-known Hartree-Fock wave functions for atoms are the
ones published by Clementi and Roetti in 1974, their wave functions
have significant errors (up to millihartree), especially for the
heavier atoms. Instead, the wave functions included in this library
for H-Xe and their cations and anions are from
Koga et al, "Analytical Hartree–Fock wave functions subject to cusp
and asymptotic constraints: He to Xe, Li+ to Cs+, H− to I−",
Int. J. Quantum Chem. 71, 491 (1999),
doi: 10.1002/(SICI)1097-461X(1999)71:6<491::AID-QUA6>3.0.CO;2-T
which yield total eneregies that deviate from the exact numerical
solution only by microhartrees. Included are also non-relativistic
Hartree-Fock wave functions for heavier atoms from
Koga et al, "Analytical Hartree–Fock wave functions for the atoms Cs
to Lr", Theor. Chem. Acc. 104, 411 (2000), doi: 10.1007/s002140000150
which have an average error of 0.19 mEh compared to exact, fully
numerical calculations.
The core functionality of the library comes from the Atom() class
which implements the get_densities(r) function to return the
above-mentioned quantities at a given numpy array of radial points r.
The Hartree-Fock orbital data are stored in the static AtomData class,
which implements dictionaries of parameters for various elements that
are accessed by upper-case element symbols.
The GridGenerator gives a simple implementation of an efficient
Gauss-Legendre quadrature grid.
A color dictionary following colors used by Jmol (roughly CPK) by
element symbol. Should be accessed using get_colors_for_elements()
which returns a list of RGB for an input list of string element labels
or Atom objects.
Module is Python 2 and 3 compatible.
Examples of use can be seen in test_densities() function.
"""
def test_densities():
"""
Simple unit testing routine.
All functions were validated against the reference FORTRAN code to good agreement.
However, it is possible to test the gradient by finite difference, so lets do that.
The Laplacian does not seem so well behaved (but agrees with the reference).
This should not need re-testing as it is independent of input data, but it can't hurt.
As long as reference total electron number and kinetic energies are included
for the atom specification we can check the Slater orbital parameters that have been
entered by integrating the density and kinetic energy density. Good agreement
indicates parameters have been entered correctly.
Whilst the threshold value on these integrations may seem a bit high,
painstaking cross checking shows that even with correct values
and converged grids, this is as close as we get. Most small changes seem to break it
so we can have confidence.
"""
actual, r, wt = GridGenerator.make_grid(400, 'ahlrichsm3')
grid = 4*pi*wt
data = AtomData()
print("\nINTEGRATED DENSITY TEST")
print("=======================")
for a in list(data.nuclear_charge.keys()):
atom = Atom(a)
Nel = data.electron_count[a]
d0, d1, g0, g1, t0, t1, l0, l1 = atom.get_densities(r)
# Count electrons per spin channel
s_occ = AtomData.s_occ.get(a, [0, 0])
p_occ = AtomData.p_occ.get(a, [0, 0])
d_occ = AtomData.d_occ.get(a, [0, 0])
f_occ = AtomData.f_occ.get(a, [0, 0])
nela = np.sum(s_occ[0])+np.sum(p_occ[0])+np.sum(d_occ[0])+np.sum(f_occ[0])
nelb = np.sum(s_occ[1])+np.sum(p_occ[1])+np.sum(d_occ[1])+np.sum(f_occ[1])
assert(nela+nelb == Nel)
id0 = np.dot(d0, grid)
id1 = np.dot(d1, grid)
diff_0 = id0 - nela
percent_diff_0 = 100*diff_0/nela
# Check to catch for Hydrogen having no beta electrons
if nelb > 0.0:
diff_1 = id1 - nelb
percent_diff_1 = 100*diff_1/nelb
else:
diff_1 = 0.0
percent_diff_1 = 0.0
print("{:>3} - N_0 = ({:4.1f}) {:+2.6e}%, N_1 = ({:4.1f}) {:+2.6e}%, {:}".format(a, id0, percent_diff_0, id1, percent_diff_1, "PASSED" if max(abs(diff_0), abs(diff_1)) < 1e-4 else "FAILED - "))
print("\nINTEGRATED KINETIC TEST")
print("=======================")
for a in list(data.ke_test.keys()):
atom = Atom(a)
t_bm = data.ke_test[a]
d0, d1, g0, g1, t0, t1, l0, l1 = atom.get_densities(r)
it0 = np.dot(t0, grid)
it1 = np.dot(t1, grid)
itot = it0 + it1
diff = itot - t_bm
print("{:>3} - T = {:+.6e}%, {:}".format(a, 100*diff/t_bm, "PASSED" if abs(100*diff/t_bm) < 1e-2 else "FAILED - "))
# The integral of the Laplacian over all space should be 0. Check that.
print("\nINTEGRATED LAPLACIAN TEST")
print("=========================")
for a in list(AtomData.ke_test.keys()):
atom = Atom(a)
d0, d1, g0, g1, t0, t1, l0, l1 = atom.get_densities(r)
il0 = np.dot(grid, l0)
il1 = np.dot(grid, l1)
print("{:>3} - L_0 = {:+.6e}, L_1 = {:+.6e}, {:}".format(a, il0, il1, "PASSED" if max(abs(il0), abs(il1)) < 1e-6 else "FAILED - "))
print("\nFINITE DIFFERENCE GRADIENT TEST")
print("===============================")
print("Testing gradient evaluation function against finite difference estimate...")
ne = Atom("Ne") # Let's use "the guvnor"
# We only need to test a few points around the core
fdh = 1e-8
fdr = np.arange(0.9, 0.9+fdh*10, fdh)
d0, d1, g0, g1, t0, t1, l0, l1 = ne.get_densities(fdr)
# First the first central difference
fdiff = (d0[2:]-d0[:-2])/(2*fdh) # Construct the central difference
if np.allclose(fdiff, g0[1:-1], atol=1e-1): # finite difference is not perfect, so lenient tolerance
print("Gradient: PASSED")
else:
print("Gradient: FAILED -")
print("\nELEMENT COLOR FUNCTIONS TEST")
print("===========================")
test_obj = [Atom("H"), Atom("C"), Atom("O")]
test_str = ["H", "C", "O"]
ref = np.array([[1., 1., 1.], [0.565, 0.565, 0.565], [1. , 0.051, 0.051]])
if np.allclose( np.array(get_colors_for_elements(test_obj)), ref):
print("\nColor from objects: PASSED")
else:
print("\nColor from objects: FAILED -")
if np.allclose( np.array(get_colors_for_elements(test_str)), ref):
print("Color from strings: PASSED")
else:
print("Color from strings: FAILED -")
if HAVE_LIBXC:
test_functional='GGA_X_PBE'
print("\nATOMIC EXCHANGE ENERGIES WITH {}".format(test_functional))
print("============================================")
for a in list(data.ke_test.keys()):
atom = Atom(a)
nE, vrho, vsigma, vtau, vlapl = atom.libxc_eval(r, functional=test_functional, restricted=False)
Exc = (np.dot(nE, grid)).item()
print('{:3s} {:.10f}'.format(a, Exc))
else:
print("\nNot doing energy calculations due to lack of libxc.\n")
class Atom:
"""
Object initialised to contain the relevant data
for the element provided to the constructor.
The returned Atom object then has the core method
get_densities(r)
that takes a distance (or numpy array of distances)
and returns the corresponding densities.
"""
def __init__(self, element):
# Currently only deal with atomic symbols.
# Eventually add dictionary to support full names and numbers?
self.element = element
try:
u_atom = element.upper()
d_atom = element.lower()
# set nuclear charge
self.nuclear_charge = AtomData.nuclear_charge[u_atom]
# We assume that every atom has S electrons
# So access with index notation to raise KeyError on missing atom.
self.s_exp = AtomData.s_exp[u_atom]
self.s_coef = AtomData.s_coef[u_atom]
self.s_n = AtomData.s_n[u_atom]
self.s_occ = AtomData.s_occ[u_atom]
# It is possible to be missing P and D electrons, so use get with defaults
self.p_exp = AtomData.p_exp.get(u_atom, None)
self.p_coef = AtomData.p_coef.get(u_atom, None)
self.p_n = AtomData.p_n.get(u_atom, None)
self.p_occ = AtomData.p_occ.get(u_atom, [0, 0])
self.d_exp = AtomData.d_exp.get(u_atom, None)
self.d_coef = AtomData.d_coef.get(u_atom, None)
self.d_n = AtomData.d_n.get(u_atom, None)
self.d_occ = AtomData.d_occ.get(u_atom, [0, 0])
self.f_exp = AtomData.f_exp.get(u_atom, None)
self.f_coef = AtomData.f_coef.get(u_atom, None)
self.f_n = AtomData.f_n.get(u_atom, None)
self.f_occ = AtomData.f_occ.get(u_atom, [0, 0])
except KeyError:
raise KeyError('Error: Atom data for "{:}" missing'.format(element))
def get_densities(self, r):
"""
Core function returning all densities for the initialised atom.
Makes use of numpy vectorisation to efficiently compute an array
of r.
In:
r space coordinate (distance from nucleus) [numpy array or scalar]
Out (Throughout X denotes spin 0/1 (up/down). All numpy arrays.):
den_X Electron density
grd_X Gradient of electron density
tau_X Orbital kinetic energy density
lap_X Laplacian of electron density
"""
# Handle scalar r values as single element array
if np.isscalar(r):
r = np.array([r])
assert np.min(r) > 0, "Error: distances must be non-zero and positive."
if self.s_exp is not None:
oS, doS, ddoS = self.get_orbitals(self.s_n, self.s_exp, self.s_coef, r)
den_0 = np.sum(self.s_occ[0][:,None]*oS**2, axis=0)
den_1 = np.sum(self.s_occ[1][:,None]*oS**2, axis=0)
grd_0 = np.sum(self.s_occ[0][:,None]*(oS*doS), axis=0)
grd_1 = np.sum(self.s_occ[1][:,None]*(oS*doS), axis=0)
tau_0 = np.sum(self.s_occ[0][:,None]*doS**2, axis=0)
tau_1 = np.sum(self.s_occ[1][:,None]*doS**2, axis=0)
lap_s = oS*ddoS + doS**2 + 2*oS*doS/r
lap_0 = np.sum(self.s_occ[0][:,None]*lap_s, axis=0)
lap_1 = np.sum(self.s_occ[1][:,None]*lap_s, axis=0)
else:
# Otherwise supply zeros in place
den_0 = np.zeros(r.shape)
den_1 = np.zeros(r.shape)
grd_0 = np.zeros(r.shape)
grd_1 = np.zeros(r.shape)
tau_0 = np.zeros(r.shape)
tau_1 = np.zeros(r.shape)
lap_0 = np.zeros(r.shape)
lap_1 = np.zeros(r.shape)
# Check if atom has occupied P orbitals
if self.p_exp is not None:
oP, doP, ddoP = self.get_orbitals(self.p_n, self.p_exp, self.p_coef, r)
den_0 += np.sum(self.p_occ[0][:,None]*oP**2, axis=0)
den_1 += np.sum(self.p_occ[1][:,None]*oP**2, axis=0)
grd_0 += np.sum(self.p_occ[0][:,None]*oP*doP, axis=0)
grd_1 += np.sum(self.p_occ[1][:,None]*oP*doP, axis=0)
tau_0 += np.sum(self.p_occ[0][:,None]*(doP**2 + 2*(oP/r)**2), axis=0)
tau_1 += np.sum(self.p_occ[1][:,None]*(doP**2 + 2*(oP/r)**2), axis=0)
lap_p = oP*ddoP + doP**2 + 2*oP*doP/r
lap_0 += np.sum(self.p_occ[0][:,None]*lap_p, axis=0)
lap_1 += np.sum(self.p_occ[1][:,None]*lap_p, axis=0)
# Check if atom has occupied D orbitals
if self.d_exp is not None:
oD, doD, ddoD = self.get_orbitals(self.d_n, self.d_exp, self.d_coef, r)
den_0 += np.sum(self.d_occ[0][:,None]*oD**2, axis=0)
den_1 += np.sum(self.d_occ[1][:,None]*oD**2, axis=0)
grd_0 += np.sum(self.d_occ[0][:,None]*oD*doD, axis=0)
grd_1 += np.sum(self.d_occ[1][:,None]*oD*doD, axis=0)
tau_0 += np.sum(self.d_occ[0][:,None]*(doD**2 + 6*(oD/r)**2), axis=0)
tau_1 += np.sum(self.d_occ[1][:,None]*(doD**2 + 6*(oD/r)**2), axis=0)
lap_d = oD*ddoD + doD**2 + 2*oD*doD/r
lap_0 += np.sum(self.d_occ[0][:,None]*lap_d, axis=0)
lap_1 += np.sum(self.d_occ[1][:,None]*lap_d, axis=0)
# Check if atom has occupied F orbitals
if self.f_exp is not None:
oF, doF, ddoF = self.get_orbitals(self.f_n, self.f_exp, self.f_coef, r)
den_0 += np.sum(self.f_occ[0][:,None]*oF**2, axis=0)
den_1 += np.sum(self.f_occ[1][:,None]*oF**2, axis=0)
grd_0 += np.sum(self.f_occ[0][:,None]*oF*doF, axis=0)
grd_1 += np.sum(self.f_occ[1][:,None]*oF*doF, axis=0)
tau_0 += np.sum(self.f_occ[0][:,None]*(doF**2 + 12*(oF/r)**2), axis=0)
tau_1 += np.sum(self.f_occ[1][:,None]*(doF**2 + 12*(oF/r)**2), axis=0)
lap_f = oF*ddoF + doF**2 + 2*oF*doF/r
lap_0 += np.sum(self.f_occ[0][:,None]*lap_f, axis=0)
lap_1 += np.sum(self.f_occ[1][:,None]*lap_f, axis=0)
# Take care of scaling
den_0 /= 4*pi
den_1 /= 4*pi
grd_0 /= 2*pi
grd_1 /= 2*pi
tau_0 /= 8*pi
tau_1 /= 8*pi
lap_0 /= 2*pi
lap_1 /= 2*pi
return den_0, den_1, grd_0, grd_1, tau_0, tau_1, lap_0, lap_1
def get_range(self):
"""
Evaluates the maximum range of the basis, at which the exponential factor for the smallest exponent is numerically zero.
"""
min_exp = np.finfo(np.core.numerictypes.float64).max
if self.s_exp is not None:
min_exp = min(min_exp, min(self.s_exp))
if self.p_exp is not None:
min_exp = min(min_exp, min(self.p_exp))
if self.d_exp is not None:
min_exp = min(min_exp, min(self.d_exp))
if self.f_exp is not None:
min_exp = min(min_exp, min(self.f_exp))
max_range = -log(np.finfo(np.core.numerictypes.float64).tiny)/min_exp
return max_range
def get_orbitals(self, q_numbers, exponents, coefficients, r):
"""
Evaluates Slater orbitals at position r.
IN:
q_numbers = array of principal quantum numbers
exponents = array of exponents
coefficients = linear combination coefficients for orbitals
r = space coordinates
OUT:
of = orbital values at each r dimension (orbitals, r)
dof = first orbital derivatives at r
ddof = second orbital derivatives
"""
# Begin by calculating the values of the Slater function, and its first 2 derivatives
# for the given exponents.
# gives (n, r) dimension [where n is number of functions]
f = self.G(q_numbers, exponents, r)
df = self.DG(q_numbers, exponents, r, f)
ddf = self.DDG(q_numbers, exponents, r, f)
# Now construct the orbitals by multiplication with orbital coefficients
of = np.einsum('ij,jk->ik', coefficients, f) # (i=orbital, j=function, k=r)
dof = np.einsum('ij,jk->ik', coefficients, df)
ddof = np.einsum('ij,jk->ik', coefficients, ddf)
return of, dof, ddof
def G(self, n, zeta, r):
"""Evaluates the radial Slater orbital R(r) = N r^{n-1} exp(-zeta r)
arguments:
n: principal quantum number
zeta: exponent
r: distance from nucleus
"""
# Principal quantum number factors for STO normalization
FACTORS = np.array([factorial(2*nn) for nn in range(1,max(n)+1)])**(-0.5)
n_facs = FACTORS[n - 1]
try:
c = n_facs*(2.0*zeta)**(n + 0.5)
except ValueError:
print("Exponents and principal number factors are different shapes.")
print("Did you typo a ',' for a decimal point? e.g. '1,23456' for '1.23456'")
raise ValueError("exponent or principal number error")
rn = np.power.outer(r, (n - 1))
es = np.einsum('j,ij->ji', c, rn)
pw = np.exp(-np.outer(zeta, r))
return es*pw
def DG(self, n, e, r, f):
"""Evaluates the first derivative R'(r) of the radial Slater orbital
R(r) = N r^{n-1} exp(-zeta r) as R'(r) = [(n-1)/r - zeta] R(r).
arguments:
n: principal quantum number
zeta: exponent
r: distance from nucleus
f: undifferentiated function
"""
pre = -e[:, None] + np.divide.outer((n - 1), r)
return pre*f
def DDG(self, n, e, r, f):
"""Evaluates the second derivative R''(r) of the radial Slater orbital
R(r) = N r^{n-1} exp(-zeta r) as
R''(r) = {[(n-1)/r - zeta]^2 - (n-1)/r^2} R(r)
arguments:
n: principal quantum number
zeta: exponent
r: distance from nucleus
f: undifferentiated function
"""
pre = (-e[:, None] + np.divide.outer((n - 1), r))**2
pre -= np.divide.outer((n - 1), r**2)
return pre*f
def get_nuclear_potential(self, r):
"""
Returns the -Z/r potential energy curve for the atom for all points r
"""
return -self.nuclear_charge/r
def get_gaussian_nuclear_potential(self, r, gamma=0.2):
"""
Returns a Gaussian approximation to the nuclear potential as defined
in
F. Brockherde, L. Vogt, L. Li, M. E. Tuckerman, K. Burke, and K. R. Müller, Nat. Commun. 8, (2017).
DOI: 10.1038/s41467-017-00839-3
v(r) = Z*exp(-r/(2*gamma**2))
where Z is the nuclear charge and gamma is a width parameter chosen as 0.2 in the reference.
"""
return -self.nuclear_charge*np.exp(-r**2/(2*gamma**2))
def get_color(self):
"""
Returns RGB color of element for plotting.
"""
return COLOR_DICT[self.element]
def libxc_eval(self, r, functional='gga_x_pbe', restricted=False, density_threshold=None, sigma_threshold=None, nan_check=False, ext_params=None):
'''Evaluates a functional with the atomic density data using libxc'''
d0, d1, g0, g1, t0, t1, l0, l1 = self.get_densities(r)
if not HAVE_LIBXC:
raise ImportError('Cannot evaluate functional since pylibxc could not be imported.')
func = pylibxc.LibXCFunctional(functional, "unpolarized" if restricted else "polarized")
# Did we get a threshold?
if density_threshold is not None:
func.set_dens_threshold(density_threshold)
if sigma_threshold is not None:
func.set_sigma_threshold(sigma_threshold)
# Did we get external parameters?
if ext_params is not None:
func.set_ext_params(ext_params)
# Create input
inp = {}
if restricted:
inp["rho"] = d0+d1
inp["sigma"] = np.multiply(g0+g1,g0+g1)
inp["lapl"]= l0+l1
inp["tau"]= t0+t1
else:
rho_array = np.zeros((d0.size,2), dtype='float64')
sigma_array = np.zeros((d0.size,3), dtype='float64')
lapl_array = np.zeros((d0.size,2), dtype='float64')
tau_array = np.zeros((d0.size,2), dtype='float64')
rho_array[:,0]=d0
rho_array[:,1]=d1
sigma_array[:,0]=np.multiply(g0, g0)
sigma_array[:,1]=np.multiply(g0, g1)
sigma_array[:,2]=np.multiply(g1, g1)
lapl_array[:,0]=l0
lapl_array[:,1]=l1
tau_array[:,0]=t0
tau_array[:,1]=t1
inp["rho"] = rho_array
inp["sigma"] = sigma_array
inp["lapl"]= lapl_array
inp["tau"]= tau_array
# Compute functional
ret = func.compute(inp)
# Get energy density per particle
zk = ret.get("zk", np.zeros_like(d0))
if zk.shape != d0.shape:
zk = np.reshape(zk, d0.shape)
# Energy density
nE = np.multiply(zk,d0+d1)
# First derivatives
vrho = ret["vrho"]
vsigma = ret.get("vsigma", np.zeros_like(inp["sigma"]))
vtau = ret.get("vtau", np.zeros_like(inp["tau"]))
vlapl = ret.get("vlapl", np.zeros_like(inp["lapl"]))
# Earlier versions of PyLibXC return the wrong shape, so reshape
# just to be sure
vrho = np.reshape(vrho, inp["rho"].shape)
vsigma = np.reshape(vsigma, inp["sigma"].shape)
vlapl = np.reshape(vlapl, inp["lapl"].shape)
vtau = np.reshape(vtau, inp["tau"].shape)
if nan_check:
# Indices of NaNs
nanidx = np.isnan(nE)
if len(nE[nanidx]):
print('nanidx len={} with {} nans'.format(nanidx.size, len(nE[nanidx])))
#print('NaN densities\n{}'.format((d0+d1)[nanidx]))
for i in range(len(nanidx)):
if nanidx[i]:
if restricted:
print('NaN at rho= {:e} sigmaa= {:e} lapl= {: e} tau={:e}'.format(rho_array[i],sigma_array[i], lapl_array[i], tau_array[i]))
else:
print('NaN at rhoa= {:e} rhob= {:e} sigmaaa= {:e} sigmaab= {: e} sigmabb= {:e} lapla= {: e} laplb= {: e} taua= {:e} taub={:e}'.format(d0[i],d1[i],sigma_array[i,0], sigma_array[i,1], sigma_array[i,2], l0[i], l1[i], t0[i], t1[i]))
return nE, vrho, vsigma, vtau, vlapl
def xcfun_eval(self, r, functional='PBEX', restricted=False):
'''Evaluates a functional with the atomic density data using XCFun'''
d0, d1, g0, g1, t0, t1, l0, l1 = self.get_densities(r)
if not HAVE_XCFUN:
raise ImportError('Cannot evaluate functional since xcfun could not be imported.')
# Initialize functional
func = xcfun.Functional(functional)
# Create input
if restricted:
density = d0+d1
densgrad = np.zeros((r.size, 3))
densgrad[:,0] = g0 + g1
res = func.eval_energy_n(density, densgrad)
else:
density = np.zeros((r.size, 2))
density[:,0] = d0
density[:,1] = d1
densgrad = np.zeros((r.size, 3, 2))
densgrad[:,0,0] = g0
densgrad[:,0,1] = g1
res = func.eval_energy_ab(density, densgrad)
return res
def static_init(cls):
if getattr(cls, "static_init", None):
cls.static_init()
return cls
@static_init
class AtomData:
"""
Class encapsulating raw data for all atoms.
Primarily used to populate Atom class on instantiation.
No need to give entries if the atom has none of that orbital, e.g. P and D for Li
In these cases the default None and [0, 0] occupation should be used.
Achieved by .get() accessing
"""
# Testing data
nuclear_charge = {
'H' : 1.0,
'HE' : 2.0,
'LI' : 3.0,
'BE' : 4.0,
'B' : 5.0,
'C' : 6.0,
'N' : 7.0,
'O' : 8.0,
'F' : 9.0,
'NE' : 10.0,
'NA' : 11.0,
'MG' : 12.0,
'AL' : 13.0,
'SI' : 14.0,
'P' : 15.0,
'S' : 16.0,
'CL' : 17.0,
'AR' : 18.0,
'K' : 19.0,
'CA' : 20.0,
'SC' : 21.0,
'TI' : 22.0,
'V' : 23.0,
'CR' : 24.0,
'MN' : 25.0,
'FE' : 26.0,
'CO' : 27.0,
'NI' : 28.0,
'CU' : 29.0,
'ZN' : 30.0,
'GA' : 31.0,
'GE' : 32.0,
'AS' : 33.0,
'SE' : 34.0,
'BR' : 35.0,
'KR' : 36.0,
'RB' : 37.0,
'SR' : 38.0,
'Y' : 39.0,
'ZR' : 40.0,
'NB' : 41.0,
'MO' : 42.0,
'TC' : 43.0,
'RU' : 44.0,
'RH' : 45.0,
'PD' : 46.0,
'AG' : 47.0,
'CD' : 48.0,
'IN' : 49.0,
'SN' : 50.0,
'SB' : 51.0,
'TE' : 52.0,
'I' : 53.0,
'XE' : 54.0,
'CS' : 55.0,
'BA' : 56.0,
'LA' : 57.0,
'CE' : 58.0,
'PR' : 59.0,
'ND' : 60.0,
'PM' : 61.0,
'SM' : 62.0,
'EU' : 63.0,
'GD' : 64.0,
'TB' : 65.0,
'DY' : 66.0,
'HO' : 67.0,
'ER' : 68.0,
'TM' : 69.0,
'YB' : 70.0,
'LU' : 71.0,
'HF' : 72.0,
'TA' : 73.0,
'W' : 74.0,
'RE' : 75.0,
'OS' : 76.0,
'IR' : 77.0,
'PT' : 78.0,
'AU' : 79.0,
'HG' : 80.0,
'TL' : 81.0,
'PB' : 82.0,
'BI' : 83.0,
'PO' : 84.0,
'AT' : 85.0,
'RN' : 86.0,
'FR' : 87.0,
'RA' : 88.0,
'AC' : 89.0,
'TH' : 90.0,
'PA' : 91.0,
'U' : 92.0,
'NP' : 93.0,
'PU' : 94.0,
'AM' : 95.0,
'CM' : 96.0,
'BK' : 97.0,
'CF' : 98.0,
'ES' : 99.0,
'FM' : 100.0,
'MD' : 101.0,
'NO' : 102.0,
'LR' : 103.0
}
electron_count = {}
ke_test = {}
s_exp = {}
s_coef = {}
s_n = {}
s_occ = {}
p_exp = {}
p_coef = {}
p_n = {}
p_occ = {}
d_exp = {}
d_coef = {}
d_n = {}
d_occ = {}
f_exp = {}
f_coef = {}
f_n = {}
f_occ = {}
@classmethod
def add_entry(self, a, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc):
'''Adds an entry to the database'''
self.ke_test[a] = Ekin
for ishell in range(len(ams)):
if ams[ishell] == 0:
self.s_exp[a] = xs[ishell]
self.s_coef[a] = np.transpose(cs[ishell])
self.s_n[a] = ns[ishell]
self.s_occ[a] = socc
elif ams[ishell] == 1:
self.p_exp[a] = xs[ishell]
self.p_coef[a] = np.transpose(cs[ishell])
self.p_n[a] = ns[ishell]
self.p_occ[a] = pocc
elif ams[ishell] == 2:
self.d_exp[a] = xs[ishell]
self.d_coef[a] = np.transpose(cs[ishell])
self.d_n[a] = ns[ishell]
self.d_occ[a] = docc
elif ams[ishell] == 3:
self.f_exp[a] = xs[ishell]
self.f_coef[a] = np.transpose(cs[ishell])
self.f_n[a] = ns[ishell]
self.f_occ[a] = focc
else:
raise ValueError('Angular momentum too large')
@classmethod
def static_init(self):
'''Initialize the data storage by reading in the tabulated wave functions'''
# Find where the current file is
curpath = os.path.dirname(os.path.abspath(__file__))
# Add neutral atoms
for a in self.nuclear_charge.keys():
# Light atoms
infile='{}/literature_data/k99l/neutral/{}'.format(curpath, a.lower())
if os.path.isfile(infile):
self.electron_count[a]=self.nuclear_charge[a]
Etot, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc = parse(infile)
self.add_entry(a, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc)
# Heavy atoms
infile='{}/literature_data/k00heavy/{}'.format(curpath, a.lower())
if os.path.isfile(infile):
self.electron_count[a]=self.nuclear_charge[a]
Etot, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc = parse(infile)
self.add_entry(a, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc)
# Add cations and anions
neutral_atoms = self.nuclear_charge.copy()
for a in neutral_atoms.keys():
infile='{}/literature_data/k99l/cation/{}.cat'.format(curpath, a.lower())
if os.path.isfile(infile):
cata="{}+".format(a)
self.nuclear_charge[cata]=neutral_atoms[a]
self.electron_count[cata]=neutral_atoms[a]-1
Etot, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc = parse(infile)
self.add_entry(cata, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc)
infile='{}/literature_data/k99l/anion/{}.an'.format(curpath, a.lower())
if os.path.isfile(infile):
ana="{}-".format(a)
self.nuclear_charge[ana]=neutral_atoms[a]
self.electron_count[ana]=neutral_atoms[a]+1
Etot, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc = parse(infile)
self.add_entry(ana, Ekin, ams, ns, xs, cs, socc, pocc, docc, focc)
class GridGenerator:
"""
Static class encapsulating grid generation methods.
Public facing method is make_grid(n)
Then spherical integration can be achieved using
np.sum(f(r)*4*pi*weight)
"""
@staticmethod
def make_grid(n, method='ahlrichsm3', R=1.0, quad='chebyshev2_mod'):
"""
Generate a radial integration grid containing n points.
n: desired grid points
method:
'ahlrichsm3' : parameter-free Ahlrichs M3 grid (default)
'ahlrichsm4' : Ahlrichs M4 grid with alpha=0.6
'handy' : Murray-Handy-Laming quadrature
'muraknowles' : Mura-Knowles quadrature
R: atomic size adjustment parameter, default R=1.0
quad: quadrature scheme, 'chebyshev2_mod' (default)
OUT:
n: number of grid points (note: may be different from that requested)
r: numpy array of coordinates
wt: numpy array of weights
"""
if method=='ahlrichsm3' or method == 'krack':
n, r, wt = GridGenerator.radial_ahlrichs(n, alpha=0.0, quad=quad)
elif method == 'ahlrichsm4':
n, r, wt = GridGenerator.radial_ahlrichs(n, alpha=0.6, quad=quad)
# The Laguerre quadrature doesn't appear to work for the wanted numbers of quadrature points
# elif method == 'laguerre':
# n, r, wt = GridGenerator.radial_laguerre(n)
elif method == 'handy':
n, r, wt = GridGenerator.radial_handy(n, quad=quad)
elif method == 'muraknowles':
n, r, wt = GridGenerator.radial_muraknowles(n, quad=quad)
elif method == 'becke':
n, r, wt = GridGenerator.radial_becke(n, quad=quad)
else:
raise ValueError('Unknown grid {}'.format(method))
# Eliminate any zero weights
r = r[wt!=0.0]
wt = wt[wt!=0.0]
n = len(r)
# Atomic scaling adjustment
r = r*R
wt = wt*R**3
return n, r, wt
@staticmethod
def chebyshev1(n):
"""Gauss-Chebyshev quadrature of the first kind for
calculating \int_{-1}^{1} f(x) dx = \sum_i w_i f(x_i).
Returns n, x, w.
See Abramowitz-Stegun p. 889
"""
# Index vector
ivec = np.asarray([ n+1-i for i in range(1,n+1)])
# Angles for sine and cosine
angles = (2*ivec-1)*pi/(2*n)
# Integration nodes
x = np.cos(angles)
# Integration weights
w = pi/n*np.sqrt(1-np.square(x))
return n, x, w
@staticmethod
def chebyshev2(n):
"""Gauss-Chebyshev quadrature of the second kind for
calculating \int_{-1}^{1} f(x) dx = \sum_i w_i f(x_i).
Returns n, x, w.
See Abramowitz-Stegun p. 889
"""
# Index vector
ivec = np.asarray([ n+1-i for i in range(1,n+1)])
# Angles for sine and cosine
angles = ivec*pi/(n+1)
# Integration nodes
x = np.cos(angles)
# Integration weights
w = pi/(n+1)*np.sqrt(1-np.square(x))
return n, x, w
@staticmethod
def chebyshev2_modified(n):
"""Modified Gauss-Chebyshev quadrature of the second kind for
calculating \int_{-1}^{1} f(x) dx = \sum_i w_i f(x_i).
Returns n, x, w.
See eqns (31)-(33) in J. M. Pérez‐Jordá, A. D. Becke, and
E. San‐Fabián, Automatic numerical integration techniques for
polyatomic molecules, J. Chem. Phys. 100, 6520 (1994);
doi:10.1063/1.467061
"""
# 1/(n+1)
oonpp=1.0/(n+1.0)
# Index vector
ivec = np.asarray([ i for i in range(1,n+1)])
# Angles for sine and cosine
angles = ivec*pi*oonpp
# Sines and cosines
sines = np.sin(angles)
cosines = np.cos(angles)
# Sine squared
sinesq = np.power(sines,2)
sinecos = np.multiply(sines,cosines)
# Integration weights
w = 16.0/3.0/(n+1.0) * np.power(sinesq,2)
# Integration nodes
x = 1.0 - 2.0*ivec*oonpp + 2/pi*np.multiply(1.0 + 2.0/3.0*sinesq, sinecos)
return n, x, w
@staticmethod
def chebyshev3(n):
"""Gauss-Chebyshev quadrature of the third kind for
calculating \int_{-1}^{1} f(x) dx = \sum_i w_i f(x_i).
Returns n, x, w.
See Abramowitz-Stegun p. 889
"""
# Index vector
ivec = np.asarray([ n+1-i for i in range(1,n+1)])
# Angles for sine and cosine
angles = 0.5*(2*ivec-1)*pi/(2*n+1)
# Original integration nodes and weights
x = np.square(np.cos(angles))
w = 2*pi/(2*n+1)*x*np.sqrt(np.divide(1-x,x))
# Transform from [0,1] to [-1,1]
x = 2*x-1
w = 2*w
return n, x, w
@staticmethod
def quadrature(n, quad='chebyshev2_mod'):
"""Quadrature rule for calculating \int_{-1}^{1} f(x) dx = \sum_i w_i f(x_i).
Returns n, x, w.
Input:
n: number of quadrature points
quad: quadrature rule
"""