Electronic repulsion integral containers.
ebcc.ham.eris.RERIs(mf, space, mo_coeff=None)
Bases: BaseERIs
, BaseRHamiltonian
Restricted ERIs container class.
Initialise the Hamiltonian.
Parameters: |
|
---|
Source code in ebcc/ham/base.py
def __init__(
self,
mf: SCF,
space: tuple[SpaceType, ...],
mo_coeff: Optional[tuple[CoeffType, ...]] = None,
) -> None:
"""Initialise the Hamiltonian.
Args:
mf: Mean-field object.
space: Space object for each index.
mo_coeff: Molecular orbital coefficients for each index.
"""
Namespace.__init__(self)
# Parameters:
self.__dict__["mf"] = mf
self.__dict__["space"] = space
self.__dict__["mo_coeff"] = mo_coeff if mo_coeff is not None else (mf.mo_coeff,) * 4
ebcc.ham.eris.RERIs.__getitem__(key)
Just-in-time getter.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/ham/eris.py
def __getitem__(self, key: str) -> NDArray[T]:
"""Just-in-time getter.
Args:
key: Key to get.
Returns:
ERIs for the given spaces.
"""
if key not in self._members.keys():
# Get the coefficients and shape
coeffs = self._get_coeffs(key)
coeffs = tuple(self._to_pyscf_backend(c) for c in coeffs)
shape = tuple(c.shape[-1] for c in coeffs)
# Transform the block
# TODO: Optimise for patially AO
if getattr(self.mf, "_eri", None) is not None:
block = pyscf.ao2mo.incore.general(self.mf._eri, coeffs, compact=False)
else:
block = pyscf.ao2mo.kernel(self.mf.mol, coeffs, compact=False)
block = np.reshape(block, shape)
# Store the block
self._members[key] = np.asarray(block, dtype=types[float])
return self._members[key]
ebcc.ham.eris.UERIs(mf, space, mo_coeff=None)
Bases: BaseERIs
, BaseUHamiltonian
Unrestricted ERIs container class.
Initialise the Hamiltonian.
Parameters: |
|
---|
Source code in ebcc/ham/base.py
def __init__(
self,
mf: SCF,
space: tuple[SpaceType, ...],
mo_coeff: Optional[tuple[CoeffType, ...]] = None,
) -> None:
"""Initialise the Hamiltonian.
Args:
mf: Mean-field object.
space: Space object for each index.
mo_coeff: Molecular orbital coefficients for each index.
"""
Namespace.__init__(self)
# Parameters:
self.__dict__["mf"] = mf
self.__dict__["space"] = space
self.__dict__["mo_coeff"] = mo_coeff if mo_coeff is not None else (mf.mo_coeff,) * 4
ebcc.ham.eris.UERIs.__getitem__(key)
Just-in-time getter.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/ham/eris.py
def __getitem__(self, key: str) -> RERIs:
"""Just-in-time getter.
Args:
key: Key to get.
Returns:
ERIs for the given spins.
"""
if key not in ("aaaa", "aabb", "bbaa", "bbbb"):
raise KeyError(f"Invalid key: {key}")
if key not in self._members:
i = "ab".index(key[0])
j = "ab".index(key[2])
self._members[key] = RERIs(
self.mf,
space=(self.space[0][i], self.space[1][i], self.space[2][j], self.space[3][j]),
mo_coeff=(
self.mo_coeff[0][i],
self.mo_coeff[1][i],
self.mo_coeff[2][j],
self.mo_coeff[3][j],
),
)
return self._members[key]
ebcc.ham.eris.GERIs(*args, **kwargs)
Bases: BaseERIs
, BaseGHamiltonian
Generalised ERIs container class.
Initialise the class.
Source code in ebcc/ham/eris.py
def __init__(self, *args: Any, **kwargs: Any) -> None:
"""Initialise the class."""
super().__init__(*args, **kwargs)
# Get the coefficients and shape
mo_a = [self._to_pyscf_backend(mo[: self.mf.mol.nao]) for mo in self.mo_coeff]
mo_b = [self._to_pyscf_backend(mo[self.mf.mol.nao :]) for mo in self.mo_coeff]
shape = tuple(mo.shape[-1] for mo in self.mo_coeff)
if len(set(shape)) != 1:
raise ValueError(
"MO coefficients must have the same number of basis functions for "
f"{self.__class__.__name__}."
)
nmo = shape[0]
if getattr(self.mf, "_eri", None) is not None:
array = pyscf.ao2mo.incore.general(self.mf._eri, mo_a)
array += pyscf.ao2mo.incore.general(self.mf._eri, mo_b)
array += pyscf.ao2mo.incore.general(self.mf._eri, mo_a[:2] + mo_b[2:])
array += pyscf.ao2mo.incore.general(self.mf._eri, mo_b[:2] + mo_a[2:])
else:
array = pyscf.ao2mo.kernel(self.mf.mol, mo_a)
array += pyscf.ao2mo.kernel(self.mf.mol, mo_b)
array += pyscf.ao2mo.kernel(self.mf.mol, mo_a[:2] + mo_b[2:])
array += pyscf.ao2mo.kernel(self.mf.mol, mo_b[:2] + mo_a[2:])
array = pyscf.ao2mo.addons.restore(1, array, nmo)
array = np.reshape(array, shape)
array = self._to_ebcc_backend(array)
# Transform to antisymmetric Physicist's notation
array = np.transpose(array, (0, 2, 1, 3)) - np.transpose(array, (0, 2, 3, 1))
# Store the array
self.__dict__["_array"] = array
ebcc.ham.eris.GERIs.__getitem__(key)
Just-in-time getter.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/ham/eris.py
def __getitem__(self, key: str) -> NDArray[T]:
"""Just-in-time getter.
Args:
key: Key to get.
Returns:
ERIs for the given spaces.
"""
if "p" in key:
raise NotImplementedError(f"AO basis not supported in {self.__class__.__name__}.")
return self._array[self._get_slices(key)]