mindquantum.simulator.Simulator
- class mindquantum.simulator.Simulator(backend, n_qubits=None, *args, seed=None, dtype=None, **kwargs)[source]
Quantum simulator that simulate quantum circuit.
- Parameters
backend (str) – which backend you want. The supported backend can be found in SUPPORTED_SIMULATOR
n_qubits (int) – number of quantum simulator. Default:
None
.seed (int) – the random seed for this simulator, if
None
, seed will generate by numpy.random.randint. Default:None
.dtype (mindquantum.dtype) – the data type of simulator. Default:
None
.
- Raises
TypeError – if backend is not str.
TypeError – if n_qubits is not int.
TypeError – if seed is not int.
ValueError – if backend is not supported.
ValueError – if n_qubits is negative.
ValueError – if seed is less than 0 or great than \(2^23 - 1\).
Examples
>>> from mindquantum.algorithm.library import qft >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 2) >>> sim.apply_circuit(qft(range(2))) >>> sim.get_qs() array([0.5+0.j, 0.5+0.j, 0.5+0.j, 0.5+0.j])
- apply_circuit(circuit, pr=None)[source]
Apply a circuit on this simulator.
- Parameters
circuit (Circuit) – The quantum circuit you want to apply on this simulator.
pr (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) – The parameter resolver for this circuit. If the circuit is not parameterized, this arg should be
None
. Default:None
.
- Returns
MeasureResult or None, if the circuit has measure gate, then return a MeasureResult, otherwise return None.
Examples
>>> import numpy as np >>> from mindquantum.core.circuit import Circuit >>> from mindquantum.core.gates import H >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 2) >>> sim.apply_circuit(Circuit().un(H, 2)) >>> sim.apply_circuit(Circuit().ry('a', 0).ry('b', 1), np.array([1.1, 2.2])) >>> sim mqvector simulator with 2 qubits (little endian). Current quantum state: -0.0721702531972066¦00⟩ -0.30090405886869676¦01⟩ 0.22178317006196263¦10⟩ 0.9246947752567126¦11⟩ >>> sim.apply_circuit(Circuit().measure(0).measure(1)) shots: 1 Keys: q1 q0│0.00 0.2 0.4 0.6 0.8 1.0 ───────────┼───────────┴───────────┴───────────┴───────────┴───────────┴ 11│▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓ │ {'11': 1}
- apply_gate(gate, pr=None, diff=False)[source]
Apply a gate on this simulator, can be a quantum gate or a measurement operator.
- Parameters
gate (BasicGate) – The gate you want to apply.
pr (Union[numbers.Number, numpy.ndarray, ParameterResolver, list]) – The parameter for parameterized gate. Default:
None
.diff (bool) – Whether to apply the derivative gate on this simulator. Default:
False
.
- Returns
int or None, if the gate if a measure gate, then return a collapsed state, Otherwise return None.
- Raises
TypeError – if gate is not a BasicGate.
ValueError – if any qubit of gate is higher than simulator qubits.
ValueError – if gate is parameterized, but no parameter supplied.
TypeError – the pr is not a ParameterResolver if gate is parameterized.
Examples
>>> import numpy as np >>> from mindquantum.core.gates import RY, Measure >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 1) >>> sim.apply_gate(RY('a').on(0), np.pi/2) >>> sim.get_qs() array([0.70710678+0.j, 0.70710678+0.j]) >>> sim.apply_gate(Measure().on(0)) 1 >>> sim.get_qs() array([0.+0.j, 1.+0.j])
- apply_hamiltonian(hamiltonian: Hamiltonian)[source]
Apply hamiltonian to a simulator, this hamiltonian can be hermitian or non hermitian.
Note
The quantum state may be not a normalized quantum state after apply hamiltonian.
- Parameters
hamiltonian (Hamiltonian) – the hamiltonian you want to apply.
Examples
>>> from mindquantum.core.circuit import Circuit >>> from mindquantum.core.operators import QubitOperator, Hamiltonian >>> from mindquantum.simulator import Simulator >>> import scipy.sparse as sp >>> sim = Simulator('mqvector', 1) >>> sim.apply_circuit(Circuit().h(0)) >>> sim.get_qs() array([0.70710678+0.j, 0.70710678+0.j]) >>> ham1 = Hamiltonian(QubitOperator('Z0')) >>> sim.apply_hamiltonian(ham1) >>> sim.get_qs() array([ 0.70710678+0.j, -0.70710678+0.j]) >>> sim.reset() >>> ham2 = Hamiltonian(sp.csr_matrix([[1, 2], [3, 4]])) >>> sim.apply_hamiltonian(ham2) >>> sim.get_qs() array([1.+0.j, 3.+0.j])
- astype(dtype, seed=None)[source]
Convert simulator to other data type.
Note
The quantum state will copied from origin simulator.
- Parameters
dtype (mindquantum.dtype) – the data type of new simulator.
seed (int) – the seed of new simulator. Default:
None
.
- copy()[source]
Copy this simulator.
- Returns
Simulator, a copy version of this simulator.
Examples
>>> from mindquantum.core.gates import RX >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 1) >>> sim.apply_gate(RX(1).on(0)) >>> sim2 = sim.copy() >>> sim2.apply_gate(RX(-1).on(0)) >>> sim2 mqvector simulator with 1 qubit (little endian). Current quantum state: 1¦0⟩
- property dtype
Get data type of simulator.
- entropy()[source]
Calculate the von Neumann entropy of current quantum state.
Definition of von Neumann entropy \(S\) shown as below.
\[S(\rho) = -\text{tr}(\rho \ln \rho)\]where \(\rho\) is density matrix.
- Returns
numbers.Number, the von Neumann entropy of current quantum state.
Examples
>>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqmatrix', 1) >>> sim.set_qs([[0.5, 0], [0, 0.5]]) >>> sim.entropy() 0.6931471805599453
- get_expectation(hamiltonian, circ_right=None, circ_left=None, simulator_left=None, pr=None)[source]
Get expectation of the given hamiltonian. The hamiltonian could be non hermitian.
This method is designed to calculate the expectation shown as below.
\[E = \left<\varphi\right|U_l^\dagger H U_r \left|\psi\right>\]where \(U_l\) is circ_left, \(U_r\) is circ_right, \(H\) is hams and \(\left|\psi\right>\) is the current quantum state of this simulator, and \(\left|\varphi\right>\) is the quantum state of simulator_left.
- Parameters
hamiltonian (Hamiltonian) – The hamiltonian you want to get expectation.
circ_right (Circuit) – The \(U_r\) circuit described above. If it is
None
, we will use empty circuit. Default:None
.circ_left (Circuit) – The \(U_l\) circuit described above. If it is
None
, then it will be the same ascirc_right
. Default:None
.simulator_left (Simulator) – The simulator that contains \(\left|\varphi\right>\). If
None
, then \(\left|\varphi\right>\) is assumed to be equals to \(\left|\psi\right>\). Default:None
.pr (Union[Dict[str, numbers.Number], ParameterResolver]) – the variable value of circuit. Default:
None
.
- Returns
numbers.Number, the expectation value.
Examples
>>> from mindquantum.core.circuit import Circuit >>> from mindquantum.core.operators import QubitOperator, Hamiltonian >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 1) >>> sim.apply_circuit(Circuit().ry(1.2, 0)) >>> ham = Hamiltonian(QubitOperator('Z0')) >>> sim.get_expectation(ham) (0.36235775447667357+0j) >>> sim.get_expectation(ham, Circuit().rx('a', 0), Circuit().ry(2.3, 0), pr={'a': 2.4}) (-0.25463350745693886+0.8507316752782879j) >>> sim1, sim2 = Simulator('mqvector', 1), Simulator('mqvector', 1) >>> sim1.apply_circuit(Circuit().ry(1.2, 0).rx(2.4, 0)) >>> sim2.apply_circuit(Circuit().ry(1.2, 0).ry(2.3, 0)) >>> sim1.apply_hamiltonian(ham) >>> from mindquantum.simulator import inner_product >>> inner_product(sim2, sim1) (-0.25463350745693886+0.8507316752782879j)
- get_expectation_with_grad(hams, circ_right, circ_left=None, simulator_left=None, parallel_worker=None, pr_shift=False)[source]
Get a function that return the forward value and gradient w.r.t circuit parameters.
This method is designed to calculate the expectation and its gradient shown as below.
\[E = \left<\varphi\right|U_l^\dagger H U_r \left|\psi\right>\]where \(U_l\) is circ_left, \(U_r\) is circ_right, \(H\) is hams and \(\left|\psi\right>\) is the current quantum state of this simulator, and \(\left|\varphi\right>\) is the quantum state of simulator_left.
- Parameters
hams (Union[
Hamiltonian
, List[Hamiltonian
]]) – AHamiltonian
or a list ofHamiltonian
that need to get expectation.circ_right (
Circuit
) – The \(U_r\) circuit described above.circ_left (
Circuit
) – The \(U_l\) circuit described above. By default, this circuit will benone
, and in this situation, \(U_l\) will be equals to \(U_r\). Default:None
.simulator_left (
Simulator
) – The simulator that contains \(\left|\varphi\right>\). IfNone
, then \(\left|\varphi\right>\) is assumed to be equals to \(\left|\psi\right>\). Default:None
.parallel_worker (int) – The parallel worker numbers. The parallel workers can handle batch in parallel threads. Default:
None
.pr_shift (bool) – Whether or not to use parameter-shift rule. Only available in "mqvector" simulator. It will be enabled automatically when circuit contains noise channel. Noted that not every gate uses the same shift value π/2, so the gradient of FSim gate and parameterized custom gate will be calculated by finite difference method with gap 0.001. Default:
False
.
- Returns
GradOpsWrapper, a grad ops wrapper than contains information to generate this grad ops.
Examples
>>> import numpy as np >>> from mindquantum.core.circuit import Circuit >>> from mindquantum.core.operators import QubitOperator, Hamiltonian >>> from mindquantum.simulator import Simulator >>> circ = Circuit().ry('a', 0) >>> ham = Hamiltonian(QubitOperator('Z0')) >>> sim = Simulator('mqvector', 1) >>> grad_ops = sim.get_expectation_with_grad(ham, circ) >>> grad_ops(np.array([1.0])) (array([[0.54030231+0.j]]), array([[[-0.84147098+0.j]]])) >>> sim1 = Simulator('mqvector', 1) >>> prep_circ = Circuit().h(0) >>> ansatz = Circuit().ry('a', 0).rz('b', 0).ry('c', 0) >>> sim1.apply_circuit(prep_circ) >>> sim2 = Simulator('mqvector', 1) >>> ham = Hamiltonian(QubitOperator("")) >>> grad_ops = sim2.get_expectation_with_grad(ham, ansatz, Circuit(), simulator_left=sim1) >>> f, g = grad_ops(np.array([7.902762e-01, 2.139225e-04, 7.795934e-01])) >>> f array([[0.99999989-7.52279618e-05j]])
- get_partial_trace(obj_qubits)[source]
Calculate the partial trace of current density matrix.
- Parameters
obj_qubits (Union[int, list[int]]) – Specific which qubits (subsystems) to trace over.
- Returns
numpy.ndarray, the partial trace of current density matrix.
Examples
>>> from mindquantum.core.circuit import Circuit >>> from mindquantum.simulator import Simulator >>> circ = Circuit().h(0).x(1, 0) >>> sim = Simulator('mqmatrix', 2) >>> sim.apply_circuit(circ) >>> mat = sim.get_partial_trace(0) >>> mat array([[0.5-0.j, 0. -0.j], [0. +0.j, 0.5-0.j]])
- get_pure_state_vector()[source]
Get state vector if current density matrix is pure.
The relation between density matrix \(\rho\) and state vector \(\left| \psi \right>\) shown as below.
\[\rho = \left| \psi \right>\!\left< \psi \right|\]Note that the state vector \(\left| \psi \right>\) may have an arbitrary global phase \(e^{i\phi}\).
- Returns
numpy.ndarray, a state vector calculated from current density matrix.
Examples
>>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqmatrix', 1) >>> sim.set_qs([[0.5, 0.5], [0.5, 0.5]]) >>> sim.get_pure_state_vector() array([0.70710678+0.j, 0.70710678+0.j])
- get_qs(ket=False)[source]
Get current quantum state of this simulator.
- Parameters
ket (bool) – Whether to return the quantum state in ket format or not. Default:
False
.- Returns
numpy.ndarray, the current quantum state.
Examples
>>> from mindquantum.algorithm.library import qft >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 2) >>> sim.apply_circuit(qft(range(2))) >>> sim.get_qs() array([0.5+0.j, 0.5+0.j, 0.5+0.j, 0.5+0.j])
- property n_qubits
Get simulator qubit.
- Returns
int, the qubit number of simulator.
- purity()[source]
Calculate the purity of current quantum state.
Definition of purity \(\gamma\) shown as below.
\[\gamma \equiv \text{tr}(\rho^2)\]where \(\rho\) is density matrix.
- Returns
numbers.Number, the purity of current quantum state.
Examples
>>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqmatrix', 1) >>> sim.set_qs([[0.5, 0], [0, 0.5]]) >>> sim.purity() 0.5
- reset()[source]
Reset simulator to zero state.
Examples
>>> from mindquantum.algorithm.library import qft >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 2) >>> sim.apply_circuit(qft(range(2))) >>> sim.reset() >>> sim.get_qs() array([1.+0.j, 0.+0.j, 0.+0.j, 0.+0.j])
- sampling(circuit, pr=None, shots=1, seed=None)[source]
Sample the measure qubit in circuit.
Sampling does not change the original quantum state of this simulator. The sampling results are represented in little-endian order by default (e.g., '01' means q1=0, q0=1). If big-endian order is needed, use MeasureResult.reverse_endian() method.
- Parameters
circuit (Circuit) – The circuit that you want to evolve and sample.
pr (Union[None, dict, ParameterResolver]) – The parameter resolver for this circuit, if this circuit is a parameterized circuit. Default:
None
.shots (int) – How many shots you want to sample this circuit. Default:
1
.seed (int) – Random seed for random sampling. If
None
, seed will be a random int number. Default:None
.
- Returns
MeasureResult, the measure result of sampling. The bit strings in the result are in little-endian order.
Examples
>>> from mindquantum.core.circuit import Circuit >>> from mindquantum.core.gates import Measure >>> from mindquantum.simulator import Simulator >>> circ = Circuit().ry('a', 0).ry('b', 1) >>> circ += Measure('q0_0').on(0) >>> circ += Measure('q0_1').on(0) >>> circ += Measure('q1').on(1) >>> sim = Simulator('mqvector', circ.n_qubits) >>> res = sim.sampling(circ, {'a': 1.1, 'b': 2.2}, shots=100, seed=42) >>> res shots: 100 Keys: q1 q0_1 q0_0│0.00 0.122 0.245 0.367 0.49 0.612 ──────────────────┼───────────┴───────────┴───────────┴───────────┴───────────┴ 000│▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒ │ 011│▒▒▒▒▒▒▒▒▒ │ 100│▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓▓ │ 111│▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒▒ │ {'000': 18, '011': 9, '100': 49, '111': 24}
- set_qs(quantum_state)[source]
Set quantum state for this simulation.
- Parameters
quantum_state (numpy.ndarray) – the quantum state that you want.
Examples
>>> import numpy as np >>> from mindquantum.simulator import Simulator >>> sim = Simulator('mqvector', 1) >>> sim.get_qs() array([1.+0.j, 0.+0.j]) >>> sim.set_qs(np.array([1, 1])) >>> sim.get_qs() array([0.70710678+0.j, 0.70710678+0.j])