Base classes for ebcc.opt
.
ebcc.opt.base.BaseOptions(e_tol=1e-08, t_tol=1e-08, max_iter=20, diis_space=9, diis_min_space=1, damping=0.0)
dataclass
Bases: _BaseOptions
Options for Brueckner-orbital calculations.
Parameters: |
|
---|
ebcc.opt.base.BaseBruecknerEBCC(cc, options=None, **kwargs)
Bases: ABC
Base class for 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.base.BaseBruecknerEBCC.spin_type: str
property
Get the spin type.
ebcc.opt.base.BaseBruecknerEBCC.name: str
property
Get the name of the method.
ebcc.opt.base.BaseBruecknerEBCC.kernel()
Run the Bruckner-orbital coupled cluster calculation.
Returns: |
|
---|
Source code in ebcc/opt/base.py
def kernel(self) -> float:
"""Run the Bruckner-orbital coupled cluster calculation.
Returns:
Correlation energy.
"""
timer = util.Timer()
# Make sure the initial CC calculation is converged:
if not self.cc.converged:
with lib.temporary_env(self.cc, log=NullLogger()):
self.cc.kernel()
# Set up DIIS:
damping = self.Damping(options=self.options)
# Initialise coefficients:
mo_coeff_new: NDArray[T] = np.copy(np.asarray(self.cc.mo_coeff, dtype=types[float]))
mo_coeff_ref: NDArray[T] = np.copy(np.asarray(self.cc.mo_coeff, dtype=types[float]))
mo_coeff_ref = self.mo_to_correlated(mo_coeff_ref)
u_tot = None
self.cc.log.output("Solving for Brueckner orbitals.")
self.cc.log.debug("")
self.log.info(
f"{ANSI.B}{'Iter':>4s} {'Energy (corr.)':>16s} {'Energy (tot.)':>18s} "
f"{'Conv.':>8s} {'Δ(Energy)':>13s} {'|T1|':>13s}{ANSI.R}"
)
self.log.info(
f"%4d %16.10f %18.10f {[ANSI.r, ANSI.g][self.cc.converged]}%8r{ANSI.R}",
0,
self.cc.e_corr,
self.cc.e_tot,
self.cc.converged,
)
converged = False
for niter in range(1, self.options.max_iter + 1):
# Update rotation matrix:
u, u_tot = self.get_rotation_matrix(u_tot=u_tot, damping=damping)
# Update MO coefficients:
mo_coeff_new = self.update_coefficients(u_tot, mo_coeff_new, mo_coeff_ref)
# Transform mean-field and amplitudes:
self.mf.mo_coeff = numpy.asarray(mo_coeff_new)
self.mf.e_tot = self.mf.energy_tot()
amplitudes = self.transform_amplitudes(u)
# Run CC calculation:
e_prev = self.cc.e_tot
with lib.temporary_env(self.cc, log=NullLogger()):
self.cc.__class__.__init__(
self.cc,
self.mf,
log=self.cc.log,
ansatz=self.cc.ansatz,
space=self.cc.space,
omega=self.cc.omega,
g=self.cc.bare_g,
G=self.cc.bare_G,
options=self.cc.options,
)
self.cc.amplitudes = amplitudes
self.cc.kernel()
de = abs(e_prev - self.cc.e_tot)
dt = self.get_t1_norm()
# Log the iteration:
converged_e = bool(de < self.options.e_tol)
converged_t = bool(dt < self.options.t_tol)
self.log.info(
f"%4s %16.10f %18.10f {[ANSI.r, ANSI.g][int(converged)]}%8r{ANSI.R}"
f" {[ANSI.r, ANSI.g][int(converged_e)]}%13.3e{ANSI.R}"
f" {[ANSI.r, ANSI.g][int(converged_t)]}%13.3e{ANSI.R}",
niter,
self.cc.e_corr,
self.cc.e_tot,
self.cc.converged,
de,
dt,
)
# Check for convergence:
converged = converged_e and converged_t
if converged:
self.log.debug("")
self.log.output(f"{ANSI.g}Converged.{ANSI.R}")
break
else:
self.log.debug("")
self.log.warning(f"{ANSI.r}Failed to converge.{ANSI.R}")
self.cc.log.debug("")
self.cc.log.output("E(corr) = %.10f", self.cc.e_corr)
self.cc.log.output("E(tot) = %.10f", self.cc.e_tot)
self.cc.log.debug("")
self.cc.log.debug("Time elapsed: %s", timer.format_time(timer()))
self.cc.log.debug("")
return self.cc.e_corr
ebcc.opt.base.BaseBruecknerEBCC.get_rotation_matrix(u_tot=None, damping=None, t1=None)
abstractmethod
Update the rotation matrix.
Also returns the total rotation matrix.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
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.
"""
pass
ebcc.opt.base.BaseBruecknerEBCC.transform_amplitudes(u, amplitudes=None)
abstractmethod
Transform the amplitudes into the Brueckner orbital basis.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
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.
"""
pass
ebcc.opt.base.BaseBruecknerEBCC.get_t1_norm(amplitudes=None)
abstractmethod
Get the norm of the T1 amplitude.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
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.
"""
pass
ebcc.opt.base.BaseBruecknerEBCC.mo_to_correlated(mo_coeff)
abstractmethod
Transform the MO coefficients into the correlated basis.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
def mo_to_correlated(self, mo_coeff: Any) -> Any:
"""Transform the MO coefficients into the correlated basis.
Args:
mo_coeff: MO coefficients.
Returns:
Correlated slice of MO coefficients.
"""
pass
ebcc.opt.base.BaseBruecknerEBCC.mo_update_correlated(mo_coeff, mo_coeff_corr)
abstractmethod
Update the correlated slice of a set of MO coefficients.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
def mo_update_correlated(self, mo_coeff: Any, mo_coeff_corr: Any) -> Any:
"""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.
"""
pass
ebcc.opt.base.BaseBruecknerEBCC.update_coefficients(u_tot, mo_coeff_new, mo_coeff_ref)
abstractmethod
Update the MO coefficients.
Parameters: |
|
---|
Returns: |
|
---|
Source code in ebcc/opt/base.py
@abstractmethod
def update_coefficients(
self, u_tot: SpinArrayType, mo_coeff_new: Any, mo_coeff_ref: Any
) -> Any:
"""Update the MO coefficients.
Args:
u_tot: Total rotation matrix.
mo_coeff_new: New MO coefficients.
mo_coeff_ref: Reference MO coefficients.
Returns:
Updated MO coefficients.
"""
pass