Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Some changes to facilitate newly implemented CCSDt #148

Merged
merged 9 commits into from
Jan 19, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .github/workflows/ci.yaml
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ jobs:
python -m pip install wheel --user
python -m pip install setuptools --upgrade --user
python -m pip install https://github.com/BoothGroup/dyson/archive/master.zip
python -m pip install -e .[dmet,ebcc] --user
python -m pip install .[dmet,ebcc] --user
- name: Run unit tests
run: |
python -m pip install pytest pytest-cov --user
Expand Down
28 changes: 18 additions & 10 deletions examples/ewf/molecules/26-ebcc-solvers.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@

mol = pyscf.gto.Mole()
mol.atom = ["N 0 0 0", "N 0 0 2"]
mol.basis = "aug-cc-pvdz"
mol.basis = "cc-pvdz"
mol.output = "pyscf.out"
mol.build()

Expand All @@ -24,19 +24,15 @@
casscf = pyscf.mcscf.CASSCF(mf, 8, 10)
casscf.kernel()

# Both these alternative specifications will always use an ebcc solver.
# emb = vayesta.ewf.EWF(mf, solver=f'EB{ansatz}', solver_options=dict(solve_lambda=False))
# emb = vayesta.ewf.EWF(mf, solver='ebcc', solver_options=dict(solve_lambda=False, ansatz=ansatz))

def get_emb_result(ansatz, bathtype="full"):
# Uses fastest available solver for given ansatz; PySCF if available, otherwise ebcc.
emb = vayesta.ewf.EWF(
mf, solver=ansatz, bath_options=dict(bathtype=bathtype), solver_options=dict(solve_lambda=False)
)
# Both these alternative specifications will always use an ebcc solver.
# Note that the capitalization of the solver name other than the ansatz is arbitrary.
# emb = vayesta.ewf.EWF(mf, solver=f'EB{ansatz}', bath_options=dict(bathtype=bathtype),
# solver_options=dict(solve_lambda=False))
# emb = vayesta.ewf.EWF(mf, solver='ebcc', bath_options=dict(bathtype=bathtype),
# solver_options=dict(solve_lambda=False, ansatz=ansatz))

with emb.iao_fragmentation() as f:
with f.rotational_symmetry(2, "y", center=(0, 0, 1)):
f.add_atomic_fragment(0)
Expand All @@ -45,12 +41,24 @@ def get_emb_result(ansatz, bathtype="full"):


e_ccsd = get_emb_result("CCSD", "full")
e_ccsdt = get_emb_result("CCSDT", "dmet")
e_CCSDT = get_emb_result("CCSDT", "dmet")
e_ccsdtprime = get_emb_result("CCSDt'", "full")


# CCSDt will set DMET orbitals as active space if no CAS is specified. To use a CAS, we must specify it explicitly.

# embfull = vayesta.ewf.EWF(mf, solver='ebcc', bath_options=dict(bathtype='full'), solver_options=dict(solve_lambda=False))
# for ffull, fact in zip( embfull.loop(), emb.loop() ):
# f.set_cas(c_occ=ffull.cluster.c_active_occ, c_vir=ffull.cluster.c_active_vir)
# emb.kernel()
e_CCSDt_ = get_emb_result("CCSDt", "dmet")
e_CCSDt = get_emb_result("CCSDt", "full")

