mindquantum.simulator.mqchem.MQChemSimulator

class mindquantum.simulator.mqchem.MQChemSimulator(n_qubits, n_electrons, seed=None, dtype='double')[源代码]

基于组态相互作用(CI)方法的量子化学模拟器。

该模拟器通过在CI向量空间(由固定电子数定义的完整希尔伯特空间的子空间)中工作,针对化学问题进行了优化。 对于典型的化学模拟,这种方法比使用完整的态矢量模拟器具有显著的性能优势。

该模拟器设计用于与 UCCExcitationGateCIHamiltonian 协同工作。它提供了应用UCC线路、 计算哈密顿量期望值以及计算变分量子算法(如VQE)所需梯度的方法。

默认情况下,模拟器在Hartree-Fock态下初始化,该状态是量子化学计算的典型参考态。

参数:
  • n_qubits (int) - 系统中的总量子比特数(自旋轨道数)。

  • n_electrons (int) - 电子数,定义了CI空间的维度。

  • seed (int) - 此模拟器的随机种子。如果为 None,将生成一个随机种子。默认值:None

  • dtype (str) - 模拟的数据类型,可以是 "float""double"。默认值:"double"

样例:

>>> 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, pr=None)[源代码]

将量子线路应用于当前模拟器状态。

说明

线路中只有 UCCExcitationGate 的实例会被应用;所有其他类型的门都将被忽略。

参数:
  • circuit (Union[Circuit, Iterable[UCCExcitationGate]]) - 要应用的量子线路或UCC门的可迭代对象。

  • pr (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) - 用于替换参数值的参数解析器。如果为 None,则直接使用门中定义的参数。默认值:None

样例:

>>> 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, pr=None)[源代码]

将单个UCC激发门应用于当前模拟器状态。

参数:
  • gate (UCCExcitationGate) - 要应用的UCC激发门。

  • pr (Union[ParameterResolver, dict, numpy.ndarray, list, numbers.Number]) - 用于替换参数值的参数解析器。如果为 None,则直接使用门中定义的参数。默认值:None

异常:

样例:

>>> 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)[源代码]

计算哈密顿量相对于当前状态的期望值。

计算 \(\langle\psi|H|\psi\rangle\),其中 \(|\psi\rangle\) 是当前的CI态矢量,\(H\) 是CI哈密顿量。

参数:
  • ham (CIHamiltonian) - 要计算期望值的哈密顿量。

返回:

float,实数期望值。

异常:

样例:

>>> 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, circuit)[源代码]

生成一个计算期望值及其梯度的函数。

该方法实现了伴随微分法,以计算期望值 \(\langle\psi(\theta)|H|\psi(\theta)\rangle\) 相对于UCC ansatz线路参数 \(\theta\) 的梯度。 状态被制备为 \(|\psi(\theta)\rangle = U(\theta)|\psi_0\rangle\),其中 \(|\psi_0\rangle\) 是模拟器的当前状态。

参数:
  • ham (CIHamiltonian) - 哈密顿量 \(H\)

  • circuit (Union[Circuit, Iterable[UCCExcitationGate]]) - 参数化的UCC线路 \(U(\theta)\)。该线路必须具有用于梯度计算的参数。

返回:

Callable,一个接受参数值NumPy数组 x 并返回元组 (expectation, gradient) 的函数。expectation 是浮点数期望值,gradient 是一个NumPy数组,包含相对于 x 中每个参数的导数。参数的顺序由 circuit.params_name 决定。

异常:

样例:

>>> 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=False)[源代码]

获取模拟器当前的量子态。

虽然模拟器内部将状态存储为紧凑的CI向量,但此方法以完整的 \(2^{n_{qubits}}\) 维计算基中的密集态矢量形式返回状态。

参数:
  • ket (bool) - 如果为 True,则以狄拉克(ket)符号的字符串形式返回量子态。如果为 False,则以NumPy数组形式返回状态。默认值:False

返回:

Union[numpy.ndarray, str],量子态矢量,以NumPy数组或ket符号字符串形式表示。

异常:
  • TypeError - 如果 ket 不是布尔值。

样例:

>>> 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()[源代码]

将模拟器的状态重置为Hartree-Fock(HF)态。

Hartree-Fock态是无相互作用费米子系统的基态,其中 n_electrons 个最低能量的自旋轨道被占据。在计算基中,这对应于状态 \(|11...100...0\rangle\)

样例:

>>> 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)[源代码]

从稀疏表示设置CI向量。

参数:
  • qs_rep (List[Tuple[int, complex]]) - 一个元组列表,其中每个元组 (mask, amplitude) 定义了基态的振幅。mask 是表示计算基态(斯莱特行列式)的整数,amplitude 是其对应的复振幅。qs_rep 中的所有基态的布居数必须等于 n_electrons

样例:

>>> 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⟩