k-UpCCGSD¶
- class tencirchem.KUPCCGSD(mol: Mole, active_space: Tuple[int, int] | None = None, mo_coeff: ndarray | None = None, k: int = 3, n_tries: int = 1, engine: str | None = None, run_hf: bool = True, run_mp2: bool = True, run_ccsd: bool = True, run_fci: bool = True)[source]¶
Bases:
UCC
Run \(k\)-UpCCGSD calculation. The interfaces are similar to
UCCSD
.- __init__(mol: Mole, active_space: Tuple[int, int] | None = None, mo_coeff: ndarray | None = None, k: int = 3, n_tries: int = 1, engine: str | None = None, run_hf: bool = True, run_mp2: bool = True, run_ccsd: bool = True, run_fci: bool = True)[source]¶
Initialize the class with molecular input.
- Parameters:
mol (Mole) – The molecule as PySCF
Mole
object.active_space (Tuple[int, int], optional) – Active space approximation. The first integer is the number of electrons and the second integer is the number or spatial-orbitals. Defaults to None.
mo_coeff (np.ndarray, optional) – Molecule coefficients. If provided then RHF is skipped. Can be used in combination with the
init_state
attribute. Defaults to None which means RHF orbitals are used.k (int, optional) – The number of layers in the ansatz. Defaults to 3
n_tries (int, optional) – The number of different initial points used for VQE calculation. For large circuits usually a lot of runs are required for good accuracy. Defaults to 1.
engine (str, optional) – The engine to run the calculation. See Engines for details.
run_hf (bool, optional) – Whether run HF for molecule orbitals. Defaults to
True
.run_mp2 (bool, optional) – Whether run MP2 for energy reference. Defaults to
True
.run_ccsd (bool, optional) – Whether run CCSD for energy reference. Defaults to
True
.run_fci (bool, optional) – Whether run FCI for energy reference. Defaults to
True
.
See also
- apply_excitation(state: Any, ex_op: Tuple, engine: str | None = None) Any ¶
Apply a given excitation operator to a given state.
- Parameters:
state (Tensor) – The input state in statevector or CI vector form.
ex_op (tuple of ints) – The excitation operator.
engine (str, optional) – The engine to use. Defaults to
None
, which usesself.engine
.
- Returns:
tensor – The resulting tensor.
- Return type:
Tensor
Examples
>>> from tencirchem import UCC >>> from tencirchem.molecule import h2 >>> ucc = UCC(h2) >>> ucc.apply_excitation([1, 0, 0, 0], (3, 1, 0, 2)) array([0, 0, 0, 1])
- classmethod as_pyscf_solver(config_function: Callable | None = None, **kwargs)¶
Converts the
UCC
class to a PySCF FCI solver.- Parameters:
config_function (callable) – A function to configure the
UCC
instance. Accepts theUCC
instance and modifies it inplace beforekernel()
is called.kwargs – Other arguments to be passed to the
__init__()
function such asengine
.
- Return type:
FCISolver
Examples
>>> from pyscf.mcscf import CASSCF >>> from tencirchem import UCCSD >>> from tencirchem.molecule import h8 >>> # normal PySCF workflow >>> hf = h8.HF() >>> round(hf.kernel(), 8) -4.14961853 >>> casscf = CASSCF(hf, 2, 2) >>> # set the FCI solver for CASSCF to be UCCSD >>> casscf.fcisolver = UCCSD.as_pyscf_solver() >>> round(casscf.kernel()[0], 8) -4.16647335
- civector(params: Any | None = None, engine: str | None = None) Any ¶
Evaluate the configuration interaction (CI) vector.
- Parameters:
params (Tensor, optional) – The circuit parameters. Defaults to None, which uses the optimized parameter and
kernel()
must be called before.engine (str, optional) – The engine to use. Defaults to
None
, which usesself.engine
.
- Returns:
civector – Corresponding CI vector
- Return type:
Tensor
See also
statevector
Evaluate the circuit state vector.
energy
Evaluate the total energy.
energy_and_grad
Evaluate the total energy and parameter gradients.
Examples
>>> from tencirchem import UCCSD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> uccsd.civector([0, 0]) # HF state array([1., 0., 0., 0.])
- property civector_size: int¶
The size of the CI vector.
- property e_kupccgsd¶
Returns \(k\)-UpCCGSD energy
- property e_ucc: float¶
Returns UCC energy
- embed_rdm_cas(rdm_cas)¶
Embed CAS RDM into RDM of the whole system
- energy(params: Any | None = None, engine: str | None = None) float ¶
Evaluate the total energy.
- Parameters:
params (Tensor, optional) – The circuit parameters. Defaults to None, which uses the optimized parameter and
kernel()
must be called before.engine (str, optional) – The engine to use. Defaults to
None
, which usesself.engine
.
- Returns:
energy – Total energy
- Return type:
float
See also
civector
Get the configuration interaction (CI) vector.
statevector
Evaluate the circuit state vector.
energy_and_grad
Evaluate the total energy and parameter gradients.
Examples
>>> from tencirchem import UCCSD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> round(uccsd.energy([0, 0]), 8) # HF state -1.11670614
- energy_and_grad(params: Any | None = None, engine: str | None = None) Tuple[float, Any] ¶
Evaluate the total energy and parameter gradients.
- Parameters:
params (Tensor, optional) – The circuit parameters. Defaults to None, which uses the optimized parameter and
kernel()
must be called before.engine (str, optional) – The engine to use. Defaults to
None
, which usesself.engine
.
- Returns:
energy (float) – Total energy
grad (Tensor) – The parameter gradients
See also
civector
Get the configuration interaction (CI) vector.
statevector
Evaluate the circuit state vector.
energy
Evaluate the total energy.
Examples
>>> from tencirchem import UCCSD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> e, g = uccsd.energy_and_grad([0, 0]) >>> round(e, 8) -1.11670614 >>> g array([..., ...])
- classmethod from_integral(int1e: ndarray, int2e: ndarray, n_elec: int, e_core: float = 0, ovlp: ndarray | None = None, **kwargs)¶
Construct UCC classes from electron integrals.
- Parameters:
int1e (np.ndarray) – One-body integral in spatial orbital.
int2e (np.ndarray) – Two-body integral, in spatial orbital, chemists’ notation, and without considering symmetry.
n_elec (int) – The number of electrons
e_core (float, optional) – The nuclear energy or core energy if active space approximation is involved. Defaults to 0.
ovlp (np.ndarray) – The overlap integral. Defaults to
None
and identity matrix is used.kwargs – Other arguments to be passed to the
__init__()
function such asengine
.
- Returns:
ucc – A UCC instance
- Return type:
- get_addr(bitstring: str) int ¶
Get the address (index) of a CI bitstring in the CI vector.
- Parameters:
bitstring (str) – The bitstring such as
"0101"
.- Returns:
address – The bitstring address.
- Return type:
int
Examples
>>> from tencirchem import UCCSD, PUCCD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> uccsd.get_addr("0101") # the HF state 0 >>> uccsd.get_addr("1010") 3 >>> puccd = PUCCD(h2) >>> puccd.get_addr("01") # the HF state 0 >>> puccd.get_addr("10") 1
- get_ci_strings(strs2addr: bool = False) ndarray ¶
Get the CI bitstrings for all configurations in the CI vector.
- Parameters:
strs2addr (bool, optional.) – Whether return the reversed mapping for one spin sector. Defaults to
False
.- Returns:
cistrings (np.ndarray) – The CI bitstrings.
string_addr (np.ndarray) – The address of the string in one spin sector. Returned when
strs2addr
is set toTrue
.
Examples
>>> from tencirchem import UCCSD, PUCCD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> uccsd.get_ci_strings() array([ 5, 6, 9, 10], dtype=uint64) >>> [f"{bin(i)[2:]:0>4}" for i in uccsd.get_ci_strings()] ['0101', '0110', '1001', '1010'] >>> uccsd.get_ci_strings(True)[1] # only one spin sector array([0, 0, 1, 0], dtype=uint64)
- get_circuit(params: Any | None = None, decompose_multicontrol: bool = False, trotter: bool = False) Circuit ¶
Get the circuit as TensorCircuit
Circuit
object- Parameters:
params (Tensor, optional) – The circuit parameters. Defaults to None, which uses the optimized parameter. If
kernel()
is not called before, the initial guess is used.decompose_multicontrol (bool, optional) – Whether decompose the Multicontrol gate in the circuit into CNOT gates. Defaults to False.
trotter (bool, optional) – Whether Trotterize the UCC factor into Pauli strings. Defaults to False.
- Returns:
circuit – The quantum circuit.
- Return type:
tc.Circuit
- get_energy_dataframe() DataFrame ¶
Returns energy information dataframe
- get_ex1_ops(t1: ndarray | None = None) Tuple[List[Tuple], List[int], ndarray] [source]¶
Get generalized one-body excitation operators.
- Parameters:
t1 (np.ndarray, optional) – Not used. Kept for consistency with the parent method.
- Returns:
ex_op (List[Tuple]) – The excitation operators. Each operator is represented by a tuple of ints.
param_ids (List[int]) – The mapping from excitations to parameters.
init_guess (np.ndarray) – The initial guess for the parameters.
See also
get_ex2_ops
Get generalized paired two-body excitation operators.
get_ex_ops
Get one-body and two-body excitation operators for \(k\)-UpCCGSD ansatz.
- get_ex2_ops(t2: ndarray | None = None) Tuple[List[Tuple], List[int], ndarray] [source]¶
Get generalized paired two-body excitation operators.
- Parameters:
t2 (np.ndarray, optional) – Not used. Kept for consistency with the parent method.
- Returns:
ex_op (List[Tuple]) – The excitation operators. Each operator is represented by a tuple of ints.
param_ids (List[int]) – The mapping from excitations to parameters.
init_guess (np.ndarray) – The initial guess for the parameters.
See also
get_ex1_ops
Get one-body excitation operators.
get_ex_ops
Get one-body and two-body excitation operators for \(k\)-UpCCGSD ansatz.
- get_ex_ops(t1: ndarray | None = None, t2: ndarray | None = None) Tuple[List[Tuple], List[int], ndarray] [source]¶
Get one-body and two-body excitation operators for \(k\)-UpCCGSD ansatz. The excitations are generalized and two-body excitations are restricted to paired ones. Initial guesses are generated randomly.
- Parameters:
t1 (np.ndarray, optional) – Not used. Kept for consistency with the parent method.
t2 (np.ndarray, optional) – Not used. Kept for consistency with the parent method.
- Returns:
ex_op (List[Tuple]) – The excitation operators. Each operator is represented by a tuple of ints.
param_ids (List[int]) – The mapping from excitations to parameters.
init_guess (np.ndarray) – The initial guess for the parameters.
See also
get_ex1_ops
Get generalized one-body excitation operators.
get_ex2_ops
Get generalized paired two-body excitation operators.
Examples
>>> from tencirchem import KUPCCGSD >>> from tencirchem.molecule import h2 >>> kupccgsd = KUPCCGSD(h2) >>> ex_op, param_ids, init_guess = kupccgsd.get_ex_ops() >>> ex_op [(1, 3, 2, 0), (3, 2), (1, 0), (1, 3, 2, 0), (3, 2), (1, 0), (1, 3, 2, 0), (3, 2), (1, 0)] >>> param_ids [0, 1, 1, 2, 3, 3, 4, 5, 5] >>> init_guess array([...])
- get_excitation_dataframe() DataFrame ¶
- get_init_state_dataframe(coeff_epsilon: float = 1e-12) DataFrame ¶
Returns initial state information dataframe.
- Parameters:
coeff_epsilon (float, optional) – The threshold to screen out states with small coefficients. Defaults to 1e-12.
- Return type:
pd.DataFrame
See also
init_state
The circuit initial state before applying the excitation operators.
Examples
>>> from tencirchem import UCC >>> from tencirchem.molecule import h2 >>> ucc = UCC(h2) >>> ucc.init_state = [0.707, 0, 0, 0.707] >>> ucc.get_init_state_dataframe() configuration coefficient 0 0101 0.707 1 1010 0.707
- get_opt_function(with_time: bool = False) Callable | Tuple[Callable, float] ¶
Returns the cost function in SciPy format for optimization. The gradient is included. Basically a wrapper to
energy_and_grad()
.- Parameters:
with_time (bool, optional) – Whether return staging time. Defaults to False.
- Returns:
opt_function (Callable) – The optimization cost function in SciPy format.
time (float) – Staging time. Returned when
with_time
is set toTrue
.
- property h_fermion_op: FermionOperator¶
Hamiltonian as openfermion.FermionOperator
- property h_qubit_op: QubitOperator¶
Hamiltonian as openfermion.QubitOperator, mapped by Jordan-Wigner transformation.
- property init_state: Any¶
The circuit initial state before applying the excitation operators. Usually RHF.
See also
get_init_state_dataframe
Returns initial state information dataframe.
- kernel()[source]¶
The kernel to perform the VQE algorithm. The L-BFGS-B method in SciPy is used for optimization and configuration is possible by setting the
self.scipy_minimize_options
attribute.- Returns:
e – The optimized energy
- Return type:
float
- make_rdm1(statevector: Any | None = None, basis: str = 'AO') ndarray ¶
Evaluate the spin-traced one-body reduced density matrix (1RDM).
\[\textrm{1RDM}[p,q] = \langle p_{\alpha}^\dagger q_{\alpha} \rangle + \langle p_{\beta}^\dagger q_{\beta} \rangle\]If active space approximation is employed, returns the full RDM of all orbitals.
- Parameters:
statevector (Tensor, optional) – Custom system state. Could be CI vector or state vector. Defaults to None, which uses the optimized state by
civector()
.basis (str, optional) – One of
"AO"
or"MO"
. Defaults to"AO"
, which is for consistency with PySCF.
- Returns:
rdm1 – The spin-traced one-body RDM.
- Return type:
np.ndarray
See also
make_rdm2
Evaluate the spin-traced two-body reduced density matrix (2RDM).
- make_rdm2(statevector: Any | None = None, basis: str = 'AO') ndarray ¶
Evaluate the spin-traced two-body reduced density matrix (2RDM).
\[\begin{split}\begin{aligned} \textrm{2RDM}[p,q,r,s] & = \langle p_{\alpha}^\dagger r_{\alpha}^\dagger s_{\alpha} q_{\alpha} \rangle + \langle p_{\beta}^\dagger r_{\alpha}^\dagger s_{\alpha} q_{\beta} \rangle \\ & \quad + \langle p_{\alpha}^\dagger r_{\beta}^\dagger s_{\beta} q_{\alpha} \rangle + \langle p_{\beta}^\dagger r_{\beta}^\dagger s_{\beta} q_{\beta} \rangle \end{aligned}\end{split}\]If active space approximation is employed, returns the full RDM of all orbitals.
- Parameters:
statevector (Tensor, optional) – Custom system state. Could be CI vector or state vector. Defaults to None, which uses the optimized state by
civector()
.basis (str, optional) – One of
"AO"
or"MO"
. Defaults to"AO"
, which is for consistency with PySCF.
- Returns:
rdm2 – The spin-traced two-body RDM.
- Return type:
np.ndarray
See also
make_rdm1
Evaluate the spin-traced one-body reduced density matrix (1RDM).
Examples
>>> import numpy as np >>> from tencirchem import UCC >>> from tencirchem.molecule import h2 >>> ucc = UCC(h2) >>> state = [1, 0, 0, 0] ## HF state >>> rdm1 = ucc.make_rdm1(state, basis="MO") >>> rdm2 = ucc.make_rdm2(state, basis="MO") >>> e_hf = ucc.int1e.ravel() @ rdm1.ravel() + 1/2 * ucc.int2e.ravel() @ rdm2.ravel() >>> np.testing.assert_allclose(e_hf + ucc.e_nuc, ucc.e_hf, atol=1e-10)
- property n_params: int¶
The number of parameter in the ansatz/circuit.
- property no: int¶
The number of occupied orbitals.
- property nv: int¶
The number of virtual (unoccupied orbitals).
- property param_ids: List[int]¶
The mapping from excitations operators to parameters.
- property param_to_ex_ops¶
- property params: Any¶
The circuit parameters.
- print_ansatz()¶
- print_circuit()¶
Prints the circuit information. If you wish to print the circuit diagram, use
get_circuit()
and then calldraw()
such asprint(ucc.get_circuit().draw())
.
- print_energy()¶
- print_excitations()¶
- print_init_state()¶
- print_summary(include_circuit: bool = False)¶
Print a summary of the class.
- Parameters:
include_circuit (bool) – Whether include the circuit section.
- statevector(params: Any | None = None, engine: str | None = None) Any ¶
Evaluate the circuit state vector.
- Parameters:
params (Tensor, optional) – The circuit parameters. Defaults to None, which uses the optimized parameter and
kernel()
must be called before.engine (str, optional) – The engine to use. Defaults to
None
, which usesself.engine
.
- Returns:
statevector – Corresponding state vector
- Return type:
Tensor
See also
civector
Evaluate the configuration interaction (CI) vector.
energy
Evaluate the total energy.
energy_and_grad
Evaluate the total energy and parameter gradients.
Examples
>>> from tencirchem import UCCSD >>> from tencirchem.molecule import h2 >>> uccsd = UCCSD(h2) >>> uccsd.statevector([0, 0]) # HF state array([0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.])
- property statevector_size: int¶
The size of the statevector.