def test_explicit_rhf_outside_solver_off_diagonal_blocks(): mol = molecules.molecule_water_sto3g() mol.build() mf = pyscf.scf.RHF(mol) mf.kernel() mocoeffs = mf.mo_coeff moenergies = mf.mo_energy ao2mo = AO2MOpyscf(mocoeffs, mol.verbose, mol) ao2mo.perform_rhf_full() tei_mo = ao2mo.tei_mo[0] C = mocoeffs E = np.diag(moenergies) occupations = utils.occupations_from_pyscf_mol(mol, mocoeffs) nocc, nvirt, _, _ = occupations A = eqns.form_rpa_a_matrix_mo_singlet_full(E, tei_mo, nocc) B = eqns.form_rpa_b_matrix_mo_singlet_full(tei_mo, nocc) G = np.block([[A, B], [B, A]]) assert G.shape == (2 * nocc * nvirt, 2 * nocc * nvirt) G_inv = np.linalg.inv(G) components = 3 integrals_dipole_ao = mol.intor("cint1e_r_sph", comp=components) integrals_dipole_mo_ai = [] for component in range(components): integrals_dipole_mo_ai_component = np.dot( C[:, nocc:].T, np.dot(integrals_dipole_ao[component, ...], C[:, :nocc])).reshape(-1, order="F") integrals_dipole_mo_ai.append(integrals_dipole_mo_ai_component) integrals_dipole_mo_ai = np.stack(integrals_dipole_mo_ai, axis=0).T integrals_dipole_mo_ai_super = np.concatenate( (integrals_dipole_mo_ai, -integrals_dipole_mo_ai), axis=0) rhsvecs = integrals_dipole_mo_ai_super rspvecs = np.dot(G_inv, rhsvecs) polarizability = 4 * np.dot(rhsvecs.T, rspvecs) / 2 # pylint: disable=bad-whitespace result__0_00 = np.array([[7.93556221, 0.0, 0.0], [0.0, 3.06821077, 0.0], [0.0, 0.0, 0.05038621]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(polarizability, result__0_00, rtol=rtol, atol=atol)
def test_integrals_pyscf(): mol = molecules.molecule_water_sto3g() mol.build() integral_generator = integrals.IntegralsPyscf(mol) np.testing.assert_equal(mol.intor("cint1e_r_sph", comp=3), integral_generator.integrals(integrals.DIPOLE)) np.testing.assert_equal( mol.intor("cint1e_cg_irxp_sph", comp=3), integral_generator.integrals(integrals.ANGMOM_COMMON_GAUGE), )
def test_explicit_uhf(): mol = molecules.molecule_water_sto3g() mol.charge = 1 mol.spin = 1 mol.build() mf = pyscf.scf.UHF(mol) mf.kernel() C = np.stack(mf.mo_coeff, axis=0) E_a = np.diag(mf.mo_energy[0]) E_b = np.diag(mf.mo_energy[1]) assert E_a.shape == E_b.shape E = np.stack((E_a, E_b), axis=0) integrals_dipole_ao = mol.intor("cint1e_r_sph", comp=3) occupations = utils.occupations_from_pyscf_mol(mol, C) solver = iterators.ExactInv(C, E, occupations) ao2mo = AO2MOpyscf(C, mol.verbose, mol) ao2mo.perform_uhf_full() solver.tei_mo = ao2mo.tei_mo solver.tei_mo_type = "full" driver = cphf.CPHF(solver) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = integrals_dipole_ao driver.add_operator(operator_dipole) driver.set_frequencies() driver.run(solver_type="exact", hamiltonian="rpa", spin="singlet") assert len(driver.frequencies) == len(driver.results) == 1 res = driver.results[0] print(res) atol = 1.0e-5 rtol = 0.0 np.testing.assert_allclose(res, ref_water_cation_UHF_HF_STO3G, rtol=rtol, atol=atol)
def test_explicit_uhf() -> None: mol = molecules.molecule_water_sto3g() mol.charge = 1 mol.spin = 1 mol.build() mf = pyscf.scf.UHF(mol) mf.kernel() C = np.stack(mf.mo_coeff, axis=0) C_a = C[0, ...] C_b = C[1, ...] E_a = np.diag(mf.mo_energy[0]) E_b = np.diag(mf.mo_energy[1]) assert C_a.shape == C_b.shape assert E_a.shape == E_b.shape E = np.stack((E_a, E_b), axis=0) norb = C_a.shape[1] nocc_a, nocc_b = mol.nelec nvirt_a, nvirt_b = norb - nocc_a, norb - nocc_b occupations = [nocc_a, nvirt_a, nocc_b, nvirt_b] C_occ_alph = C_a[:, :nocc_a] C_virt_alph = C_a[:, nocc_a:] C_occ_beta = C_b[:, :nocc_b] C_virt_beta = C_b[:, nocc_b:] C_ovov_aaaa = (C_occ_alph, C_virt_alph, C_occ_alph, C_virt_alph) C_ovov_aabb = (C_occ_alph, C_virt_alph, C_occ_beta, C_virt_beta) C_ovov_bbaa = (C_occ_beta, C_virt_beta, C_occ_alph, C_virt_alph) C_ovov_bbbb = (C_occ_beta, C_virt_beta, C_occ_beta, C_virt_beta) C_oovv_aaaa = (C_occ_alph, C_occ_alph, C_virt_alph, C_virt_alph) C_oovv_bbbb = (C_occ_beta, C_occ_beta, C_virt_beta, C_virt_beta) tei_mo_ovov_aaaa = pyscf.ao2mo.general(mol, C_ovov_aaaa, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nvirt_a, nocc_a, nvirt_a) tei_mo_ovov_aabb = pyscf.ao2mo.general(mol, C_ovov_aabb, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nvirt_a, nocc_b, nvirt_b) tei_mo_ovov_bbaa = pyscf.ao2mo.general(mol, C_ovov_bbaa, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nvirt_b, nocc_a, nvirt_a) tei_mo_ovov_bbbb = pyscf.ao2mo.general(mol, C_ovov_bbbb, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nvirt_b, nocc_b, nvirt_b) tei_mo_oovv_aaaa = pyscf.ao2mo.general(mol, C_oovv_aaaa, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nocc_a, nvirt_a, nvirt_a) tei_mo_oovv_bbbb = pyscf.ao2mo.general(mol, C_oovv_bbbb, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nocc_b, nvirt_b, nvirt_b) integrals_dipole_ao = mol.intor("cint1e_r_sph", comp=3) solver = solvers.ExactInv(C, E, occupations) solver.tei_mo = ( tei_mo_ovov_aaaa, tei_mo_ovov_aabb, tei_mo_ovov_bbaa, tei_mo_ovov_bbbb, tei_mo_oovv_aaaa, tei_mo_oovv_bbbb, ) solver.tei_mo_type = AO2MOTransformationType.partial driver = cphf.CPHF(solver) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = integrals_dipole_ao driver.add_operator(operator_dipole) driver.set_frequencies() driver.run( hamiltonian=Hamiltonian.RPA, spin=Spin.singlet, program=Program.PySCF, program_obj=mol, ) assert len(driver.frequencies) == len(driver.results) == 1 res = driver.results[0] print(res) atol = 1.0e-5 rtol = 0.0 np.testing.assert_allclose(res, ref_water_cation_UHF_HF_STO3G, rtol=rtol, atol=atol)
def test_explicit_uhf_outside_solver() -> None: mol = molecules.molecule_water_sto3g() mol.charge = 1 mol.spin = 1 mol.build() mf = pyscf.scf.UHF(mol) mf.kernel() C_a = mf.mo_coeff[0] C_b = mf.mo_coeff[1] E_a = np.diag(mf.mo_energy[0]) E_b = np.diag(mf.mo_energy[1]) assert C_a.shape == C_b.shape assert E_a.shape == E_b.shape norb = C_a.shape[1] nocc_a, nocc_b = mol.nelec nvirt_a, nvirt_b = norb - nocc_a, norb - nocc_b occupations = [nocc_a, nvirt_a, nocc_b, nvirt_b] C_occ_alph = C_a[:, :nocc_a] C_virt_alph = C_a[:, nocc_a:] C_occ_beta = C_b[:, :nocc_b] C_virt_beta = C_b[:, nocc_b:] C_ovov_aaaa = (C_occ_alph, C_virt_alph, C_occ_alph, C_virt_alph) C_ovov_aabb = (C_occ_alph, C_virt_alph, C_occ_beta, C_virt_beta) C_ovov_bbaa = (C_occ_beta, C_virt_beta, C_occ_alph, C_virt_alph) C_ovov_bbbb = (C_occ_beta, C_virt_beta, C_occ_beta, C_virt_beta) C_oovv_aaaa = (C_occ_alph, C_occ_alph, C_virt_alph, C_virt_alph) C_oovv_bbbb = (C_occ_beta, C_occ_beta, C_virt_beta, C_virt_beta) tei_mo_ovov_aaaa = pyscf.ao2mo.general(mol, C_ovov_aaaa, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nvirt_a, nocc_a, nvirt_a) tei_mo_ovov_aabb = pyscf.ao2mo.general(mol, C_ovov_aabb, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nvirt_a, nocc_b, nvirt_b) tei_mo_ovov_bbaa = pyscf.ao2mo.general(mol, C_ovov_bbaa, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nvirt_b, nocc_a, nvirt_a) tei_mo_ovov_bbbb = pyscf.ao2mo.general(mol, C_ovov_bbbb, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nvirt_b, nocc_b, nvirt_b) tei_mo_oovv_aaaa = pyscf.ao2mo.general(mol, C_oovv_aaaa, aosym="s4", compact=False, verbose=5).reshape( nocc_a, nocc_a, nvirt_a, nvirt_a) tei_mo_oovv_bbbb = pyscf.ao2mo.general(mol, C_oovv_bbbb, aosym="s4", compact=False, verbose=5).reshape( nocc_b, nocc_b, nvirt_b, nvirt_b) A_s_ss_a = eqns.form_rpa_a_matrix_mo_singlet_ss_partial( E_a, tei_mo_ovov_aaaa, tei_mo_oovv_aaaa) A_s_os_a = eqns.form_rpa_a_matrix_mo_singlet_os_partial(tei_mo_ovov_aabb) B_s_ss_a = eqns.form_rpa_b_matrix_mo_singlet_ss_partial(tei_mo_ovov_aaaa) B_s_os_a = eqns.form_rpa_b_matrix_mo_singlet_os_partial(tei_mo_ovov_aabb) A_s_ss_b = eqns.form_rpa_a_matrix_mo_singlet_ss_partial( E_b, tei_mo_ovov_bbbb, tei_mo_oovv_bbbb) A_s_os_b = eqns.form_rpa_a_matrix_mo_singlet_os_partial(tei_mo_ovov_bbaa) B_s_ss_b = eqns.form_rpa_b_matrix_mo_singlet_ss_partial(tei_mo_ovov_bbbb) B_s_os_b = eqns.form_rpa_b_matrix_mo_singlet_os_partial(tei_mo_ovov_bbaa) # Since the "triplet" part contains no Coulomb contribution, and # (xx|yy) is only in the Coulomb part, there is no ss/os # separation for the triplet part. G_aa = np.block([[A_s_ss_a, B_s_ss_a], [B_s_ss_a, A_s_ss_a]]) G_bb = np.block([[A_s_ss_b, B_s_ss_b], [B_s_ss_b, A_s_ss_b]]) G_ab = np.block([[A_s_os_a, B_s_os_a], [B_s_os_a, A_s_os_a]]) G_ba = np.block([[A_s_os_b, B_s_os_b], [B_s_os_b, A_s_os_b]]) G_aa_inv = np.linalg.inv(G_aa) G_bb_inv = np.linalg.inv(G_bb) nov_aa = nocc_a * nvirt_a nov_bb = nocc_b * nvirt_b assert G_aa_inv.shape == (2 * nov_aa, 2 * nov_aa) assert G_bb_inv.shape == (2 * nov_bb, 2 * nov_bb) # Form the operator-independent part of the response vectors. left_a = np.linalg.inv(G_aa - np.dot(G_ab, np.dot(G_bb_inv, G_ba))) left_b = np.linalg.inv(G_bb - np.dot(G_ba, np.dot(G_aa_inv, G_ab))) integrals_dipole_ao = mol.intor("cint1e_r_sph", comp=3) integrals_dipole_mo_ai_a = [] integrals_dipole_mo_ai_b = [] for comp in range(3): integrals_dipole_mo_ai_comp_a = np.dot( C_a[:, nocc_a:].T, np.dot(integrals_dipole_ao[comp, ...], C_a[:, :nocc_a])) integrals_dipole_mo_ai_comp_b = np.dot( C_b[:, nocc_b:].T, np.dot(integrals_dipole_ao[comp, ...], C_b[:, :nocc_b])) integrals_dipole_mo_ai_comp_a = np.reshape( integrals_dipole_mo_ai_comp_a, -1, order="F") integrals_dipole_mo_ai_comp_b = np.reshape( integrals_dipole_mo_ai_comp_b, -1, order="F") integrals_dipole_mo_ai_a.append(integrals_dipole_mo_ai_comp_a) integrals_dipole_mo_ai_b.append(integrals_dipole_mo_ai_comp_b) integrals_dipole_mo_ai_a = np.stack(integrals_dipole_mo_ai_a, axis=0).T integrals_dipole_mo_ai_b = np.stack(integrals_dipole_mo_ai_b, axis=0).T integrals_dipole_mo_ai_a_super = np.concatenate( (integrals_dipole_mo_ai_a, -integrals_dipole_mo_ai_a), axis=0) integrals_dipole_mo_ai_b_super = np.concatenate( (integrals_dipole_mo_ai_b, -integrals_dipole_mo_ai_b), axis=0) # Form the operator-dependent part of the response vectors. right_a = integrals_dipole_mo_ai_a_super - (np.dot( G_ab, np.dot(G_bb_inv, integrals_dipole_mo_ai_b_super))) right_b = integrals_dipole_mo_ai_b_super - (np.dot( G_ba, np.dot(G_aa_inv, integrals_dipole_mo_ai_a_super))) # The total response vector for each spin is the product of the # operator-independent (left) and operator-dependent (right) # parts. rspvec_a = np.dot(left_a, right_a) rspvec_b = np.dot(left_b, right_b) res_a = np.dot(integrals_dipole_mo_ai_a_super.T, rspvec_a) / 2 res_b = np.dot(integrals_dipole_mo_ai_b_super.T, rspvec_b) / 2 res_u = 2 * (res_a + res_b) print(res_u) atol = 1.0e-5 rtol = 0.0 np.testing.assert_allclose(res_u, ref_water_cation_UHF_HF_STO3G, rtol=rtol, atol=atol)
def test_explicit_uhf_from_rhf_outside_solver(): mol = molecules.molecule_water_sto3g() mol.build() mf = pyscf.scf.RHF(mol) mf.kernel() mocoeffs = mf.mo_coeff moenergies = mf.mo_energy ao2mo = AO2MOpyscf(mocoeffs, mol.verbose, mol) ao2mo.perform_rhf_full() tei_mo = ao2mo.tei_mo[0] C_a = mocoeffs C_b = C_a.copy() E_a = np.diag(moenergies) # E_b = E_a.copy() occupations = utils.occupations_from_pyscf_mol(mol, mocoeffs) nocc_a, nvirt_a, nocc_b, nvirt_b = occupations # Same-spin and opposite-spin contributions should add together # properly for restricted wavefunction. A_s = eqns.form_rpa_a_matrix_mo_singlet_full(E_a, tei_mo, nocc_a) A_s_ss = eqns.form_rpa_a_matrix_mo_singlet_ss_full(E_a, tei_mo, nocc_a) A_s_os = eqns.form_rpa_a_matrix_mo_singlet_os_full(tei_mo, nocc_a, nocc_b) np.testing.assert_allclose(A_s, A_s_ss + A_s_os) B_s = eqns.form_rpa_b_matrix_mo_singlet_full(tei_mo, nocc_a) B_s_ss = eqns.form_rpa_b_matrix_mo_singlet_ss_full(tei_mo, nocc_a) B_s_os = eqns.form_rpa_b_matrix_mo_singlet_os_full(tei_mo, nocc_a, nocc_b) np.testing.assert_allclose(B_s, B_s_ss + B_s_os) # Since the "triplet" part contains no Coulomb contribution, and # (xx|yy) is only in the Coulomb part, there is no ss/os # separation for the triplet part. G_r = np.block([[A_s, B_s], [B_s, A_s]]) G_aa = np.block([[A_s_ss, B_s_ss], [B_s_ss, A_s_ss]]) G_bb = G_aa.copy() G_ab = np.block([[A_s_os, B_s_os], [B_s_os, A_s_os]]) G_ba = G_ab.copy() np.testing.assert_allclose(G_r, (G_aa + G_ab)) G_r_inv = np.linalg.inv(G_r) G_aa_inv = np.linalg.inv(G_aa) G_bb_inv = np.linalg.inv(G_bb) assert G_r_inv.shape == (2 * nocc_a * nvirt_a, 2 * nocc_a * nvirt_a) assert G_aa_inv.shape == (2 * nocc_a * nvirt_a, 2 * nocc_a * nvirt_a) assert G_bb_inv.shape == (2 * nocc_b * nvirt_b, 2 * nocc_b * nvirt_b) # Form the operator-independent part of the response vectors. left_alph = np.linalg.inv(G_aa - np.dot(G_ab, np.dot(G_bb_inv, G_ba))) left_beta = np.linalg.inv(G_bb - np.dot(G_ba, np.dot(G_aa_inv, G_ab))) integrals_dipole_ao = mol.intor("cint1e_r_sph", comp=3) integrals_dipole_mo_ai_r = [] integrals_dipole_mo_ai_a = [] integrals_dipole_mo_ai_b = [] for comp in range(3): integrals_dipole_mo_ai_comp_r = np.dot( C_a[:, nocc_a:].T, np.dot(integrals_dipole_ao[comp, ...], C_a[:, :nocc_a]) ) integrals_dipole_mo_ai_comp_a = np.dot( C_a[:, nocc_a:].T, np.dot(integrals_dipole_ao[comp, ...], C_a[:, :nocc_a]) ) integrals_dipole_mo_ai_comp_b = np.dot( C_b[:, nocc_b:].T, np.dot(integrals_dipole_ao[comp, ...], C_b[:, :nocc_b]) ) integrals_dipole_mo_ai_comp_r = np.reshape(integrals_dipole_mo_ai_comp_r, -1, order="F") integrals_dipole_mo_ai_comp_a = np.reshape(integrals_dipole_mo_ai_comp_a, -1, order="F") integrals_dipole_mo_ai_comp_b = np.reshape(integrals_dipole_mo_ai_comp_b, -1, order="F") integrals_dipole_mo_ai_r.append(integrals_dipole_mo_ai_comp_r) integrals_dipole_mo_ai_a.append(integrals_dipole_mo_ai_comp_a) integrals_dipole_mo_ai_b.append(integrals_dipole_mo_ai_comp_b) integrals_dipole_mo_ai_r = np.stack(integrals_dipole_mo_ai_r, axis=0).T integrals_dipole_mo_ai_a = np.stack(integrals_dipole_mo_ai_a, axis=0).T integrals_dipole_mo_ai_b = np.stack(integrals_dipole_mo_ai_b, axis=0).T integrals_dipole_mo_ai_r_super = np.concatenate( (integrals_dipole_mo_ai_r, -integrals_dipole_mo_ai_r), axis=0 ) integrals_dipole_mo_ai_a_super = np.concatenate( (integrals_dipole_mo_ai_a, -integrals_dipole_mo_ai_a), axis=0 ) integrals_dipole_mo_ai_b_super = np.concatenate( (integrals_dipole_mo_ai_b, -integrals_dipole_mo_ai_b), axis=0 ) # Form the operator-dependent part of the response vectors. right_alph = integrals_dipole_mo_ai_a_super - ( np.dot(G_ab, np.dot(G_bb_inv, integrals_dipole_mo_ai_b_super)) ) right_beta = integrals_dipole_mo_ai_b_super - ( np.dot(G_ba, np.dot(G_aa_inv, integrals_dipole_mo_ai_a_super)) ) rspvec_r = np.dot(G_r_inv, integrals_dipole_mo_ai_r_super) # The total response vector for each spin is the product of the # operator-independent (left) and operator-dependent (right) # parts. rspvec_a = np.dot(left_alph, right_alph) rspvec_b = np.dot(left_beta, right_beta) res_r = 4 * np.dot(integrals_dipole_mo_ai_r_super.T, rspvec_r) / 2 res_a = np.dot(integrals_dipole_mo_ai_a_super.T, rspvec_a) / 2 res_b = np.dot(integrals_dipole_mo_ai_b_super.T, rspvec_b) / 2 res_u = 2 * (res_a + res_b) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(res_u, res_r, rtol=rtol, atol=atol) print(res_r) print(res_u)
def test_final_result_rhf_h2o_sto3g_tda_triplet() -> None: """Test correctness of the final result for water/STO-3G with the TDA approximation/CIS for triplet response induced by the dipole length operator computed with quantities from disk. """ hamiltonian = Hamiltonian.TDA spin = Spin.triplet C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = np.asarray([5, 2, 5, 2], dtype=int) stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = solvers.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = AO2MOTransformationType.full driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(hamiltonian=hamiltonian, spin=spin, program=None, program_obj=None) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[14.64430714, 0.0, 0.0], [0.0, 8.80921432, 0.0], [0.0, 0.0, 0.06859496]]) result__0_02 = np.array([[14.68168443, 0.0, 0.0], [0.0, 8.83562647, 0.0], [0.0, 0.0, 0.0689291]]) result__0_06 = np.array([[14.98774296, 0.0, 0.0], [0.0, 9.0532224, 0.0], [0.0, 0.0, 0.07172414]]) result__0_10 = np.array([[15.63997724, 0.0, 0.0], [0.0, 9.52504267, 0.0], [0.0, 0.0, 0.07805428]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability( Program.PySCF, mol, cphf.CPHF(solvers.ExactInv(C, E, occupations)), C, E, occupations, frequencies=frequencies, ) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol)
def test_final_result_rhf_h2o_sto3g_rpa_singlet() -> None: """Test correctness of the final result for water/STO-3G with full RPA for singlet response induced by the dipole length operator (the electric polarizability) computed with quantities from disk. """ hamiltonian = Hamiltonian.RPA spin = Spin.singlet C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = np.asarray([5, 2, 5, 2], dtype=int) stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = solvers.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = AO2MOTransformationType.full driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(hamiltonian=hamiltonian, spin=spin, program=None, program_obj=None) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[7.93556221, 0.0, 0.0], [0.0, 3.06821077, 0.0], [0.0, 0.0, 0.05038621]]) result__0_02 = np.array([[7.94312371, 0.0, 0.0], [0.0, 3.07051688, 0.0], [0.0, 0.0, 0.05054685]]) result__0_06 = np.array([[8.00414009, 0.0, 0.0], [0.0, 3.08913608, 0.0], [0.0, 0.0, 0.05186977]]) result__0_10 = np.array([[8.1290378, 0.0, 0.0], [0.0, 3.12731363, 0.0], [0.0, 0.0, 0.05473482]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) # Reminder: there's no call to do SCF here because we already have # the MO coefficients. mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability( Program.PySCF, mol, cphf.CPHF(solvers.ExactInv(C, E, occupations)), C, E, occupations, frequencies=frequencies, ) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol)
def test_final_result_rhf_h2o_sto3g_rpa_triplet() -> None: """Test correctness of the final result for water/STO-3G with full RPA for triplet response induced by the dipole length operator computed with quantities from disk. """ hamiltonian = Hamiltonian.RPA spin = Spin.triplet C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = np.asarray([5, 2, 5, 2], dtype=int) stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = solvers.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = AO2MOTransformationType.full driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(hamiltonian=hamiltonian, spin=spin, program=None, program_obj=None) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[26.59744305, 0.0, 0.0], [0.0, 18.11879557, 0.0], [0.0, 0.0, 0.07798969]]) result__0_02 = np.array([[26.68282287, 0.0, 0.0], [0.0, 18.19390051, 0.0], [0.0, 0.0, 0.07837521]]) result__0_06 = np.array([[27.38617401, 0.0, 0.0], [0.0, 18.81922578, 0.0], [0.0, 0.0, 0.08160226]]) result__0_10 = np.array([[28.91067234, 0.0, 0.0], [0.0, 20.21670386, 0.0], [0.0, 0.0, 0.08892512]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability( Program.PySCF, mol, cphf.CPHF(solvers.ExactInv(C, E, occupations)), C, E, occupations, frequencies=frequencies, ) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol)
def test_final_result_rhf_h2o_sto3g_tda_triplet(): hamiltonian = "tda" spin = "triplet" C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = [5, 2, 5, 2] stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = iterators.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = "full" driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(solver_type="exact", hamiltonian=hamiltonian, spin=spin) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[14.64430714, 0.0, 0.0], [0.0, 8.80921432, 0.0], [0.0, 0.0, 0.06859496]]) result__0_02 = np.array([[14.68168443, 0.0, 0.0], [0.0, 8.83562647, 0.0], [0.0, 0.0, 0.0689291]]) result__0_06 = np.array([[14.98774296, 0.0, 0.0], [0.0, 9.0532224, 0.0], [0.0, 0.0, 0.07172414]]) result__0_10 = np.array([[15.63997724, 0.0, 0.0], [0.0, 9.52504267, 0.0], [0.0, 0.0, 0.07805428]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability(Program.PySCF, mol, C, E, occupations, frequencies) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol) return
def test_final_result_rhf_h2o_sto3g_tda_singlet(): hamiltonian = "tda" spin = "singlet" C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = [5, 2, 5, 2] stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = iterators.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = "full" driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(solver_type="exact", hamiltonian=hamiltonian, spin=spin) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[8.89855952, 0.0, 0.0], [0.0, 4.00026556, 0.0], [0.0, 0.0, 0.0552774]]) result__0_02 = np.array([[8.90690928, 0.0, 0.0], [0.0, 4.00298342, 0.0], [0.0, 0.0, 0.05545196]]) result__0_06 = np.array([[8.97427725, 0.0, 0.0], [0.0, 4.02491517, 0.0], [0.0, 0.0, 0.05688918]]) result__0_10 = np.array([[9.11212633, 0.0, 0.0], [0.0, 4.06981937, 0.0], [0.0, 0.0, 0.05999934]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability(Program.PySCF, mol, C, E, occupations, frequencies) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol) return
def test_final_result_rhf_h2o_sto3g_rpa_singlet(): hamiltonian = "rpa" spin = "singlet" C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = [5, 2, 5, 2] stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = iterators.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = "full" driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(solver_type="exact", hamiltonian=hamiltonian, spin=spin) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[7.93556221, 0.0, 0.0], [0.0, 3.06821077, 0.0], [0.0, 0.0, 0.05038621]]) result__0_02 = np.array([[7.94312371, 0.0, 0.0], [0.0, 3.07051688, 0.0], [0.0, 0.0, 0.05054685]]) result__0_06 = np.array([[8.00414009, 0.0, 0.0], [0.0, 3.08913608, 0.0], [0.0, 0.0, 0.05186977]]) result__0_10 = np.array([[8.1290378, 0.0, 0.0], [0.0, 3.12731363, 0.0], [0.0, 0.0, 0.05473482]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) # Reminder: there's no call to do SCF here because we already have # the MO coefficients. mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability(Program.PySCF, mol, C, E, occupations, frequencies) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol) return
def test_final_result_rhf_h2o_sto3g_rpa_triplet(): hamiltonian = "rpa" spin = "triplet" C = utils.fix_mocoeffs_shape(utils.np_load(REFDIR / "C.npz")) E = utils.fix_moenergies_shape(utils.np_load(REFDIR / "F_MO.npz")) TEI_MO = utils.np_load(REFDIR / "TEI_MO.npz") # nocc_alph, nvirt_alph, nocc_beta, nvirt_beta occupations = [5, 2, 5, 2] stub = "h2o_sto3g_" dim = occupations[0] + occupations[1] mat_dipole_x = utils.parse_int_file_2(REFDIR / f"{stub}mux.dat", dim) mat_dipole_y = utils.parse_int_file_2(REFDIR / f"{stub}muy.dat", dim) mat_dipole_z = utils.parse_int_file_2(REFDIR / f"{stub}muz.dat", dim) solver = iterators.ExactInv(C, E, occupations) solver.tei_mo = (TEI_MO, ) solver.tei_mo_type = "full" driver = cphf.CPHF(solver) ao_integrals_dipole = np.stack((mat_dipole_x, mat_dipole_y, mat_dipole_z), axis=0) operator_dipole = operators.Operator(label="dipole", is_imaginary=False, is_spin_dependent=False) operator_dipole.ao_integrals = ao_integrals_dipole driver.add_operator(operator_dipole) frequencies = (0.0, 0.02, 0.06, 0.1) driver.set_frequencies(frequencies) driver.run(solver_type="exact", hamiltonian=hamiltonian, spin=spin) assert len(driver.results) == len(frequencies) result__0_00 = np.array([[26.59744305, 0.0, 0.0], [0.0, 18.11879557, 0.0], [0.0, 0.0, 0.07798969]]) result__0_02 = np.array([[26.68282287, 0.0, 0.0], [0.0, 18.19390051, 0.0], [0.0, 0.0, 0.07837521]]) result__0_06 = np.array([[27.38617401, 0.0, 0.0], [0.0, 18.81922578, 0.0], [0.0, 0.0, 0.08160226]]) result__0_10 = np.array([[28.91067234, 0.0, 0.0], [0.0, 20.21670386, 0.0], [0.0, 0.0, 0.08892512]]) atol = 1.0e-8 rtol = 0.0 np.testing.assert_allclose(driver.results[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(driver.results[3], result__0_10, rtol=rtol, atol=atol) mol = molecules.molecule_water_sto3g() mol.build() polarizability = electric.Polarizability(Program.PySCF, mol, C, E, occupations, frequencies) polarizability.form_operators() polarizability.run(hamiltonian=hamiltonian, spin=spin) polarizability.form_results() np.testing.assert_allclose(polarizability.polarizabilities[0], result__0_00, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[1], result__0_02, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[2], result__0_06, rtol=rtol, atol=atol) np.testing.assert_allclose(polarizability.polarizabilities[3], result__0_10, rtol=rtol, atol=atol) return
def test_jk_pyscf(): mol = molecules.molecule_water_sto3g() mol.build() jk_generator = integrals.JKPyscf(mol)