print("E(HF)= %+16.8f Ha" % mf.e_tot)
print("E(CASCI)= %+16.8f Ha" % casci.e_tot)
print("E(CASSCF)= %+16.8f Ha" % casscf.e_tot)
print("E(CCSD, complete)= %+16.8f Ha" % e_ccsd)
print("E(emb. CCSDT, DMET CAS)= %+16.8f Ha" % e_ccsdt)
print("E(emb. CCSDT, DMET CAS)= %+16.8f Ha" % e_CCSDT)
print("E(emb. CCSDt, DMET CAS)= %+16.8f Ha" % e_CCSDt_)
print("E(emb. CCSDt', complete+DMET active space)= %+16.8f Ha" % e_ccsdtprime)
print("E(emb. CCSDt, complete+DMET active space)= %+16.8f Ha" % e_CCSDt)
2 changes: 1 addition & 1 deletion vayesta/core/qemb/qemb.py
Original file line number Diff line number Diff line change
Expand Up @@ -1574,7 +1574,7 @@ def _check_fragmentation(self, complete_occupied=True, complete_virtual=True, to
for s in range(nspin):
nmo_s = tspin(self.nmo, s)
nelec_s = tspin(nelec, s)
fragments = self.get_fragments(active=True, flags=dict(is_secfrag=False))
fragments = self.get_fragments(contributes=True, flags=dict(is_secfrag=False))
if not fragments:
return False
c_frags = np.hstack([tspin(x.c_frag, s) for x in fragments])
Expand Down
3 changes: 2 additions & 1 deletion vayesta/ewf/fragment.py
Original file line number Diff line number Diff line change
Expand Up @@ -306,7 +306,8 @@ def get_solver_options(self, solver):
has_actspace = (
(solver == "TCCSD")
or ("CCSDt'" in solver)
or (solver.upper() == "EBCC" and self.opts.solver_options["ansatz"] == "CCSDt'")
or ("CCSDt" in solver)
or (solver.upper() == "EBCC" and self.opts.solver_options["ansatz"] in ["CCSDt", "CCSDt'"])
)
if has_actspace:
# Set CAS orbitals
Expand Down
11 changes: 9 additions & 2 deletions vayesta/solver/ebcc.py
Original file line number Diff line number Diff line change
Expand Up @@ -61,7 +61,7 @@ def add_nonull_opt(d, key, newkey):
d[newkey] = val

opts = {}
for key, newkey in zip(["max_cycle", "conv_tol", "conv_tol_normt"], ["max_cycle", "e_tol", "t_tol"]):
for key, newkey in zip(["max_cycle", "conv_tol", "conv_tol_normt"], ["max_iter", "e_tol", "t_tol"]):
add_nonull_opt(opts, key, newkey)
return opts

Expand Down Expand Up @@ -225,7 +225,14 @@ class Options(UEBCC_Solver.Options, EB_REBCC_Solver.Options):

def get_couplings(self):
# EBCC wants contribution g_{xpq} p^\\dagger q b; need to transpose to get this contribution.
return tuple([x.transpose(0, 2, 1) for x in self.hamil.couplings])
#return tuple([x.transpose(0, 2, 1) for x in self.hamil.couplings])
# EBCC now wants an array, with spin as the first index
assert(np.allclose(self.hamil.couplings[0].shape, self.hamil.couplings[1].shape))
sh = self.hamil.couplings[0].shape
g = np.zeros((2, sh[0], sh[2], sh[1]), dtype=self.hamil.couplings[0].dtype)
g[0,:,:,:] = self.hamil.couplings[0].transpose(0,2,1)
g[1,:,:,:] = self.hamil.couplings[1].transpose(0,2,1)
return g

def construct_wavefunction(self, mycc, mo, mbos=None):
self.wf = EBCC_WaveFunction(
Expand Down
5 changes: 3 additions & 2 deletions vayesta/solver/hamiltonian.py
Original file line number Diff line number Diff line change
Expand Up @@ -360,8 +360,9 @@ def pad(a, diag_val):

clusmf.get_hcore = lambda *args, **kwargs: heff
if overwrite_fock:
clusmf.get_fock = lambda *args, **kwargs: pad_to_match(
clusmf.get_fock = lambda *args, **kwargs: np.asarray(pad_to_match(
self.get_fock(with_vext=True, use_seris=not force_bare_eris), dummy_energy
)
)
clusmf.get_veff = lambda *args, **kwargs: np.array(clusmf.get_fock(*args, **kwargs)) - np.array(
clusmf.get_hcore()
Expand Down Expand Up @@ -639,7 +640,7 @@ def get_fock(self, with_vext=True, use_seris=True, with_exxdiv=False):
+ einsum("iipq->pq", gab[oa, oa]) # Coulomb
- einsum("ipqi->pq", gbb[ob, :, :, ob])
) # Exchange
fock = ((fock[0] + dfa), (fock[1] + dfb))
fock = np.asarray(((fock[0] + dfa), (fock[1] + dfb)))
Copy link
Contributor

@maxnus maxnus Nov 8, 2023

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Can this can be packed into a single array, if fock[0] and fock[1] could have different shapes?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It can not be if the shapes are different. Wouldn't Fock matrices for alpha and beta always have the same shape?


return fock

Expand Down
3 changes: 3 additions & 0 deletions vayesta/tests/solver/test_ebcc.py
Original file line number Diff line number Diff line change
Expand Up @@ -103,6 +103,9 @@ def test_rccsdtprime_water_sto3g_dmet(self):
def test_uccsdtprime_water_sto3g_dmet(self):
return self._test("water_sto3g", "uhf", "CCSDt'", "CCSDT", bathtype="dmet", setcas=False)

def test_uCCSDt_water_sto3g_dmet(self):
return self._test("water_sto3g", "uhf", "CCSDt", "CCSDT", bathtype="dmet", setcas=False)

@pytest.mark.slow
def test_rccsdtprime_h4_sto3g_setcas_full(self):
return self._test("h4_sto3g", "rhf", "CCSDt'", "CCSDT", bathtype="full", setcas=True)
Expand Down