Unrestricted Brueckner-orbital coupled cluster.

ebcc.opt.ubrueckner.BruecknerUEBCC(cc, options=None, **kwargs)

Bases: BaseBruecknerEBCC

Unrestricted Brueckner-orbital coupled cluster.

Initialise the Brueckner EBCC object.

Parameters:
  • cc (BaseEBCC) –

    Parent EBCC object.

  • options (Optional[BaseOptions], default: None ) –

    Options for the EOM calculation.

  • **kwargs (Any, default: {} ) –

    Additional keyword arguments used to update options.

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.ubrueckner.BruecknerUEBCC.get_rotation_matrix(u_tot=None, damping=None, t1=None)

Update the rotation matrix.

Also returns the total rotation matrix.

Parameters:
  • u_tot (Optional[SpinArrayType], default: None ) –

    Total rotation matrix.

  • damping (Optional[BaseDamping], default: None ) –

    Damping object.

  • t1 (Optional[SpinArrayType], default: None ) –

    T1 amplitude.

Returns:
  • tuple[SpinArrayType, SpinArrayType]

    Rotation matrix and total rotation matrix.

Source code in ebcc/opt/ubrueckner.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 = util.Namespace(
            aa=np.eye(self.cc.space[0].ncorr, dtype=types[float]),
            bb=np.eye(self.cc.space[1].ncorr, dtype=types[float]),
        )

    t1_block: Namespace[NDArray[T]] = util.Namespace()
    zocc = np.zeros((self.cc.space[0].ncocc, self.cc.space[0].ncocc))
    zvir = np.zeros((self.cc.space[0].ncvir, self.cc.space[0].ncvir))
    t1_block.aa = np.block([[zocc, -t1.aa], [np.transpose(t1.aa), zvir]])
    zocc = np.zeros((self.cc.space[1].ncocc, self.cc.space[1].ncocc))
    zvir = np.zeros((self.cc.space[1].ncvir, self.cc.space[1].ncvir))
    t1_block.bb = np.block([[zocc, -t1.bb], [np.transpose(t1.bb), zvir]])

    u = util.Namespace(aa=scipy.linalg.expm(t1_block.aa), bb=scipy.linalg.expm(t1_block.bb))

    u_tot.aa = u_tot.aa @ u.aa
    u_tot.bb = u_tot.bb @ u.bb
    if np.linalg.det(u_tot.aa) < 0:
        u_tot.aa = _put(
            u_tot.aa, np.ix_(np.arange(u_tot.aa.shape[0]), np.array([0])), -u_tot.aa[:, 0]
        )
    if np.linalg.det(u_tot.bb) < 0:
        u_tot.bb = _put(
            u_tot.bb, np.ix_(np.arange(u_tot.aa.shape[0]), np.array([0])), -u_tot.bb[:, 0]
        )

    a = np.concatenate(
        [np.ravel(scipy.linalg.logm(u_tot.aa)), np.ravel(scipy.linalg.logm(u_tot.bb))], axis=0
    )
    a: NDArray[T] = np.asarray(np.real(a), dtype=types[float])
    if damping is not None:
        xerr = np.concatenate([t1.aa.ravel(), t1.bb.ravel()])
        a = damping(a, error=xerr)

    u_tot.aa = scipy.linalg.expm(np.reshape(a[: u_tot.aa.size], u_tot.aa.shape))
    u_tot.bb = scipy.linalg.expm(np.reshape(a[u_tot.aa.size :], u_tot.bb.shape))

    return u, u_tot

ebcc.opt.ubrueckner.BruecknerUEBCC.transform_amplitudes(u, amplitudes=None)

Transform the amplitudes into the Brueckner orbital basis.

Parameters:
  • u (SpinArrayType) –

    Rotation matrix.

  • amplitudes (Optional[Namespace[SpinArrayType]], default: None ) –

    Cluster amplitudes.

Returns:
  • Namespace[SpinArrayType]

    Transformed cluster amplitudes.

