Generalised Brueckner-orbital coupled cluster.
ebcc.opt.gbrueckner.BruecknerGEBCC(cc, options=None, **kwargs)
Bases: BaseBruecknerEBCC
Generalised Brueckner-orbital coupled cluster.
Initialise the Brueckner EBCC object.
Parameters: |
|
---|
Source code in ebcc/opt/base.py
def __init__(
self,
cc: BaseEBCC,
options: Optional[BaseOptions] = None,
**kwargs: Any,
) -> None:
r"""Initialise the Brueckner EBCC object.
Args:
cc: Parent `EBCC` object.
options: Options for the EOM calculation.
**kwargs: Additional keyword arguments used to update `options`.
"""
# Options:
if options is None:
options = self.Options()
self.options = options
for key, val in kwargs.items():
setattr(self.options, key, val)
# Parameters:
self.cc = cc
self.mf = cc.mf
self.space = cc.space
self.log = cc.log
# Attributes:
self.converged = False
# Logging:
init_logging(cc.log)
cc.log.info(f"\n{ANSI.B}{ANSI.U}{self.name}{ANSI.R}")
cc.log.debug(f"{ANSI.B}{'*' * len(self.name)}{ANSI.R}")
cc.log.debug("")
cc.log.info(f"{ANSI.B}Options{ANSI.R}:")
cc.log.info(f" > e_tol: {ANSI.y}{self.options.e_tol}{ANSI.R}")
cc.log.info(f" > t_tol: {ANSI.y}{self.options.t_tol}{ANSI.R}")
cc.log.info(f" > max_iter: {ANSI.y}{self.options.max_iter}{ANSI.R}")
cc.log.info(f" > diis_space: {ANSI.y}{self.options.diis_space}{ANSI.R}")
cc.log.info(f" > diis_min_space: {ANSI.y}{self.options.diis_min_space}{ANSI.R}")
cc.log.info(f" > damping: {ANSI.y}{self.options.damping}{ANSI.R}")
cc.log.debug("")
ebcc.opt.gbrueckner.BruecknerGEBCC.get_rotation_matrix(u_tot=None, damping=None, t1=None)
Update the rotation matrix.
Also returns the total rotation matrix.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def get_rotation_matrix(
self,
u_tot: Optional[SpinArrayType] = None,
damping: Optional[BaseDamping] = None,
t1: Optional[SpinArrayType] = None,
) -> tuple[SpinArrayType, SpinArrayType]:
"""Update the rotation matrix.
Also returns the total rotation matrix.
Args:
u_tot: Total rotation matrix.
damping: Damping object.
t1: T1 amplitude.
Returns:
Rotation matrix and total rotation matrix.
"""
if t1 is None:
t1 = self.cc.t1
if u_tot is None:
u_tot = np.eye(self.cc.space.ncorr, dtype=types[float])
zocc = np.zeros((self.cc.space.ncocc, self.cc.space.ncocc))
zvir = np.zeros((self.cc.space.ncvir, self.cc.space.ncvir))
t1_block: NDArray[T] = np.block([[zocc, -t1], [np.transpose(t1), zvir]])
u = scipy.linalg.expm(t1_block)
u_tot = u_tot @ u
if np.linalg.det(u_tot) < 0:
u_tot = _put(u_tot, np.ix_(np.arange(u_tot.shape[0]), np.array([0])), -u_tot[:, 0])
a: NDArray[T] = np.asarray(np.real(scipy.linalg.logm(u_tot)), dtype=types[float])
if damping is not None:
a = damping(a, error=t1)
u_tot = scipy.linalg.expm(a)
return u, u_tot
ebcc.opt.gbrueckner.BruecknerGEBCC.transform_amplitudes(u, amplitudes=None)
Transform the amplitudes into the Brueckner orbital basis.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def transform_amplitudes(
self, u: SpinArrayType, amplitudes: Optional[Namespace[SpinArrayType]] = None
) -> Namespace[SpinArrayType]:
"""Transform the amplitudes into the Brueckner orbital basis.
Args:
u: Rotation matrix.
amplitudes: Cluster amplitudes.
Returns:
Transformed cluster amplitudes.
"""
if not amplitudes:
amplitudes = self.cc.amplitudes
nocc = self.cc.space.ncocc
ci = u[:nocc, :nocc]
ca = u[nocc:, nocc:]
# Transform T amplitudes:
for name, key, n in self.cc.ansatz.fermionic_cluster_ranks(spin_type=self.spin_type):
args: list[Union[tuple[int, ...], NDArray[T]]] = [
self.cc.amplitudes[name],
tuple(range(n * 2)),
]
for i in range(n):
args += [ci, (i, i + n * 2)]
for i in range(n):
args += [ca, (i + n, i + n * 3)]
args += [tuple(range(n * 2, n * 4))]
self.cc.amplitudes[name] = util.einsum(*args)
# Transform S amplitudes:
for name, key, n in self.cc.ansatz.bosonic_cluster_ranks(spin_type=self.spin_type):
raise util.ModelNotImplemented # TODO
# Transform U amplitudes:
for name, key, nf, nb in self.cc.ansatz.coupling_cluster_ranks(spin_type=self.spin_type):
raise util.ModelNotImplemented # TODO
return self.cc.amplitudes
ebcc.opt.gbrueckner.BruecknerGEBCC.get_t1_norm(amplitudes=None)
Get the norm of the T1 amplitude.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def get_t1_norm(self, amplitudes: Optional[Namespace[SpinArrayType]] = None) -> T:
"""Get the norm of the T1 amplitude.
Args:
amplitudes: Cluster amplitudes.
Returns:
Norm of the T1 amplitude.
"""
if not amplitudes:
amplitudes = self.cc.amplitudes
weight: T = np.linalg.norm(amplitudes["t1"])
return weight
ebcc.opt.gbrueckner.BruecknerGEBCC.mo_to_correlated(mo_coeff)
Transform the MO coefficients into the correlated basis.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def mo_to_correlated(self, mo_coeff: NDArray[T]) -> NDArray[T]:
"""Transform the MO coefficients into the correlated basis.
Args:
mo_coeff: MO coefficients.
Returns:
Correlated slice of MO coefficients.
"""
return mo_coeff[:, self.cc.space.correlated]
ebcc.opt.gbrueckner.BruecknerGEBCC.mo_update_correlated(mo_coeff, mo_coeff_corr)
Update the correlated slice of a set of MO coefficients.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def mo_update_correlated(self, mo_coeff: NDArray[T], mo_coeff_corr: NDArray[T]) -> NDArray[T]:
"""Update the correlated slice of a set of MO coefficients.
Args:
mo_coeff: MO coefficients.
mo_coeff_corr: Correlated slice of MO coefficients.
Returns:
Updated MO coefficients.
"""
mo_coeff = _put(
mo_coeff,
np.ix_(np.arange(mo_coeff.shape[0]), self.cc.space.correlated), # type: ignore
mo_coeff_corr,
)
return mo_coeff
ebcc.opt.gbrueckner.BruecknerGEBCC.update_coefficients(u_tot, mo_coeff, mo_coeff_ref)
Update the MO coefficients.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/gbrueckner.py
def update_coefficients(
self, u_tot: SpinArrayType, mo_coeff: NDArray[T], mo_coeff_ref: NDArray[T]
) -> NDArray[T]:
"""Update the MO coefficients.
Args:
u_tot: Total rotation matrix.
mo_coeff: New MO coefficients.
mo_coeff_ref: Reference MO coefficients.
Returns:
Updated MO coefficients.
"""
mo_coeff_new_corr = util.einsum("pi,ij->pj", mo_coeff_ref, u_tot)
mo_coeff_new = self.mo_update_correlated(mo_coeff, mo_coeff_new_corr)
return mo_coeff_new