/
test_OneParticleOperator.py
331 lines (280 loc) · 12.5 KB
/
test_OneParticleOperator.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
#!/usr/bin/env python3
## vi: tabstop=4 shiftwidth=4 softtabstop=4 expandtab
## ---------------------------------------------------------------------
##
## Copyright (C) 2018 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 unittest
import numpy as np
from numpy.testing import assert_array_almost_equal_nulp, assert_equal
from adcc import OneParticleOperator, zeros_like
from adcc.OneParticleOperator import product_trace
from adcc.testdata.cache import cache
from pytest import approx
class TestOneParticleOperator(unittest.TestCase):
def test_to_ndarray(self):
mp2diff = adcc.LazyMp(cache.refstate["h2o_sto3g"]).mp2_diffdm
dm_oo = mp2diff["o1o1"].to_ndarray()
dm_ov = mp2diff["o1v1"].to_ndarray()
dm_vv = mp2diff["v1v1"].to_ndarray()
dm_o = np.hstack((dm_oo, dm_ov))
dm_v = np.hstack((dm_ov.transpose(), dm_vv))
dm_full = np.vstack((dm_o, dm_v))
np.testing.assert_almost_equal(dm_full, mp2diff.to_ndarray(),
decimal=12)
def test_to_ndarray_nosym(self):
ref = cache.refstate["h2o_sto3g"]
dm = OneParticleOperator(ref.mospaces, is_symmetric=False)
dm["o1o1"].set_random()
dm["o1v1"].set_random()
dm["v1o1"].set_random()
dm["v1v1"].set_random()
dm_oo = dm["o1o1"].to_ndarray()
dm_ov = dm["o1v1"].to_ndarray()
dm_vo = dm["v1o1"].to_ndarray()
dm_vv = dm["v1v1"].to_ndarray()
dm_o = np.hstack((dm_oo, dm_ov))
dm_v = np.hstack((dm_vo, dm_vv))
dm_full = np.vstack((dm_o, dm_v))
assert_array_almost_equal_nulp(dm_full, dm.to_ndarray())
def test_product_trace_symmetric(self):
ref = cache.refstate["h2o_sto3g"]
dipx_mo = ref.operators.electric_dipole[0]
mp2diff_mo = adcc.LazyMp(ref).mp2_diffdm
mp2diff_ao = mp2diff_mo.to_ao_basis(ref)
mp2a = mp2diff_ao[0].to_ndarray()
mp2b = mp2diff_ao[1].to_ndarray()
dipx_ao = ref.operators.provider_ao.electric_dipole[0]
dipx_ref = np.sum(mp2a * dipx_ao) + np.sum(mp2b * dipx_ao)
oo = np.sum(
mp2diff_mo["o1o1"].to_ndarray() * dipx_mo["o1o1"].to_ndarray()
)
ov = 2.0 * np.sum(
mp2diff_mo["o1v1"].to_ndarray() * dipx_mo["o1v1"].to_ndarray()
)
vv = np.sum(
mp2diff_mo["v1v1"].to_ndarray() * dipx_mo["v1v1"].to_ndarray()
)
dipx_np = oo + ov + vv
assert dipx_np == approx(product_trace(mp2diff_mo, dipx_mo))
assert product_trace(mp2diff_mo, dipx_mo) == approx(dipx_ref)
assert product_trace(dipx_mo, mp2diff_mo) == approx(dipx_ref)
def test_product_trace_nonsymmetric(self):
ref = cache.refstate["cn_sto3g"]
dipx_mo = ref.operators.electric_dipole[0]
mp2diff_mo = adcc.LazyMp(ref).mp2_diffdm
mp2diff_nosym = OneParticleOperator(ref.mospaces, is_symmetric=False)
mp2diff_nosym.set_block("o1o1", mp2diff_mo["o1o1"])
mp2diff_nosym.set_block("o1v1", mp2diff_mo["o1v1"])
mp2diff_nosym.set_block("v1v1", mp2diff_mo["v1v1"])
mp2diff_nosym.set_block("v1o1",
zeros_like(mp2diff_mo["o1v1"].transpose()))
mp2diff_ao = mp2diff_nosym.to_ao_basis(ref)
mp2a = mp2diff_ao[0].to_ndarray()
mp2b = mp2diff_ao[1].to_ndarray()
dipx_ao = ref.operators.provider_ao.electric_dipole[0]
dipx_ref = np.sum(mp2a * dipx_ao) + np.sum(mp2b * dipx_ao)
oo = np.sum(
mp2diff_nosym["o1o1"].to_ndarray() * dipx_mo["o1o1"].to_ndarray()
)
ov = np.sum(
mp2diff_nosym["o1v1"].to_ndarray() * dipx_mo["o1v1"].to_ndarray()
)
vo = np.sum(
mp2diff_nosym["v1o1"].to_ndarray() * dipx_mo["o1v1"].to_ndarray().T
)
vv = np.sum(
mp2diff_nosym["v1v1"].to_ndarray() * dipx_mo["v1v1"].to_ndarray()
)
dipx_np = oo + ov + vo + vv
assert dipx_np == approx(product_trace(mp2diff_nosym, dipx_mo))
assert product_trace(mp2diff_nosym, dipx_mo) == approx(dipx_ref)
assert product_trace(dipx_mo, mp2diff_nosym) == approx(dipx_ref)
def test_product_trace_both_nonsymmetric(self):
ref = cache.refstate["cn_sto3g"]
dipx_mo = ref.operators.electric_dipole[0]
mp2diff_mo = adcc.LazyMp(ref).mp2_diffdm
mp2diff_nosym = OneParticleOperator(ref.mospaces, is_symmetric=False)
dipx_nosym = OneParticleOperator(ref.mospaces, is_symmetric=False)
mp2diff_nosym.set_block("o1o1", mp2diff_mo["o1o1"])
mp2diff_nosym.set_block("o1v1", mp2diff_mo["o1v1"])
mp2diff_nosym.set_block("v1v1", mp2diff_mo["v1v1"])
mp2diff_nosym.set_block("v1o1",
zeros_like(mp2diff_mo["o1v1"].transpose()))
mp2diff_ao = mp2diff_nosym.to_ao_basis(ref)
dipx_nosym.set_block("o1o1", dipx_mo["o1o1"])
dipx_nosym.set_block("o1v1", dipx_mo["o1v1"])
dipx_nosym.set_block("v1v1", dipx_mo["v1v1"])
dipx_nosym.set_block("v1o1", zeros_like(dipx_mo["o1v1"].transpose()))
dipx_ao = dipx_nosym.to_ao_basis(ref)
mp2a = mp2diff_ao[0].to_ndarray()
mp2b = mp2diff_ao[1].to_ndarray()
dipxa = dipx_ao[0].to_ndarray()
dipxb = dipx_ao[1].to_ndarray()
dipx_ref = np.sum(mp2a * dipxa) + np.sum(mp2b * dipxb)
oo = np.sum(
mp2diff_nosym["o1o1"].to_ndarray() * dipx_nosym["o1o1"].to_ndarray()
)
ov = np.sum(
mp2diff_nosym["o1v1"].to_ndarray() * dipx_nosym["o1v1"].to_ndarray()
)
vo = np.sum(
mp2diff_nosym["v1o1"].to_ndarray() * dipx_nosym["v1o1"].to_ndarray()
)
vv = np.sum(
mp2diff_nosym["v1v1"].to_ndarray() * dipx_nosym["v1v1"].to_ndarray()
)
dipx_np = oo + ov + vo + vv
assert dipx_np == approx(product_trace(mp2diff_nosym, dipx_nosym))
assert product_trace(mp2diff_nosym, dipx_nosym) == approx(dipx_ref)
assert product_trace(dipx_nosym, mp2diff_nosym) == approx(dipx_ref)
#
# Test operators
#
def test_copy(self):
ref = cache.refstate["h2o_sto3g"]
mp2diff = adcc.LazyMp(ref).mp2_diffdm
cpy = mp2diff.copy()
assert cpy.blocks == mp2diff.blocks
assert cpy.blocks_nonzero == mp2diff.blocks_nonzero
assert cpy.reference_state == mp2diff.reference_state
assert cpy.mospaces == mp2diff.mospaces
assert cpy.cartesian_transform == mp2diff.cartesian_transform
for b in mp2diff.blocks:
assert cpy.is_zero_block(b) == mp2diff.is_zero_block(b)
if not mp2diff.is_zero_block(b):
assert_equal(cpy.block(b).to_ndarray(),
mp2diff.block(b).to_ndarray())
assert cpy.block(b) is not mp2diff.block(b)
def test_add(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=False)
a["o1o1"].set_random()
a["v1o1"].set_random()
a["v1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=True)
b["o1o1"].set_random()
b["v1v1"].set_random()
assert_array_almost_equal_nulp((a + b).to_ndarray(),
a.to_ndarray() + b.to_ndarray())
assert_array_almost_equal_nulp((b + a).to_ndarray(),
a.to_ndarray() + b.to_ndarray())
def test_add_nosym(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=True)
a["o1o1"].set_random()
a["o1v1"].set_random()
a["v1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=False)
b["o1o1"].set_random()
b["o1v1"].set_random()
b["v1o1"].set_random()
b["v1v1"].set_random()
assert_array_almost_equal_nulp((a + b).to_ndarray(),
a.to_ndarray() + b.to_ndarray())
assert_array_almost_equal_nulp((b + a).to_ndarray(),
b.to_ndarray() + a.to_ndarray())
def test_iadd(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=False)
a["o1o1"].set_random()
a["v1o1"].set_random()
a["o1v1"].set_random()
a["v1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=False)
b["o1o1"].set_random()
b["o1v1"].set_random()
b["v1v1"].set_random()
ref = a.to_ndarray() + b.to_ndarray()
a += b
assert_array_almost_equal_nulp(a.to_ndarray(), ref)
def test_sub(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=True)
a["o1o1"].set_random()
a["o1v1"].set_random()
a["v1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=True)
b["o1o1"].set_random()
b["o1v1"].set_random()
b["v1v1"].set_random()
assert_array_almost_equal_nulp((a - b).to_ndarray(),
a.to_ndarray() - b.to_ndarray())
assert_array_almost_equal_nulp((b - a).to_ndarray(),
b.to_ndarray() - a.to_ndarray())
assert_array_almost_equal_nulp((a - b).to_ndarray(),
(a + (-1 * b)).to_ndarray())
assert_array_almost_equal_nulp((b - a).to_ndarray(),
(b + (-1 * a)).to_ndarray())
def test_sub_nosym(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=True)
a["o1o1"].set_random()
a["o1v1"].set_random()
a["v1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=False)
b["o1o1"].set_random()
b["o1v1"].set_random()
b["v1o1"].set_random()
b["v1v1"].set_random()
assert_array_almost_equal_nulp((a - b).to_ndarray(),
a.to_ndarray() - b.to_ndarray())
assert_array_almost_equal_nulp((b - a).to_ndarray(),
b.to_ndarray() - a.to_ndarray())
assert_array_almost_equal_nulp((a - b).to_ndarray(),
(a + (-1.0 * b)).to_ndarray())
assert_array_almost_equal_nulp((b - a).to_ndarray(),
(b + (-1.0 * a)).to_ndarray())
def test_isub(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=True)
a["o1o1"].set_random()
a["o1v1"].set_random()
b = OneParticleOperator(ref.mospaces, is_symmetric=True)
b["o1v1"].set_random()
b["v1v1"].set_random()
ref = a.to_ndarray() - b.to_ndarray()
a -= b
assert_array_almost_equal_nulp(a.to_ndarray(), ref)
def test_mul(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=True)
a["o1o1"].set_random()
a["o1v1"].set_random()
assert_array_almost_equal_nulp((1.2 * a).to_ndarray(),
1.2 * a.to_ndarray())
def test_rmul(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=False)
a["o1o1"].set_random()
a["v1o1"].set_random()
a["o1v1"].set_random()
assert_array_almost_equal_nulp((a * -1.8).to_ndarray(),
-1.8 * a.to_ndarray())
def test_imul(self):
ref = cache.refstate["h2o_sto3g"]
a = OneParticleOperator(ref.mospaces, is_symmetric=False)
a["o1o1"].set_random()
a["v1o1"].set_random()
a["o1v1"].set_random()
a["v1v1"].set_random()
ref = 12 * a.to_ndarray()
a *= 12
assert_array_almost_equal_nulp(a.to_ndarray(), ref)