Source code in ebcc/opt/ubrueckner.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[0].ncocc, self.cc.space[1].ncocc)
    ci = {"a": u.aa[: nocc[0], : nocc[0]], "b": u.bb[: nocc[1], : nocc[1]]}
    ca = {"a": u.aa[nocc[0] :, nocc[0] :], "b": u.bb[nocc[1] :, nocc[1] :]}

    # Transform T amplitudes:
    for name, key, n in self.cc.ansatz.fermionic_cluster_ranks(spin_type=self.spin_type):
        for comb in util.generate_spin_combinations(n, unique=True):
            args = [getattr(self.cc.amplitudes[name], comb), tuple(range(n * 2))]
        for i in range(n):
            args += [ci[comb[i]], (i, i + n * 2)]
        for i in range(n):
            args += [ca[comb[i + n]], (i + n, i + n * 3)]
        args += [tuple(range(n * 2, n * 4))]
        setattr(self.cc.amplitudes[name], comb, 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.ubrueckner.BruecknerUEBCC.get_t1_norm(amplitudes=None)

Get the norm of the T1 amplitude.

Parameters:
  • amplitudes (Optional[Namespace[SpinArrayType]], default: None ) –

    Cluster amplitudes.

Returns:
  • T

    Norm of the T1 amplitude.

Source code in ebcc/opt/ubrueckner.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_a = np.linalg.norm(amplitudes["t1"].aa)
    weight_b = np.linalg.norm(amplitudes["t1"].bb)
    weight: T = (weight_a**2 + weight_b**2) ** 0.5
    return weight

ebcc.opt.ubrueckner.BruecknerUEBCC.mo_to_correlated(mo_coeff)

Transform the MO coefficients into the correlated basis.

Parameters:
  • mo_coeff (tuple[NDArray[T], NDArray[T]]) –

    MO coefficients.

Returns:
  • tuple[NDArray[T], NDArray[T]]

    Correlated slice of MO coefficients.

Source code in ebcc/opt/ubrueckner.py
def mo_to_correlated(
    self, mo_coeff: tuple[NDArray[T], NDArray[T]]
) -> tuple[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[0][:, self.cc.space[0].correlated],
        mo_coeff[1][:, self.cc.space[1].correlated],
    )

ebcc.opt.ubrueckner.BruecknerUEBCC.mo_update_correlated(mo_coeff, mo_coeff_corr)

Update the correlated slice of a set of MO coefficients.

Parameters:
  • mo_coeff (tuple[NDArray[T], NDArray[T]]) –

    MO coefficients.

  • mo_coeff_corr (tuple[NDArray[T], NDArray[T]]) –

    Correlated slice of MO coefficients.

Returns:
  • tuple[NDArray[T], NDArray[T]]

    Updated MO coefficients.

Source code in ebcc/opt/ubrueckner.py
def mo_update_correlated(
    self,
    mo_coeff: tuple[NDArray[T], NDArray[T]],
    mo_coeff_corr: tuple[NDArray[T], NDArray[T]],
) -> tuple[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.
    """
    space = self.cc.space
    mo_coeff = (
        _put(
            mo_coeff[0],
            np.ix_(np.arange(mo_coeff[0].shape[0]), space[0].correlated),  # type: ignore
            mo_coeff_corr[0],
        ),
        _put(
            mo_coeff[1],
            np.ix_(np.arange(mo_coeff[1].shape[0]), space[1].correlated),  # type: ignore
            mo_coeff_corr[1],
        ),
    )
    return mo_coeff

ebcc.opt.ubrueckner.BruecknerUEBCC.update_coefficients(u_tot, mo_coeff, mo_coeff_ref)

Update the MO coefficients.

Parameters:
  • u_tot (SpinArrayType) –

    Total rotation matrix.

  • mo_coeff (tuple[NDArray[T], NDArray[T]]) –

    New MO coefficients.

  • mo_coeff_ref (tuple[NDArray[T], NDArray[T]]) –

    Reference MO coefficients.

Returns:
  • tuple[NDArray[T], NDArray[T]]

    Updated MO coefficients.

Source code in ebcc/opt/ubrueckner.py
def update_coefficients(
    self,
    u_tot: SpinArrayType,
    mo_coeff: tuple[NDArray[T], NDArray[T]],
    mo_coeff_ref: tuple[NDArray[T], NDArray[T]],
) -> tuple[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[0], u_tot.aa),
        util.einsum("pi,ij->pj", mo_coeff_ref[1], u_tot.bb),
    )
    mo_coeff_new = self.mo_update_correlated(mo_coeff, mo_coeff_new_corr)
    return mo_coeff_new