Fock matrix containers.
ebcc.ham.fock.RFock(mf, space, mo_coeff=None, g=None, shift=None, xi=None)
Bases: BaseFock
, BaseRHamiltonian
Restricted Fock matrix 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,
g: Optional[Namespace[Any]] = None,
shift: Optional[bool] = None,
xi: Optional[NDArray[float64]] = None,
) -> None:
"""Initialise the Hamiltonian.
Args:
mf: Mean-field object.
space: Space object for each index.
mo_coeff: Molecular orbital coefficients for each index.
g: Namespace containing blocks of the electron-boson coupling matrix.
shift: Shift the boson operators such that the Hamiltonian is normal-ordered with
respect to a coherent state. This removes the bosonic coupling to the static
mean-field density, introducing a constant energy shift.
xi: Shift in the bosonic operators to diagonalise the photon Hamiltonian.
"""
super().__init__(mf, space, mo_coeff=mo_coeff)
# Boson parameters:
self.__dict__["g"] = g
self.__dict__["shift"] = shift
self.__dict__["xi"] = xi
ebcc.ham.fock.RFock.__getitem__(key)
Just-in-time getter.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/ham/fock.py
def __getitem__(self, key: str) -> NDArray[T]:
"""Just-in-time getter.
Args:
key: Key to get.
Returns:
Fock matrix for the given spaces.
"""
if key not in self._members:
# Get the coefficients
coeffs = self._get_coeffs(key)
# Transform the block
fock_ao: NDArray[T] = np.asarray(self.mf.get_fock(), dtype=types[float])
if self._spin_index is not None:
fock_ao = fock_ao[self._spin_index]
block = util.einsum("pq,pi,qj->ij", fock_ao, *coeffs)
# Store the block
self._members[key] = block
if self.shift:
# Shift for bosons
xi = self.xi
g = np.copy(self.g.__getattr__(f"b{key}"))
g += np.transpose(self.g.__getattr__(f"b{key[::-1]}"), (0, 2, 1))
self._members[key] -= util.einsum("I,Ipq->pq", xi, g)
return self._members[key]
ebcc.ham.fock.UFock(mf, space, mo_coeff=None, g=None, shift=None, xi=None)
Bases: BaseFock
, BaseUHamiltonian
Unrestricted Fock matrix 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,
g: Optional[Namespace[Any]] = None,
shift: Optional[bool] = None,
xi: Optional[NDArray[float64]] = None,
) -> None:
"""Initialise the Hamiltonian.
Args:
mf: Mean-field object.
space: Space object for each index.
mo_coeff: Molecular orbital coefficients for each index.
g: Namespace containing blocks of the electron-boson coupling matrix.
shift: Shift the boson operators such that the Hamiltonian is normal-ordered with
respect to a coherent state. This removes the bosonic coupling to the static
mean-field density, introducing a constant energy shift.
xi: Shift in the bosonic operators to diagonalise the photon Hamiltonian.
"""
super().__init__(mf, space, mo_coeff=mo_coeff)
# Boson parameters:
self.__dict__["g"] = g
self.__dict__["shift"] = shift
self.__dict__["xi"] = xi
ebcc.ham.fock.UFock.__getitem__(key)
Just-in-time getter.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/ham/fock.py
def __getitem__(self, key: str) -> RFock:
"""Just-in-time getter.
Args:
key: Key to get.
Returns:
Fock matrix for the given spin.
"""
if key not in ("aa", "bb"):
raise KeyError(f"Invalid key: {key}")
if key not in self._members:
i = "ab".index(key[0])
self._members[key] = RFock(
self.mf,
space=(self.space[0][i], self.space[1][i]),
mo_coeff=(self.mo_coeff[0][i], self.mo_coeff[1][i]),
g=self.g[key] if self.g is not None else None,
shift=self.shift,
xi=self.xi,
)
self._members[key].__dict__["_spin_index"] = i
return self._members[key]
ebcc.ham.fock.GFock(mf, space, mo_coeff=None, g=None, shift=None, xi=None)
Bases: RFock
Generalised Fock matrix 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,
g: Optional[Namespace[Any]] = None,
shift: Optional[bool] = None,
xi: Optional[NDArray[float64]] = None,
) -> None:
"""Initialise the Hamiltonian.
Args:
mf: Mean-field object.
space: Space object for each index.
mo_coeff: Molecular orbital coefficients for each index.
g: Namespace containing blocks of the electron-boson coupling matrix.
shift: Shift the boson operators such that the Hamiltonian is normal-ordered with
respect to a coherent state. This removes the bosonic coupling to the static
mean-field density, introducing a constant energy shift.
xi: Shift in the bosonic operators to diagonalise the photon Hamiltonian.
"""
super().__init__(mf, space, mo_coeff=mo_coeff)
# Boson parameters:
self.__dict__["g"] = g
self.__dict__["shift"] = shift
self.__dict__["xi"] = xi