mindquantum.simulator.mqchem.MQChemSimulator
- class mindquantum.simulator.mqchem.MQChemSimulator(n_qubits, n_electrons, seed=None, dtype='double')[source]
A simulator for quantum chemistry based on the Configuration Interaction (CI) method.
This simulator is optimized for chemistry problems by working within a CI vector space, which is a subspace of the full Hilbert space defined by a fixed number of electrons. This approach is significantly more memory-efficient than using a full state-vector simulator for typical chemistry simulations.
The simulator is designed to work with
UCCExcitationGate
andCIHamiltonian
. It provides methods to apply UCC circuits, calculate Hamiltonian expectation values, and compute gradients required for variational quantum algorithms like VQE.By default, the simulator is initialized in the Hartree-Fock state, which serves as the typical reference state for quantum chemistry calculations.
- Parameters
n_qubits (int) – The total number of qubits (spin-orbitals) in the system.
n_electrons (int) – The number of electrons, which defines the dimension of the CI space.
seed (int) – The random seed for this simulator. If
None
, a random seed will be generated. Default:None
.dtype (str) – The data type for simulation, either
"float"
or"double"
. Default:"double"
.
Examples
>>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2, seed=42) >>> sim.n_qubits 4 >>> sim.n_electrons 2 >>> # The simulator starts in the Hartree-Fock state |0011> >>> print(sim.get_qs(ket=True)) 1¦0011⟩
- apply_circuit(circuit: Union[Circuit, Iterable[UCCExcitationGate]], pr: ParameterResolver = None)[source]
Apply a quantum circuit to the current simulator state.
Note
Only
UCCExcitationGate
instances within the circuit will be applied; all other gate types are ignored.- Parameters
circuit (Union[Circuit, Iterable[UCCExcitationGate]]) – The quantum circuit or an iterable of UCC gates to apply.
pr (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) – A parameter resolver to substitute parameter values. If
None
, the parameters defined in the gates are used directly. Default:None
.
Examples
>>> import numpy as np >>> from mindquantum.core.operators import FermionOperator >>> from mindquantum.core.circuit import Circuit >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> t1 = FermionOperator('2^ 0', 'a') >>> t2 = FermionOperator('3^ 1', 'b') >>> circ = Circuit([mqchem.UCCExcitationGate(t1), mqchem.UCCExcitationGate(t2)]) >>> sim.apply_circuit(circ, {'a': 0.1, 'b': 0.2}) >>> print(sim.get_qs(ket=True)) 0.97517033¦0011⟩ -0.0978434¦0110⟩ 0.19767681¦1001⟩ 0.01983384¦1100⟩
- apply_gate(gate: UCCExcitationGate, pr: ParameterResolver = None)[source]
Apply a single UCC excitation gate to the current simulator state.
- Parameters
gate (UCCExcitationGate) – The UCC excitation gate to apply.
pr (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) – A parameter resolver to substitute parameter values. If
None
, the parameter defined in the gate is used directly. Default:None
.
- Raises
TypeError – If gate is not a
UCCExcitationGate
.
Examples
>>> import numpy as np >>> from mindquantum.core.operators import FermionOperator >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> t = FermionOperator('2^ 0', 'theta') >>> gate = mqchem.UCCExcitationGate(t) >>> # Apply gate with theta = 0.1 >>> sim.apply_gate(gate, {'theta': 0.1}) >>> print(sim.get_qs(ket=True)) 0.99500417¦0011⟩ -0.09983342¦0110⟩
- get_expectation(ham)[source]
Compute the expectation value of a Hamiltonian with respect to the current state.
Calculates \(\langle\psi|H|\psi\rangle\), where \(|\psi\rangle\) is the current CI state vector and \(H\) is a CI Hamiltonian.
- Parameters
ham (CIHamiltonian) – The Hamiltonian for which to compute the expectation.
- Returns
float, The real-valued expectation value.
- Raises
TypeError – If ham is not a
CIHamiltonian
.
Examples
>>> from mindquantum.core.operators import FermionOperator >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> # Hartree-Fock state is |0011> >>> ham_op = FermionOperator('0^ 0') + FermionOperator('1^ 1') >>> ham = mqchem.CIHamiltonian(ham_op) >>> # Expectation of <HF| (n_0 + n_1) |HF> should be 1+1=2 >>> sim.get_expectation(ham) 2.0
- get_expectation_with_grad(ham: CIHamiltonian, circuit: Union[Circuit, Iterable[UCCExcitationGate]])[source]
Generate a function to compute the expectation value and its gradient.
This method implements the adjoint differentiation method to calculate the gradient of the expectation value \(\langle\psi(\theta)|H|\psi(\theta)\rangle\) with respect to the parameters \(\theta\) of a UCC ansatz circuit. The state is prepared as \(|\psi(\theta)\rangle = U(\theta)|\psi_0\rangle\), where \(|\psi_0\rangle\) is the current state of the simulator.
- Parameters
ham (CIHamiltonian) – The Hamiltonian \(H\).
circuit (Union[Circuit, Iterable[UCCExcitationGate]]) – The parameterized UCC circuit \(U(\theta)\). The circuit must have parameters for gradient calculation.
- Returns
Callable, A function that accepts a NumPy array of parameter values x and returns a tuple (expectation, gradient). expectation is the float expectation value, and gradient is a NumPy array containing the derivatives with respect to each parameter in x. The order of parameters is determined by circuit.params_name.
- Raises
TypeError – If ham is not a
CIHamiltonian
.
Examples
>>> import numpy as np >>> from mindquantum.core.operators import FermionOperator >>> from mindquantum.core.circuit import Circuit >>> from mindquantum.simulator import mqchem >>> >>> # Prepare simulator, Hamiltonian, and ansatz >>> sim = mqchem.MQChemSimulator(4, 2) >>> ham_op = FermionOperator('0^ 0', 1) + FermionOperator('1^ 1', 1) >>> ham = mqchem.CIHamiltonian(ham_op) >>> t = FermionOperator('2^ 0', 'theta') >>> ansatz = Circuit(mqchem.UCCExcitationGate(t)) >>> >>> # Get the gradient function >>> grad_fn = sim.get_expectation_with_grad(ham, ansatz) >>> >>> # Calculate expectation and gradient for theta = 0.1 >>> params = np.array([0.1]) >>> expectation, gradient = grad_fn(params) >>> print(f"Expectation: {expectation:.5f}") Expectation: 1.99003+0.00000j >>> print(f"Gradient: {gradient[0]:.5f}") Gradient: -0.19867+0.00000j
- get_qs(ket: bool = False)[source]
Get the current quantum state of the simulator.
While the simulator internally stores the state as a compact CI vector, this method returns the state as a dense state vector in the full \(2^{n_{qubits}}\) dimensional computational basis.
- Parameters
ket (bool) – If
True
, returns the quantum state in Dirac (ket) notation as a string. IfFalse
, returns the state as a NumPy array. Default:False
.- Returns
Union[numpy.ndarray, str], The quantum state vector as a NumPy array or a string in ket notation.
- Raises
TypeError – If ket is not a boolean.
Examples
>>> import numpy as np >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> # Get state vector (Hartree-Fock |0011>) >>> np.round(sim.get_qs(), 2) array([0.+0.j, 0.+0.j, 0.+0.j, 1.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j, 0.+0.j]) >>> # Get state in ket string format >>> print(sim.get_qs(ket=True)) 1¦0011⟩
- reset()[source]
Reset the simulator's state to the Hartree-Fock (HF) state.
The Hartree-Fock state is the ground state of a non-interacting fermionic system, where the n_electrons lowest-energy spin-orbitals are occupied. In the computational basis, this corresponds to the state \(|11...100...0⟩\).
Examples
>>> from mindquantum.core.operators import FermionOperator >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> t = FermionOperator('2^ 0', 0.1) >>> ucc_gate = mqchem.UCCExcitationGate(t) >>> sim.apply_gate(ucc_gate) # State is no longer HF state >>> sim.reset() >>> print(sim.get_qs(ket=True)) 1¦0011⟩
- set_qs(qs_rep)[source]
Set the CI vector from a sparse representation.
- Parameters
qs_rep (List[Tuple[int, complex]]) – A list of tuples, where each tuple (mask, amplitude) defines the amplitude of a basis state. mask is an integer representing the computational basis state (Slater determinant), and amplitude is its corresponding complex amplitude. All basis states in qs_rep must have a population count equal to n_electrons.
Examples
>>> import numpy as np >>> from mindquantum.simulator import mqchem >>> sim = mqchem.MQChemSimulator(4, 2) >>> # Prepare a state that is a superposition of |0101> and |1001> >>> # Note: masks 5 (0101) and 9 (1001) both have 2 electrons. >>> new_state_rep = [(5, 1/np.sqrt(2)), (9, 1/np.sqrt(2))] >>> sim.set_qs(new_state_rep) >>> print(sim.get_qs(ket=True)) √2/2¦0101⟩ √2/2¦1001⟩