-
Notifications
You must be signed in to change notification settings - Fork 19
/
testing.py
206 lines (179 loc) · 7.65 KB
/
testing.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
#!/usr/bin/env python3
## vi: tabstop=4 shiftwidth=4 softtabstop=4 expandtab
## ---------------------------------------------------------------------
##
## Copyright (C) 2019 by the adcc authors
##
## This file is part of adcc.
##
## adcc is free software: you can redistribute it and/or modify
## it under the terms of the GNU General Public License as published
## by the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## adcc is distributed in the hope that it will be useful,
## but WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with adcc. If not, see <http://www.gnu.org/licenses/>.
##
## ---------------------------------------------------------------------
import adcc
import numpy as np
from ..misc import assert_allclose_signfix
from collections import namedtuple
ii, jj, kk, ll = 0, 1, 2, 3
eri_chem_permutations = [(ii, jj, kk, ll), # (ij|kl)
(kk, ll, ii, jj), # (kl|ij)
(jj, ii, ll, kk), # (ji|lk)
(ll, kk, jj, ii), # (lk|ji)
(jj, ii, kk, ll), # (ji|kl)
(ll, kk, ii, jj), # (lk|ij)
(ii, jj, ll, kk), # (ij|lk)
(kk, ll, jj, ii)] # (kl|ji)
del ii, jj, kk, ll
# Spin symmetry helper
# provide allowed spin block transposition
# together with the prefactors that are needed to form the
# antisymmetrized integral from Chemists' notation ERIs.
# Example: Consider the <ab||ab> block in Physicists' notation
# <ab||ab> = <ab|ab> - <ab|ba> = (aa|bb) - (ab|ba)
# Here, the last term vanished, so this must be respected when
# computing the final integral. The first prefactor (pref1) in
# this case is 1, whereas the second prefactor (pref2) is 0 due to
# vanishing block of the antisymmetrized integral.
EriPermutationSpinAntiSymm = namedtuple('EriPermutationSpinAntiSymm',
['pref1', 'pref2', 'transposition'])
eri_phys_asymm_spin_allowed_prefactors = [
EriPermutationSpinAntiSymm(1, 1, "aaaa"),
EriPermutationSpinAntiSymm(1, 1, "bbbb"),
EriPermutationSpinAntiSymm(1, 0, "abab"),
EriPermutationSpinAntiSymm(1, 0, "baba"),
EriPermutationSpinAntiSymm(0, 1, "baab"),
EriPermutationSpinAntiSymm(0, 1, "abba"),
]
def eri_asymm_construction_test(scfres, core_orbitals=0):
hfdata = adcc.backends.import_scf_results(scfres)
refstate = adcc.ReferenceState(hfdata, core_orbitals=core_orbitals)
subspaces = refstate.mospaces.subspaces
ss_pairs = []
for i in range(len(subspaces)):
for j in range(i, len(subspaces)):
ss1, ss2 = subspaces[i], subspaces[j]
ss_pairs.append(ss1 + ss2)
# build the full ERI tensor
eri_chem = np.empty((hfdata.n_orbs, hfdata.n_orbs,
hfdata.n_orbs, hfdata.n_orbs))
sfull = slice(hfdata.n_orbs)
hfdata.fill_eri_ffff((sfull, sfull, sfull, sfull), eri_chem)
eri_phys = eri_chem.transpose(0, 2, 1, 3)
# full anti-symmetrized ERI tensor
eri_asymm = eri_phys - eri_phys.transpose(1, 0, 2, 3)
n_orbs = hfdata.n_orbs
n_alpha = hfdata.n_alpha
n_beta = hfdata.n_beta
n_orbs_alpha = hfdata.n_orbs_alpha
aro1 = slice(core_orbitals, n_alpha, 1)
bro1 = slice(
n_orbs_alpha + core_orbitals, n_orbs_alpha + n_beta, 1
)
aro2 = slice(0, core_orbitals, 1)
bro2 = slice(
n_orbs_alpha, n_orbs_alpha + core_orbitals, 1
)
arv = slice(n_alpha, n_orbs_alpha, 1)
brv = slice(n_orbs_alpha + n_beta, n_orbs, 1)
lookuptable = {
"o1": {"a": aro1, "b": bro1, },
"o2": {"a": aro2, "b": bro2, },
"v1": {"a": arv, "b": brv, }
}
n_elec = n_alpha + n_beta
n_virt_a = n_orbs_alpha - n_alpha
lookuptable_prelim = {
"o1": {
"a": slice(
0, n_alpha - core_orbitals, 1
), "b": slice(n_alpha - core_orbitals, n_elec, 1)
},
"o2": {
"a": slice(0, core_orbitals, 1),
"b": slice(core_orbitals, 2 * core_orbitals, 1)
},
"v1": {
"a": slice(0, n_virt_a, 1), "b": slice(n_virt_a, n_orbs - n_elec, 1)
}
}
# loop over all spaces and compare imported
# tensor to full tensor
for i in range(len(ss_pairs)):
for j in range(i, len(ss_pairs)):
p1, p2 = ss_pairs[i], ss_pairs[j]
s = p1 + p2
n = 2
s_clean = [s[i:i + n] for i in range(0, len(s), n)]
imported_asymm = refstate.eri(s).to_ndarray()
for allowed_spin in eri_phys_asymm_spin_allowed_prefactors:
sl = [lookuptable[x] for x in s_clean]
sl = [sl[x][y] for x, y in
enumerate(list(allowed_spin.transposition))]
sl2 = [lookuptable_prelim[x] for x in s_clean]
sl2 = [sl2[x][y] for x, y in
enumerate(list(allowed_spin.transposition))]
sl = tuple(sl)
sl2 = tuple(sl2)
np.testing.assert_almost_equal(
eri_asymm[sl], imported_asymm[sl2],
err_msg="""ERIs wrong in space {} """
"""and spin block {}""".format(s, allowed_spin)
)
def operator_import_test(scfres, ao_dict):
refstate = adcc.ReferenceState(scfres)
occa = refstate.orbital_coefficients_alpha("o1b").to_ndarray()
occb = refstate.orbital_coefficients_beta("o1b").to_ndarray()
virta = refstate.orbital_coefficients_alpha("v1b").to_ndarray()
virtb = refstate.orbital_coefficients_beta("v1b").to_ndarray()
for i, ao_component in enumerate(ao_dict):
dip_oo = np.einsum('ib,ba,ja->ij', occa, ao_component, occa)
dip_oo += np.einsum('ib,ba,ja->ij', occb, ao_component, occb)
dip_ov = np.einsum('ib,ba,ja->ij', occa, ao_component, virta)
dip_ov += np.einsum('ib,ba,ja->ij', occb, ao_component, virtb)
dip_vv = np.einsum('ib,ba,ja->ij', virta, ao_component, virta)
dip_vv += np.einsum('ib,ba,ja->ij', virtb, ao_component, virtb)
dip_mock = {"o1o1": dip_oo, "o1v1": dip_ov, "v1v1": dip_vv}
dip_imported = refstate.operators.electric_dipole[i]
for b in dip_imported.blocks:
assert_allclose_signfix(
dip_mock[b], dip_imported[b].to_ndarray(),
atol=refstate.conv_tol
)
def cached_backend_hf(backend, molecule, basis, multiplicity=1):
"""
Run the SCF for a backend and a particular test case (if not done)
and return the result.
"""
import adcc.backends
from adcc.testdata import geometry
global __cache_cached_backend_hf
def payload():
hfres = adcc.backends.run_hf(backend, xyz=geometry.xyz[molecule],
basis=basis, conv_tol=1e-13,
multiplicity=multiplicity,
conv_tol_grad=1e-11)
return adcc.backends.import_scf_results(hfres)
# For reasons not clear to me (mfh), caching does not work
# with pyscf
if backend == "pyscf":
return payload()
key = (backend, molecule, basis, str(multiplicity))
try:
return __cache_cached_backend_hf[key]
except NameError:
__cache_cached_backend_hf = {}
__cache_cached_backend_hf[key] = payload()
return __cache_cached_backend_hf[key]
except KeyError:
__cache_cached_backend_hf[key] = payload()
return __cache_cached_backend_hf[key]