mindquantum.algorithm.nisq
NISQ algorithms
- class mindquantum.algorithm.nisq.Ansatz(name, n_qubits, *args, **kwargs)[source]
Basic class for Ansatz.
- Parameters
- property circuit
Get the quantum circuit of this ansatz.
- Returns
Circuit, the quantum circuit of this ansatz.
- class mindquantum.algorithm.nisq.HardwareEfficientAnsatz(n_qubits, single_rot_gate_seq, entangle_gate=X, entangle_mapping='linear', depth=1)[source]
Hardware efficient ansatz is a kind of ansatz that can be easily implement on quantum chip.
The hardware efficient is constructed by a layer of single qubit rotation gate and a layer of two qubits entanglement gate. The single qubit rotation gate layer is constructed by one or several rotation gate that act on every qubit. The two qubits entanglement gate layer is constructed by CNOT, CZ, XX, YY, ZZ, etc. acting on entangle_mapping. For more detail, please refers https://www.nature.com/articles/nature23879.
- Parameters
n_qubits (int) – number of qubit that this ansatz act on.
single_rot_gate_seq (list[BasicGate]) – A list of parameterized rotation gate that act on each qubit.
entangle_gate (BasicGate) – The non parameterized entanglement gate. If it is a single qubit gate, than the control version will be used. Default: XGate.
entangle_mapping (Union[str, list[tuple[int]]]) – The entanglement mapping of entanglement gate. ‘linear’ means the entanglement gate will be act on every neighboring qubits. ‘all’ means the entanglemtn gate will be act on any two qbuits. Besides, you can specific which two qubits you want to do entanglement by setting the entangle_mapping to a list of two qubits tuple. Default: “linear”.
depth (int) – The depth of ansatz. Default: 1.
Examples
>>> from mindquantum.algorithm.nisq.chem import HardwareEfficientAnsatz >>> from mindquantum import RY, RZ, Z >>> hea = HardwareEfficientAnsatz(3, [RY, RZ], Z, [(0, 1), (0, 2)]) >>> hea.circuit q0: ──RY(d0_n0_0)────RZ(d0_n0_1)────●────●────RY(d1_n0_0)────RZ(d1_n0_1)── │ │ q1: ──RY(d0_n1_0)────RZ(d0_n1_1)────Z────┼────RY(d1_n1_0)────RZ(d1_n1_1)── │ q2: ──RY(d0_n2_0)────RZ(d0_n2_1)─────────Z────RY(d1_n2_0)────RZ(d1_n2_1)──
- class mindquantum.algorithm.nisq.IQPEncoding(n_feature, first_rotation_gate=RZ, second_rotation_gate=RZ, num_repeats=1)[source]
General IQP Encoding.
- Parameters
Examples
>>> from mindquantum.algorithm.library import IQPEncoding >>> iqp = IQPEncoding(3) >>> iqp q0: ──H────RZ(alpha0)────●───────────────────────────●─────────────────────────────────── │ │ q1: ──H────RZ(alpha1)────X────RZ(alpha0 * alpha1)────X────●───────────────────────────●── │ │ q2: ──H────RZ(alpha2)─────────────────────────────────────X────RZ(alpha1 * alpha2)────X── >>> iqp.circuit.params_name ['alpha0', 'alpha1', 'alpha2', 'alpha0 * alpha1', 'alpha1 * alpha2'] >>> iqp.circuit.params_name >>> a = np.array([0, 1, 2]) >>> iqp.data_preparation(a) array([0, 1, 2, 0, 2]) >>> iqp.circuit.get_qs(pr=iqp.data_preparation(a)) array([-0.28324704-0.21159186j, -0.28324704-0.21159186j, 0.31027229+0.16950252j, 0.31027229+0.16950252j, 0.02500938+0.35266773j, 0.02500938+0.35266773j, 0.31027229+0.16950252j, 0.31027229+0.16950252j])
- data_preparation(data)[source]
The IQPEncoding ansatz provide a ansatz to encode classical data into quantum state. This method will prepare the classical data into suitable dimension for IQPEncoding. Suppose the origin data has \(n\) features, then the output data will have \(2n-1\) features, with first \(n\) features keep the same and for \(m > n\),
\[\text{data}_m = \text{data}_{m - n} * \text{data}_{m - n - 1}\]- Parameters
data ([list, numpy.ndarray]) – The classical data need to encode with IQPEncoding ansatz.
- Returns
numpy.ndarray, the prepared data that is suitable for this ansatz.
- class mindquantum.algorithm.nisq.Max2SATAnsatz(clauses, depth=1)[source]
The Max-2-SAT ansatz. For more detail, please refers to https://arxiv.org/pdf/1906.11259.pdf.
\[U(\beta, \gamma) = e^{-\beta_pH_b}e^{-\gamma_pH_c} \cdots e^{-\beta_0H_b}e^{-\gamma_0H_c}H^{\otimes n}\]Where,
\[H_b = \sum_{i\in n}X_{i}, H_c = \sum_{l\in m}P(l)\]Here \(n\) is the number of Boolean variables and \(m\) is the number of total clauses and \(P(l)\) is rank-one projector.
- Parameters
Examples
>>> from mindquantum.algorithm.nisq.qaoa import Max2SATAnsatz >>> clauses = [(2, -3)] >>> max2sat = Max2SATAnsatz(clauses, 1) >>> max2sat.circuit q1: ──H─────RZ(0.5*beta_0)────●───────────────────────●────RX(alpha_0)── │ │ q2: ──H────RZ(-0.5*beta_0)────X────RZ(-0.5*beta_0)────X────RX(alpha_0)──
>>> max2sat.hamiltonian 0.25 [] + 0.25 [Z1] + -0.25 [Z1 Z2] + -0.25 [Z2] >>> sats = max2sat.get_sat(4, np.array([4, 1])) >>> sats ['001', '000', '011', '010'] >>> for i in sats: >>> print(f'sat value: {max2sat.get_sat_value([i])}') sat value: 1 sat value: 0 sat value: 2 sat value: 1
- get_sat(max_n, weight)[source]
Get the strings of this max-2-sat problem.
- Parameters
max_n (int) – how many strings you want.
weight (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) – parameter value for Max-2-SAT ansatz.
- Returns
list, a list of strings.
- get_sat_value(string)[source]
Get the sat values for given strings. The string is a str that satisfies all the clauses of the given max-2-sat problem.
- Parameters
string (str) – a string of the max-2-sat problem consided.
- Returns
int, sat_value under the given string.
- property hamiltonian
Get the hamiltonian of this max-2-sat problem.
- Returns
QubitOperator, hamiltonian of this max-2-sat problem.
- class mindquantum.algorithm.nisq.MaxCutAnsatz(graph, depth=1)[source]
The MaxCut ansatz. For more detail, please refers to https://arxiv.org/pdf/1411.4028.pdf.
\[U(\beta, \gamma) = e^{-\beta_pH_b}e^{-\gamma_pH_c} \cdots e^{-\beta_0H_b}e^{-\gamma_0H_c}H^{\otimes n}\]Where,
\[H_b = \sum_{i\in n}X_{i}, H_c = \sum_{(i,j)\in C}Z_iZ_j\]Here \(n\) is the set of nodes and \(C\) is the set of edges of the graph.
- Parameters
Examples
>>> from mindquantum.algorithm.nisq.qaoa import MaxCutAnsatz >>> graph = [(0, 1), (1, 2), (0, 2)] >>> maxcut = MaxCutAnsatz(graph, 1) >>> maxcut.circuit q0: ──H────ZZ(beta_0)──────────────────ZZ(beta_0)────RX(alpha_0)── │ │ q1: ──H────ZZ(beta_0)────ZZ(beta_0)────────┼─────────RX(alpha_0)── │ │ q2: ──H──────────────────ZZ(beta_0)────ZZ(beta_0)────RX(alpha_0)──
>>> maxcut.hamiltonian 1.5 [] + -0.5 [Z0 Z1] + -0.5 [Z0 Z2] + -0.5 [Z1 Z2] >>> maxcut.hamiltonian >>> partitions = maxcut.get_partition(5, np.array([4, 1])) >>> for i in partitions: >>> print(f'partition: left: {i[0]}, right: {i[1]}, cut value: {maxcut.get_cut_value(i)}') partition: left: [2], right: [0, 1], cut value: 2 partition: left: [0, 1], right: [2], cut value: 2 partition: left: [0], right: [1, 2], cut value: 2 partition: left: [0, 1, 2], right: [], cut value: 0 partition: left: [], right: [0, 1, 2], cut value: 0
- get_cut_value(partition)[source]
Get the cut values for given partitions. The partition is a list that contains two lists, each list contains the nodes of the given graph.
- Parameters
partition (list) – a partition of the graph consided.
- Returns
int, cut_value under the given partition.
- get_partition(max_n, weight)[source]
Get the partitions of this max-cut problem.
- Parameters
max_n (int) – how many partitions you want.
weight (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) – parameter value for max-cut ansatz.
- Returns
list, a list of partitions.
- property hamiltonian
Get the hamiltonian of this max cut problem.
- Returns
QubitOperator, hamiltonian of this max cut problem.
- class mindquantum.algorithm.nisq.QubitUCCAnsatz(n_qubits=None, n_electrons=None, occ_orb=None, vir_orb=None, generalized=False, trotter_step=1)[source]
Qubit Unitary Coupled-Cluster (qUCC) ansatz is a variant of unitary coupled-cluster ansatz which uses qubit excitation operators instead of Fermion excitation operators. The Fock space spanned by qubit excitation operators is equivalent as Fermion operators, therefore the exact wave function can be approximated using qubit excitation operators at the expense of a higher order of Trotterization.
The greatest advantange of qUCC is that the number of CNOT gates is much smaller than the original version of UCC, even using a 3rd or 4th order Trotterization. Also, the accuracy is greatly improved despite that the number of variational parameters is increased.
Note
The Hartree-Fock circuit is not included. Currently, generalized=True is not allowed since the theory needs verification. Reference: Yordan S. Yordanov et al. Phys. Rev. A, 102, 062612 (2020)
- Parameters
n_qubits (int) – The number of qubits (spin-orbitals) in the simulation. Default: None.
n_electrons (int) – The number of electrons of the given molecule. Default: None.
occ_orb (list) – Indices of manually assigned occupied spatial orbitals. Default: None.
vir_orb (list) – Indices of manually assigned virtual spatial orbitals. Default: None.
generalized (bool) – Whether to use generalized excitations which do not distinguish occupied or virtual orbitals (qUCCGSD). Currently, generalized=True is not allowed since the theory needs verification. Default: False.
trotter_step (int) – The number of Trotter steps. Default is one. It is recommended to set a value larger than or equal to 2 to achieve a good accuracy. Default: 1.
Examples
>>> from mindquantum.algorithm.nisq.chem import QubitUCCAnsatz >>> QubitUCCAnsatz().n_qubits 0 >>> qucc = QubitUCCAnsatz(4, 2, trotter_step=2) >>> qucc.circuit[:10] q0: ──X──────────●──────────X───────────────────────────────X──────────●──────────X─────── │ │ │ │ │ │ q1: ──┼──────────┼──────────┼────X──────────●──────────X────┼──────────┼──────────┼────X── │ │ │ │ │ │ │ │ │ │ q2: ──●────RY(t_0_q_s_0)────●────●────RY(t_0_q_s_1)────●────┼──────────┼──────────┼────┼── │ │ │ │ q3: ────────────────────────────────────────────────────────●────RY(t_0_q_s_2)────●────●── >>> qucc.n_qubits 4 >>> qucc_2 = QubitUCCAnsatz(occ_orb=[0, 1], vir_orb=[2]) >>> qucc_2.operator_pool [-1.0*t_0_q_s_0 [Q0^ Q4] + 1.0*t_0_q_s_0 [Q4^ Q0] , -1.0*t_0_q_s_1 [Q1^ Q4] + 1.0*t_0_q_s_1 [Q4^ Q1] , -1.0*t_0_q_s_2 [Q2^ Q4] + 1.0*t_0_q_s_2 [Q4^ Q2] , -1.0*t_0_q_s_3 [Q3^ Q4] + 1.0*t_0_q_s_3 [Q4^ Q3] , -1.0*t_0_q_s_4 [Q0^ Q5] + 1.0*t_0_q_s_4 [Q5^ Q0] , -1.0*t_0_q_s_5 [Q1^ Q5] + 1.0*t_0_q_s_5 [Q5^ Q1] , -1.0*t_0_q_s_6 [Q2^ Q5] + 1.0*t_0_q_s_6 [Q5^ Q2] , -1.0*t_0_q_s_7 [Q3^ Q5] + 1.0*t_0_q_s_7 [Q5^ Q3] , -1.0*t_0_q_d_0 [Q1^ Q0^ Q5 Q4] + 1.0*t_0_q_d_0 [Q5^ Q4^ Q1 Q0] , -1.0*t_0_q_d_1 [Q2^ Q0^ Q5 Q4] + 1.0*t_0_q_d_1 [Q5^ Q4^ Q2 Q0] , -1.0*t_0_q_d_2 [Q2^ Q1^ Q5 Q4] + 1.0*t_0_q_d_2 [Q5^ Q4^ Q2 Q1] , -1.0*t_0_q_d_3 [Q3^ Q0^ Q5 Q4] + 1.0*t_0_q_d_3 [Q5^ Q4^ Q3 Q0] , -1.0*t_0_q_d_4 [Q3^ Q1^ Q5 Q4] + 1.0*t_0_q_d_4 [Q5^ Q4^ Q3 Q1] , -1.0*t_0_q_d_5 [Q3^ Q2^ Q5 Q4] + 1.0*t_0_q_d_5 [Q5^ Q4^ Q3 Q2] ]
- class mindquantum.algorithm.nisq.Transform(operator, n_qubits=None)[source]
Class for transforms of fermionic and qubit operators. Methods jordan_wigner, parity, bravyi_kitaev, bravyi_kitaev_tree, bravyi_kitaev_superfast make transform of fermionic operators to qubit ones, they are initialized by FermionOperator, return QubitOperator. Note method reversed_jordan_wigner makes transform of qubit operator to fermionic one, it is initialized by QubitOperator, returns FermionOperator.
- Parameters
operator (Union[FermionOperator, QubitOperator]) – The input FermionOperator or QubitOperator that need to do transform.
n_qubits (int) – The total qubits of this operator. Default: None
Examples
>>> from mindquantum.core.operators import FermionOperator >>> op1 = FermionOperator('1^') >>> op1 1.0 [1^] >>> from mindquantum.algorithm.nisq.chem import Transform >>> op_transform = Transform(op1) >>> op_transform.jordan_wigner() 0.5 [Z0 X1] + -0.5j [Z0 Y1] >>> op_transform.parity() 0.5 [Z0 X1] + -0.5j [Y1] >>> op_transform.bravyi_kitaev() 0.5 [Z0 X1] + -0.5j [Y1] >>> op_transform.ternary_tree() 0.5 [X0 Z1] + -0.5j [Y0 X2] >>> op2 = FermionOperator('1^', 'a') >>> Transform(op2).jordan_wigner() 0.5*a [Z0 X1] + -0.5*I*a [Z0 Y1]
- bravyi_kitaev()[source]
Apply Bravyi-Kitaev transform.
The Bravyi-Kitaev basis is a middle between Jordan-Wigner and parity transform. That is, it balances the locality of occupation and parity information for improved simulation efficiency. In this scheme, qubits store the parity of a set of \(2^x\) orbitals, where \(x \ge 0\). A qubit of index j always stores orbital \(j\). For even values of \(j\), this is the only orbital that it stores, but for odd values of \(j\), it also stores a certain set of adjacent orbitals with index less than \(j\). For the occupation transformation, we follow the formular:
\[b_{i} = \sum{[\beta_{n}]_{i,j}} f_{j},\]where \(\beta_{n}\) is the \(N\times N\) square matrix, \(N\) is the total qubit number. The qubits index are divide into three sets, the parity set, the update set and flip set. The parity of this set of qubits has the same parity as the set of orbitals with index less than \(j\), and so we will call this set of qubit indices the “parity set” of index \(j\), or \(P(j)\).
the update set of index \(j\), or \(U(j)\) contains the set of qubits (other than qubit \(j\)) that must be updated when the occupation of orbital \(j\) This is the set of qubits in the Bravyi-Kitaev basis that store a partial sum including orbital \(j\). the flip set of index \(j\), or \(F(j)\) contains the set of BravyiKitaev qubits determines whether qubit \(j\) has the same parity or inverted parity with respect to orbital \(j\).
Please see some detail explanation in the paper (THE JOURNAL OF CHEMICAL PHYSICS 137, 224109 (2012)).
Implementation from https://arxiv.org/pdf/quant-ph/0003137.pdf and “A New Data Structure for Cumulative Frequency Tables” by Peter M. Fenwick.
- Returns
QubitOperator, qubit operator after bravyi_kitaev transformation.
- bravyi_kitaev_superfast()[source]
Apply Bravyi-Kitaev Superfast transform. Implementation from https://arxiv.org/pdf/1712.00446.pdf
Note that only hermitian operators of form
\[C + \sum_{p, q} h_{p, q} a^\dagger_p a_q + \sum_{p, q, r, s} h_{p, q, r, s} a^\dagger_p a^\dagger_q a_r a_s\]where \(C\) is a constant, be transformed.
- Returns
QubitOperator, qubit operator after bravyi_kitaev_superfast.
- jordan_wigner()[source]
Apply Jordan-Wigner transform. The Jordan-Wigner transform holds the initial occupation number locally. which change the formular of fermion operator into qubit operator following the equation.
\[ \begin{align}\begin{aligned}a^\dagger_{j}\rightarrow \sigma^{-}_{j} X \prod_{i=0}^{j-1}\sigma^{Z}_{i}\\a_{j}\rightarrow \sigma^{+}_{j} X \prod_{i=0}^{j-1}\sigma^{Z}_{i},\end{aligned}\end{align} \]where the \(\sigma_{+}= \sigma^{X} + i \sigma^{Y}\) and \(\sigma_{-} = \sigma^{X} - i\sigma^{Y}\) is the Pauli spin raising and lowring operator.
- Returns
QubitOperator, qubit operator after jordan_wigner transformation.
- parity()[source]
Apply parity transform. The parity transform stores the initial occupation number nonlocally. with the formular:
\[\left|f_{M−1}, f_{M−2},\cdots, f_0\right> → \left|q_{M−1}, q_{M−2},\cdots, q_0\right>,\]where
\[q_{m} = \left|\left(\sum_{i=0}^{m-1}f_{i}\right) mod\ 2 \right>\]Basically, this formular could be written as this,
\[p_{i} = \sum{[\pi_{n}]_{i,j}} f_{j},\]where \(\pi_{n}\) is the \(N\times N\) square matrix, \(N\) is the total qubit number. The operator changes follows the following equation as:
\[ \begin{align}\begin{aligned}a^\dagger_{j}\rightarrow\frac{1}{2}\left(\prod_{i=j+1}^N \left(\sigma_i^X X\right)\right)\left( \sigma^{X}_{j}-i\sigma_j^Y\right) X \sigma^{Z}_{j-1}\\a_{j}\rightarrow\frac{1}{2}\left(\prod_{i=j+1}^N \left(\sigma_i^X X\right)\right)\left( \sigma^{X}_{j}+i\sigma_j^Y\right) X \sigma^{Z}_{j-1}\end{aligned}\end{align} \]- Returns
QubitOperator, qubits operator after parity transformation.
- reversed_jordan_wigner()[source]
Apply reversed Jordan-Wigner transform.
- Returns
FermionOperator, fermion operator after reversed_jordan_wigner transformation.
- ternary_tree()[source]
Apply Ternary tree transform. Implementation from https://arxiv.org/pdf/1910.10746.pdf.
- Returns
QubitOperator, qubit operator after ternary_tree transformation.
- class mindquantum.algorithm.nisq.UCCAnsatz(n_qubits=None, n_electrons=None, occ_orb=None, vir_orb=None, generalized=False, trotter_step=1)[source]
The unitary coupled-cluster ansatz for molecular simulations.
\[U(\vec{\theta}) = \prod_{j=1}^{N(N\ge1)}{\prod_{i=0}^{N_{j}}{\exp{(\theta_{i}\hat{\tau}_{i})}}}\]where \(\hat{\tau}\) are anti-Hermitian operators.
Note
Currently, the circuit is construncted using JW transformation. In addition, the reference state wave function (Hartree-Fock) will NOT be included.
- Parameters
n_qubits (int) – Number of qubits (spin-orbitals). Default: None.
n_electrons (int) – Number of electrons (occupied spin-orbitals). Default: None.
occ_orb (list) – Indices of manually assigned occupied spatial orbitals, for ansatz construction only. Default: None.
vir_orb (list) – Indices of manually assigned virtual spatial orbitals, for ansatz construction only. Default: None.
generalized (bool) – Whether to use generalized excitations which do not distinguish occupied or virtual orbitals (UCCGSD). Default: False.
trotter_step (int) – The order of Trotterization step. Default: 1.
Examples
>>> from mindquantum.algorithm.nisq.chem import UCCAnsatz >>> ucc = UCCAnsatz(12, 4, occ_orb=[1], ... vir_orb=[2, 3], ... generalized=True, ... trotter_step=2) >>> circuit = ucc.circuit.remove_barrier() >>> len(circuit) 3624 >>> params_list = ucc.circuit.params_name >>> len(params_list) 48 >>> circuit[-10:] q5: ──●────RX(7π/2)───────H───────●────────────────────────────●───────H────── │ │ │ q7: ──X───────H────────RX(π/2)────X────RZ(-0.5*t_1_d0_d_17)────X────RX(7π/2)──
- mindquantum.algorithm.nisq.generate_uccsd(molecular, th=0)[source]
Generate a uccsd quantum circuit based on a molecular data generated by HiQfermion or openfermion.
- Parameters
- Returns
uccsd_circuit (Circuit), the ansatz circuit generated by uccsd method.
initial_amplitudes (numpy.ndarray), the initial parameter values of uccsd circuit.
parameters_name (list[str]), the name of initial parameters.
qubit_hamiltonian (QubitOperator), the hamiltonian of the molecule.
n_qubits (int), the number of qubits in simulation.
n_electrons, the number of electrons of the molecule.
- mindquantum.algorithm.nisq.get_qubit_hamiltonian(mol)[source]
Get the qubit hamiltonian of a molecular data.
- Parameters
mol (MolecularData) – molecular data.
- Returns
QubitOperator, qubit operator of this molecular.
- mindquantum.algorithm.nisq.quccsd_generator(n_qubits=None, n_electrons=None, anti_hermitian=True, occ_orb=None, vir_orb=None, generalized=False)[source]
Generate qubit-UCCSD (qUCCSD) ansatz using qubit-excitation operators.
Note
Currently, unrestricted version is implemented, i.e., excitations from the same spatial-orbital but with different spins will use distinct variational parameters.
- Parameters
n_qubits (int) – Number of qubits (spin-orbitals). Default: None.
n_electrons (int) – Number of electrons (occupied spin-orbitals). Default: None.
anti_hermitian (bool) – Whether to subtract the hermitian conjugate to form anti-Hermitian operators. Default: True.
occ_orb (list) – Indices of manually assigned occupied spatial orbitals. Default: None.
vir_orb (list) – Indices of manually assigned virtual spatial orbitals. Default: None.
generalized (bool) – Whether to use generalized excitations which do not distinguish occupied or virtual orbitals (qUCCGSD). Default: False.
- Returns
Generator of the qUCCSD operators.
- Return type
Examples
>>> from mindquantum.algorithm.nisq.chem import quccsd_generator >>> quccsd_generator() 0 >>> quccsd_generator(4, 2) -1.0*q_s_0 [Q0^ Q2] + -1.0*q_s_2 [Q0^ Q3] + -1.0*q_d_0 [Q1^ Q0^ Q3 Q2] + -1.0*q_s_1 [Q1^ Q2] + -1.0*q_s_3 [Q1^ Q3] + 1.0*q_s_0 [Q2^ Q0] + 1.0*q_s_1 [Q2^ Q1] + 1.0*q_s_2 [Q3^ Q0] + 1.0*q_s_3 [Q3^ Q1] + 1.0*q_d_0 [Q3^ Q2^ Q1 Q0] >>> q_op = quccsd_generator(occ_orb=[0], vir_orb=[1], generalized=True) >>> q_qubit_op = q_op.to_qubit_operator() >>> print(str(q_qubit_op)[:315]) 0.125*I*q_d_4 + 0.125*I*q_d_7 + 0.125*I*q_d_9 [X0 X1 X2 Y3] + 0.125*I*q_d_4 - 0.125*I*q_d_7 - 0.125*I*q_d_9 [X0 X1 Y2 X3] + 0.25*I*q_d_12 + 0.25*I*q_d_5 + 0.5*I*q_s_0 - 0.5*I*q_s_3 [X0 Y1] + -0.125*I*q_d_4 + 0.125*I*q_d_7 - 0.125*I*q_d_9 [X0 Y1 X2 X3] + 0.125*I*q_d_4 + 0.125*I*q_d_7 - 0.125*I*q_d_9 [X0 Y1 Y2 Y3] +
- mindquantum.algorithm.nisq.uccsd0_singlet_generator(n_qubits=None, n_electrons=None, anti_hermitian=True, occ_orb=None, vir_orb=None, generalized=False)[source]
Generate UCCSD operators using CCD0 ansatz for molecular systems.
Note
Manually assigned occ_orb or vir_orb are indices of spatial orbitals instead of spin-orbitals. They will override n_electrons and n_qubits. This is to some degree similar to the active space, therefore can reduce the number of variational parameters. However, it may not reduce the number of required qubits, since Fermion excitation operators are non-local, i.e., \(a_{7}^{\dagger} a_{0}\) involves not only the 0th and 7th qubit, but also the 1st, 2nd, … 6th qubit.
- Parameters
n_qubits (int) – Number of qubits (spin-orbitals). Default: None.
n_electrons (int) – Number of electrons (occupied spin-orbitals). Default: None.
anti_hermitian (bool) – Whether to subtract the hermitian conjugate to form anti-Hermitian operators. Default: True.
occ_orb (list) – Indices of manually assigned occupied spatial orbitals. Default: None.
vir_orb (list) – Indices of manually assigned virtual spatial orbitals. Default: None.
generalized (bool) – Whether to use generalized excitations which do not distinguish occupied or virtual orbitals (UCCGSD). Default: False.
- Returns
FermionOperator, Generator of the UCCSD operators that uses CCD0 ansatz.
Examples
>>> from mindquantum.algorithm.nisq.chem.uccsd0 import uccsd0_singlet_generator >>> uccsd0_singlet_generator(4, 2) -1.0*d0_s_0 [0^ 2] + 2.0*d0_d_0 [1^ 0^ 3 2] + -1.0*d0_s_0 [1^ 3] + 1.0*d0_s_0 [2^ 0] + 1.0*d0_s_0 [3^ 1] + -2.0*d0_d_0 [3^ 2^ 1 0] >>> uccsd0_singlet_generator(4, 2, generalized=True) 1.0*d0_s_0 - 1.0*d0_s_1 [0^ 2] + 1.0*d0_d_0 [1^ 0^ 2 1] + -1.0*d0_d_0 [1^ 0^ 3 0] + -2.0*d0_d_1 [1^ 0^ 3 2] + 1.0*d0_s_0 - 1.0*d0_s_1 [1^ 3] + -1.0*d0_s_0 + 1.0*d0_s_1 [2^ 0] + -1.0*d0_d_0 [2^ 1^ 1 0] + 1.0*d0_d_2 [2^ 1^ 3 2] + 1.0*d0_d_0 [3^ 0^ 1 0] + -1.0*d0_d_2 [3^ 0^ 3 2] + -1.0*d0_s_0 + 1.0*d0_s_1 [3^ 1] + 2.0*d0_d_1 [3^ 2^ 1 0] + -1.0*d0_d_2 [3^ 2^ 2 1] + 1.0*d0_d_2 [3^ 2^ 3 0] >>> uccsd0_singlet_generator(6, 2, occ_orb=[0], vir_orb=[1]) -1.0*d0_s_0 [0^ 2] + 2.0*d0_d_0 [1^ 0^ 3 2] + -1.0*d0_s_0 [1^ 3] + 1.0*d0_s_0 [2^ 0] + 1.0*d0_s_0 [3^ 1] + -2.0*d0_d_0 [3^ 2^ 1 0]
- mindquantum.algorithm.nisq.uccsd_singlet_generator(n_qubits, n_electrons, anti_hermitian=True)[source]
Create a singlet UCCSD generator for a system with n_electrons
This function generates a FermionOperator for a UCCSD generator designed to act on a single reference state consisting of n_qubits spin orbitals and n_electrons electrons, that is a spin singlet operator, meaning it conserves spin.
- Parameters
n_qubits (int) – Number of spin-orbitals used to represent the system, which also corresponds to number of qubits in a non-compact map.
n_electrons (int) – Number of electrons in the physical system.
anti_hermitian (bool) – Flag to generate only normal CCSD operator rather than unitary variant, primarily for testing
- Returns
FermionOperator, Generator of the UCCSD operator that builds the UCCSD wavefunction.
Examples
>>> from mindquantum.algorithm.nisq.chem import uccsd_singlet_generator >>> uccsd_singlet_generator(4, 2) -s_0 [0^ 2] + -d1_0 [0^ 2 1^ 3] + -s_0 [1^ 3] + -d1_0 [1^ 3 0^ 2] + s_0 [2^ 0] + d1_0 [2^ 0 3^ 1] + s_0 [3^ 1] + d1_0 [3^ 1 2^ 0]
- mindquantum.algorithm.nisq.uccsd_singlet_get_packed_amplitudes(single_amplitudes, double_amplitudes, n_qubits, n_electrons)[source]
Convert amplitudes for use with singlet UCCSD
The output list contains only those amplitudes that are relevant to singlet UCCSD, in an order suitable for use with the function uccsd_singlet_generator.
- Parameters
single_amplitudes (numpy.ndarray) – \(N\times N\) array storing single excitation amplitudes corresponding to \(t_{i,j} * (a_i^\dagger a_j - \text{H.C.})\)
double_amplitudes (numpy.ndarray) – \(N\times N\times N\times N\) array storing double excitation amplitudes corresponding to \(t_{i,j,k,l} * (a_i^\dagger a_j a_k^\dagger a_l - \text{H.C.})\)
n_qubits (int) – Number of spin-orbitals used to represent the system, which also corresponds to number of qubits in a non-compact map.
n_electrons (int) – Number of electrons in the physical system.
- Returns
ParameterResolver, List storing the unique single and double excitation amplitudes for a singlet UCCSD operator. The ordering lists unique single excitations before double excitations.
Examples
>>> import numpy as np >>> from mindquantum.algorithm.nisq.chem import uccsd_singlet_get_packed_amplitudes >>> n_qubits, n_electrons = 4, 2 >>> np.random.seed(42) >>> ccsd_single_amps = np.random.random((4, 4)) >>> ccsd_double_amps = np.random.random((4, 4, 4, 4)) >>> uccsd_singlet_get_packed_amplitudes(ccsd_single_amps, ccsd_double_amps, ... n_qubits, n_electrons) {'s_0': 0.6011150117432088, 'd1_0': 0.7616196153287176}