E_{CCSD}>E_{FCI}$. After the computing is complete, save the result to the `molecule_file` file (`molecule_of.filename`)."
]
},
{
"cell_type": "code",
"execution_count": 4,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"H1-Li1_sto3g_singlet\n"
]
}
],
"source": [
"molecule_of.save()\n",
"molecule_file = molecule_of.filename\n",
"print(molecule_file.split('/')[-1])"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"One of the major obstacles to quantum chemistry is the volume of computation. As the system size (electron number and atomic number) increases, the time required for solving the FCI wave function and ground state energy increases by about $2^{N}$. Even for small molecules such as ethylene molecules, FCI computing is not easy. Quantum computers provide a possible solution to this problem. Research shows that quantum computers can simulate the time-dependent evolution of Hamiltonian in terms of polynomial time complexity. Compared with classical computers, quantum computers exponentially accelerate the chemical simulation on quantum processors. This tutorial introduces one of the quantum algorithms: VQE.\n",
"\n",
"## Variational Quantum Eigensolver (VQE)\n",
"\n",
"The VQE is a hybrid quantum-classical algorithm. It uses the variational principle to solve the ground state wave function. The optimization of variational parameters is carried out on the classical computer.\n",
"\n",
"### Variational Principle\n",
"\n",
"The variational principle may be expressed in the following form:\n",
"\n",
"$$\n",
"E_{0} \\le \\frac{\\langle \\Psi_{t} | \\hat{H} | \\Psi_{t} \\rangle}{\\langle \\Psi_{t} | \\Psi_{t} \\rangle}\n",
"$$\n",
"\n",
"In the preceding formula, $| \\Psi_{t} \\rangle$ indicates the probe wave function. The variational principle shows that the ground state energy obtained by any probe wave function is always greater than or equal to the real ground state energy under certain conditions. The variational principle provides a method for solving the molecular ground state Schrödinger equation. A parameterized function $f(\\theta)$ is used as an approximation of the accurate ground state wave function, and the accurate ground state energy is approximated by optimizing the parameter $\\theta$.\n",
"\n",
"### Initial State Preparation\n",
"\n",
"The $N$-electron HF wave function also has a very concise form under the quadratic quantization expression:\n",
"\n",
"$$\n",
"| \\Psi_{HF} \\rangle = \\prod^{i=0} _{N-1}{a^{\\dagger} _{i}| 0 \\rangle}\n",
"$$\n",
"\n",
"The above formula builds a bridge from quantum chemical wave function to quantum computing: $|0\\rangle$ is used to represent a non-occupied orbit, and $|1\\rangle$ is used to represent an orbit occupied by an electron. Therefore, the $N$-electron HF wave function may be mapped to a string of $M+N$ quantum bits $| 00\\dots 11\\dots \\rangle$. $M$ indicates the number of unoccupied tracks.\n",
"\n",
"The following code builds an HF initial state wave function corresponding to the LiH molecule. In Jordan-Wigner transformation, $N$ $\\text{X}$ gates are applied to $|000\\dots\\rangle$."
]
},
{
"cell_type": "code",
"execution_count": 5,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
" ┏━━━┓ \n",
"q0: ──┨╺╋╸┠───\n",
" ┗━━━┛ \n",
" ┏━━━┓ \n",
"q1: ──┨╺╋╸┠───\n",
" ┗━━━┛ \n",
" ┏━━━┓ \n",
"q2: ──┨╺╋╸┠───\n",
" ┗━━━┛ \n",
" ┏━━━┓ \n",
"q3: ──┨╺╋╸┠───\n",
" ┗━━━┛ \n"
]
}
],
"source": [
"hartreefock_wfn_circuit = Circuit([X.on(i) for i in range(molecule_of.n_electrons)])\n",
"print(hartreefock_wfn_circuit)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"We can build a probe wave function in the following form:\n",
"\n",
"$$\n",
"| \\Psi_{t} \\rangle = U(\\theta) | \\Psi_{HF} \\rangle\n",
"$$\n",
"\n",
"$U(\\theta)$ represents a unitary transformation that may be simulated by using a quantum circuit. $| \\Psi_{HF} \\rangle$ is used as an initial state, and may be conveniently prepared by using a plurality of single-bit $\\text{X}$ gates. A specific form of the $U(\\theta) | \\Psi_{HF} \\rangle$ is also referred to as wave function ansatz.\n",
"\n",
"### Wave Function Ansatz\n",
"\n",
"The coupled-cluster theory mentioned above is a very efficient wave function ansatz. To use it on a quantum computer, you need to make the following modifications:\n",
"\n",
"$$\n",
"| \\Psi_{UCC} \\rangle = \\exp{(\\hat{T} - \\hat{T}^{\\dagger})} | \\Psi_{HF} \\rangle\n",
"$$\n",
"\n",
"UCC is short for unitary coupled-cluster theory. $\\hat{T}^{\\dagger}$ represents the Hermite conjugate of $\\hat{T}$. In this way, $\\exp{(\\hat{T} - \\hat{T}^{\\dagger})}$ is the unitary operator. [Peruzzo et al.](https://doi.org/10.1038/ncomms5213) first performed chemical simulation experiments on quantum computers using VQE and unitary coupled-cluster with singles and doubles (UCCSD) in 2014. It should be noted that, by default, the parameter $\\{\\theta\\}$ in the coupled-cluster operator is a real number. There is no problem with this hypothesis in molecular systems. In periodic systems, the study of [Liu Jie et al.](https://doi.org/10.1021/acs.jctc.0c00881) suggests that a unitary coupled-cluster can result in errors due to the neglect of the complex numbers. This tutorial does not discuss the application of unitary coupled-cluster in periodic systems.\n",
"\n",
"The `generate_uccsd` function in the circuit module of MindSpore Quantum can be used to read the computing result saved in `molecule_file`, build the UCCSD wave function by one click, and obtain the corresponding quantum circuit."
]
},
{
"cell_type": "code",
"execution_count": 6,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"ccsd:-7.882352909152679.\n",
"fci:-7.882362286798728.\n"
]
}
],
"source": [
"ansatz_circuit, \\\n",
"init_amplitudes, \\\n",
"ansatz_parameter_names, \\\n",
"hamiltonian_QubitOp, \\\n",
"n_qubits, n_electrons = generate_uccsd(molecule_file, threshold=-1)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"[generate_uccsd](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.generate_uccsd.html) packs functions related to the unitary coupled-cluster, including multiple steps such as deriving a molecular Hamiltonian, building a unitary coupled-cluster ansatz operator, and extracting a coupled-cluster coefficient computed by CCSD. This function reads the molecule by entering its file path. The parameter `th` indicates the to-be-updated gradient threshold of a parameter in the quantum circuit. In the section [Building a Unitary Coupled-Cluster Ansatz Step by Step](#building-a-unitary-coupled-cluster-ansatz-step-by-step), we will demonstrate how to use the related interfaces of MindSpore Quantum to complete the steps. A complete quantum circuit includes an HF initial state and a UCCSD ansatz, as shown in the following code:"
]
},
{
"cell_type": "code",
"execution_count": 7,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
" Circuit Summary \n",
"╭──────────────────────┬──────────────────────────────────────────────╮\n",
"│ Info │ value │\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ Number of qubit │ 12 │\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ Total number of gate │ 12612 │\n",
"│ Barrier │ 2560 │\n",
"│ Noise Channel │ 0 │\n",
"│ Measurement │ 0 │\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ Parameter gate │ 640 │\n",
"│ 44 ansatz parameters │ p0, p8, p1, p9, p2, p10, p3, p11, p4, p12... │\n",
"╰──────────────────────┴──────────────────────────────────────────────╯\n",
"
\n"
],
"text/plain": [
"\u001b[1;38;2;255;0;0m Circuit Summary \u001b[0m\n",
"╭──────────────────────┬──────────────────────────────────────────────╮\n",
"│\u001b[1m \u001b[0m\u001b[1;38;2;59;59;149mInfo\u001b[0m\u001b[1m \u001b[0m\u001b[1m \u001b[0m│\u001b[1m \u001b[0m\u001b[1;38;2;59;59;149mvalue\u001b[0m\u001b[1m \u001b[0m\u001b[1m \u001b[0m│\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ \u001b[1mNumber of qubit\u001b[0m │ 12 │\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ \u001b[1mTotal number of gate\u001b[0m │ 12612 │\n",
"│ Barrier │ 2560 │\n",
"│ Noise Channel │ 0 │\n",
"│ Measurement │ 0 │\n",
"├──────────────────────┼──────────────────────────────────────────────┤\n",
"│ \u001b[1mParameter gate\u001b[0m │ 640 │\n",
"│ 44 ansatz parameters │ \u001b[38;2;72;201;176mp0, p8, p1, p9, p2, p10, p3, p11, p4, p12...\u001b[0m │\n",
"╰──────────────────────┴──────────────────────────────────────────────╯\n"
]
},
"metadata": {},
"output_type": "display_data"
},
{
"name": "stdout",
"output_type": "stream",
"text": [
"Number of parameters: 44\n"
]
}
],
"source": [
"total_circuit = hartreefock_wfn_circuit + ansatz_circuit\n",
"total_circuit.summary()\n",
"print(\"Number of parameters: %d\" % (len(ansatz_parameter_names)))"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"For the LiH molecule, the UCCSD wave function ansatz includes 44 variational parameters. The total number of quantum bit gates of the circuit is 12612, and a total of 12 quantum bits are needed for simulation.\n",
"\n",
"### VQE Procedure\n",
"\n",
"The procedure for solving the molecular ground state by using the VQE is as follows:\n",
"\n",
"1. Prepare the HF initial state: $| 00\\dots11\\dots \\rangle$.\n",
"2. Define the wave function ansatz, such as UCCSD.\n",
"3. Convert the wave function into a variational quantum circuit.\n",
"4. Initialize the variational parameters, for example, set all parameters to 0.\n",
"5. Obtain the energy $E(\\theta)$ of the molecular Hamiltonian under the set of variational parameters and the derivative $\\{ {\\partial E} / {\\partial \\theta_{i}} \\}$ of the energy about the parameters by means of multiple measurements on the quantum computer.\n",
"6. Use optimization algorithms, such as gradient descent and BFGS, to update variational parameters on classical computers.\n",
"7. Transfer the new variational parameters to the quantum circuit for updating.\n",
"8. Repeat steps 5 to 7 until the convergence criteria are met.\n",
"9. End.\n",
"\n",
"In step 5, the derivative $\\{ {\\partial E} / {\\partial \\theta_{i}} \\}$ of the energy about the parameter may be computed by using a parameter-shift rule on a quantum computer, or may be computed by simulating a parameter-shift rule or a finite difference method in a simulator. This is a relatively time-consuming process. Based on the MindSpore framework, MindSpore Quantum provides the automatic derivation function similar to machine learning, which can efficiently compute the derivatives of variational quantum circuits in simulation. The following uses MindSpore Quantum to build a parameterized UCCSD quantum circuit with an automatic derivation function:"
]
},
{
"cell_type": "code",
"execution_count": 8,
"metadata": {},
"outputs": [],
"source": [
"sim = Simulator('mqvector', total_circuit.n_qubits)\n",
"molecule_pqc = sim.get_expectation_with_grad(Hamiltonian(hamiltonian_QubitOp), total_circuit)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"You can obtain the energy $E(\\theta)=\\langle \\Psi_{UCC}(\\theta) | \\hat{H} | \\Psi_{UCC}(\\theta) \\rangle$ corresponding to the variational parameter and the derivative of each variational parameter by transferring a specific value of the parameter to `molecule_pqc`.\n",
"\n",
"For example, we can use the following code to get the expectation of hamiltonian and the corresponding gradient when initial parameters of variational quantum circuit is zero."
]
},
{
"cell_type": "code",
"execution_count": 9,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Energy: [[-7.86335762+0.j]] \n",
"shape: (1, 1) \n",
"\n",
"Gradient: [[[-5.76761758e-11+0.j -8.60518134e-02+0.j 1.19275904e-08+0.j\n",
" -4.85545124e-02+0.j -1.24342308e-14+0.j -3.92769093e-02+0.j\n",
" -3.28108247e-15+0.j -9.59481745e-02+0.j 1.24424913e-14+0.j\n",
" -3.92769093e-02+0.j -3.82408561e-15+0.j -9.59481745e-02+0.j\n",
" -8.60105773e-11+0.j -2.89649669e-02+0.j 2.97699096e-10+0.j\n",
" -4.91813233e-01+0.j -9.35293242e-04+0.j -1.63755737e-16+0.j\n",
" -3.35024598e-17+0.j 1.55629243e-17+0.j 1.59462515e-16+0.j\n",
" 2.67655514e-17+0.j -2.16577662e-17+0.j 5.06813117e-03+0.j\n",
" 1.08542346e-02+0.j -1.28614265e-02+0.j -4.90653908e-17+0.j\n",
" 1.71728860e-17+0.j 1.33973787e-01+0.j -3.03063680e-02+0.j\n",
" -3.37324552e-18+0.j 9.19976050e-19+0.j -1.05649122e-28+0.j\n",
" 5.33586109e-17+0.j -3.68357159e-17+0.j 3.33490478e-17+0.j\n",
" 1.09202635e-28+0.j 1.15391965e-16+0.j -3.03063680e-02+0.j\n",
" -5.27452935e-17+0.j 4.81071347e-17+0.j -4.22508288e-17+0.j\n",
" -2.44510199e-17+0.j -1.68035039e-03+0.j]]] \n",
"shape: (1, 1, 44)\n"
]
}
],
"source": [
"import numpy as np\n",
"\n",
"n_params = len(total_circuit.params_name)\n",
"p0 = np.zeros(n_params)\n",
"f, g = molecule_pqc(p0)\n",
"print(\"Energy: \", f, \"\\nshape: \", f.shape, '\\n')\n",
"print(\"Gradient: \", g, \"\\nshape: \", g.shape)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Throw the above calculation, we get the energy and gradient value and the user can use these data according their practical needs. Now we following the step of (5)~(7) of optimization of VQE, to optimize of variational quantum circuit. Here we use the optimizer in scipy to optimize our quantum circuit. First, we need to define the optimization function that scipy required:"
]
},
{
"cell_type": "code",
"execution_count": 10,
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"(-7.863357621536964,\n",
" array([-5.76761758e-11, -8.60518134e-02, 1.19275904e-08, -4.85545124e-02,\n",
" -1.24342308e-14, -3.92769093e-02, -3.28108247e-15, -9.59481745e-02,\n",
" 1.24424913e-14, -3.92769093e-02, -3.82408561e-15, -9.59481745e-02,\n",
" -8.60105773e-11, -2.89649669e-02, 2.97699096e-10, -4.91813233e-01,\n",
" -9.35293242e-04, -1.63755737e-16, -3.35024598e-17, 1.55629243e-17,\n",
" 1.59462515e-16, 2.67655514e-17, -2.16577662e-17, 5.06813117e-03,\n",
" 1.08542346e-02, -1.28614265e-02, -4.90653908e-17, 1.71728860e-17,\n",
" 1.33973787e-01, -3.03063680e-02, -3.37324552e-18, 9.19976050e-19,\n",
" -1.05649122e-28, 5.33586109e-17, -3.68357159e-17, 3.33490478e-17,\n",
" 1.09202635e-28, 1.15391965e-16, -3.03063680e-02, -5.27452935e-17,\n",
" 4.81071347e-17, -4.22508288e-17, -2.44510199e-17, -1.68035039e-03]))"
]
},
"execution_count": 10,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"def fun(p0, molecule_pqc, energy_list=None):\n",
" f, g = molecule_pqc(p0)\n",
" f = np.real(f)[0, 0]\n",
" g = np.real(g)[0, 0]\n",
" if energy_list is not None:\n",
" energy_list.append(f)\n",
" if len(energy_list) % 5 == 0:\n",
" print(f\"Step: {len(energy_list)},\\tenergy: {f}\")\n",
" return f, g\n",
"\n",
"fun(p0, molecule_pqc)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Here, the `fun` that we define can correctly to return the data that we need: a real energy value, and a array of gradient value with the same size of parameters. Now, we use `bfgs` optimizer in scipy to finish the optimization."
]
},
{
"cell_type": "code",
"execution_count": 11,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Step: 5,\tenergy: -7.880227726111826\n",
"Step: 10,\tenergy: -7.8818171240648365\n",
"Step: 15,\tenergy: -7.882213242905283\n",
"Step: 20,\tenergy: -7.8823453369936445\n",
"Step: 25,\tenergy: -7.8823524949903465\n",
"Step: 30,\tenergy: -7.8823526912721755\n",
"Step: 35,\tenergy: -7.882352703403696\n",
"Step: 40,\tenergy: -7.882352708341612\n"
]
}
],
"source": [
"from scipy.optimize import minimize\n",
"\n",
"energy_list = []\n",
"res = minimize(fun, p0, args=(molecule_pqc, energy_list), method='bfgs', jac=True)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"So, we finished the gradient optimization of variational quantum circuit. Here, `energy_list` is going to store the energy during optimization. Here, we briefly introduce the usage of `minimize`:\n",
"\n",
"- `fun`: The first arg is the function you want to optimize.\n",
"- `p0`: The second arg is the initial value of variables.\n",
"- `args`: The other argument of `fun` except the first argument. According the definition of `fun`, we choose `args=(molecule_pqc, energy_list)`.\n",
"- `method`: The optimization algorithm. Here we use a second order optimization algorithm `bfgs`. For more optimization algorithm, please refer: [scipy tutorial](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html).\n",
"- `jac`: To info that whether `fun` return gradient. Here we use `True`, because MindSpore Quantum can calculate the accuracy gradient value of variational quantum circuit. If use `False`, `minimize` framework will calculate the approximated gradient value base on difference method.\n",
"\n",
"`res` is the optimization result of `scipy`, including optimized parameters, the optimized value and evolution steps."
]
},
{
"cell_type": "code",
"execution_count": 12,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Ground state: \n",
"-7.882352708341612\n",
"\n",
"FCI: \n",
"-7.882362286798721\n",
"\n",
"Optimized amplitudes: \n",
"[ 2.38815317e-04 1.89073105e-03 3.52372251e-02 1.60368559e-02\n",
" -4.30846145e-09 9.09437670e-04 -5.45444468e-10 1.41627870e-02\n",
" 3.62569484e-09 9.08701303e-04 -1.05938178e-10 1.41711693e-02\n",
" -5.47703644e-04 4.26818370e-04 2.87168405e-03 5.38109720e-02\n",
" 2.34920831e-04 -9.88498886e-08 -2.51432020e-07 2.54936335e-07\n",
" -2.65409411e-08 4.64583011e-08 -4.42267421e-08 1.32809550e-05\n",
" -1.04167881e-04 7.98983208e-04 5.39567507e-10 -2.28346700e-10\n",
" -5.50005419e-02 3.09112487e-03 -4.88778663e-10 -4.89650670e-09\n",
" 4.92879467e-09 1.13389943e-07 4.45780819e-09 -1.85315776e-08\n",
" -1.42394064e-12 4.57319777e-10 3.09100040e-03 5.36042785e-08\n",
" -3.87192856e-09 9.99056234e-10 1.51453640e-10 3.72805204e-04]\n"
]
}
],
"source": [
"print(f\"Ground state: \\n{res.fun}\\n\")\n",
"print(f\"FCI: \\n-7.882362286798721\\n\")\n",
"print(f\"Optimized amplitudes: \\n{res.x}\")"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"We can see here the result of ucc method is very close to FCI method with very good accuracy.\n",
"\n",
"## Building a Unitary Coupled-Cluster Ansatz Step by Step\n",
"\n",
"\n",
"\n",
"In the preceding part, the [generate_uccsd](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.generate_uccsd.html) is used to build all the content required for designing the unitary coupled-cluster. In this section, the steps are split, we get the coupled-cluster operator, the corresponding quantum circuit and the initial guess of the variational parameters from the classical CCSD results.\n",
"First, import some extra dependencies, including the related functions in MindSpore Quantum."
]
},
{
"cell_type": "code",
"execution_count": 13,
"metadata": {},
"outputs": [],
"source": [
"from mindquantum.algorithm.nisq import Transform\n",
"from mindquantum.algorithm.nisq import get_qubit_hamiltonian\n",
"from mindquantum.algorithm.nisq import uccsd_singlet_generator, uccsd_singlet_get_packed_amplitudes\n",
"from mindquantum.core.operators import TimeEvolution"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The molecule Hamiltonian uses [get_qubit_hamiltonian](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.get_qubit_hamiltonian.html) to read the previous computing result. The result is as follows:"
]
},
{
"cell_type": "code",
"execution_count": 14,
"metadata": {},
"outputs": [],
"source": [
"hamiltonian_QubitOp = get_qubit_hamiltonian(molecule_of)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The unitary coupled-cluster operator $\\hat{T} - \\hat{T}^{\\dagger}$ can be built using [uccsd_singlet_generator](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.uccsd_singlet_generator.html). Provide the total number of quantum bits (total number of spin orbits) and the total number of electrons, and set `anti_hermitian=True`."
]
},
{
"cell_type": "code",
"execution_count": 15,
"metadata": {},
"outputs": [],
"source": [
"ucc_fermion_ops = uccsd_singlet_generator(\n",
" molecule_of.n_qubits, molecule_of.n_electrons, anti_hermitian=True)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The `ucc_fermion_ops` built in the previous step is parameterized. Use the Jordan-Wigner transformation to map the Fermi excitation operator to the Pauli operator:"
]
},
{
"cell_type": "code",
"execution_count": 16,
"metadata": {},
"outputs": [],
"source": [
"ucc_qubit_ops = Transform(ucc_fermion_ops).jordan_wigner()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Next, we need to obtain the quantum circuit corresponding to the unitary operator $\\exp{(\\hat{T} - \\hat{T}^{\\dagger})}$. [TimeEvolution](https://www.mindspore.cn/mindquantum/docs/en/master/core/operators/mindquantum.core.operators.TimeEvolution.html) can generate the circuit corresponding to $\\exp{(-i\\hat{H}t)}$, where $\\hat{H}$ is a Hermitian operator, and $t$ is a real number. Note that when [TimeEvolution](https://www.mindspore.cn/mindquantum/docs/en/master/core/operators/mindquantum.core.operators.TimeEvolution.html) is used, `ucc_qubit_ops` already contains the complex number factor $i$. Therefore, you need to divide `ucc_qubit_ops` by $i$ or extract the imaginary part of `ucc_qubit_ops`."
]
},
{
"cell_type": "code",
"execution_count": 17,
"metadata": {},
"outputs": [],
"source": [
"ansatz_circuit = TimeEvolution(ucc_qubit_ops.imag, 1.0).circuit\n",
"ansatz_parameter_names = ansatz_circuit.params_name"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"`ansatz_parameter_names` is used to record the parameter names in the circuit. So far, we have obtained the contents required by the VQE quantum circuit, including the Hamiltonian `hamiltonian_QubitOp` and the parameterized wave function ansatz `ansatz_circuit`. By referring to the preceding steps, we can obtain a complete state preparation circuit. `hartreefock_wfn_circuit` mentioned above is used as the Hartree-Fock reference state:"
]
},
{
"cell_type": "code",
"execution_count": 18,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
" Circuit Summary \n",
"╭──────────────────────┬──────────────────────────────────────────────────────────╮\n",
"│ Info │ value │\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ Number of qubit │ 12 │\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ Total number of gate │ 12612 │\n",
"│ Barrier │ 2560 │\n",
"│ Noise Channel │ 0 │\n",
"│ Measurement │ 0 │\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ Parameter gate │ 640 │\n",
"│ 44 ansatz parameters │ s_0, d1_0, s_1, d1_1, s_2, d1_2, s_3, d1_3, s_4, d1_4... │\n",
"╰──────────────────────┴──────────────────────────────────────────────────────────╯\n",
"
\n"
],
"text/plain": [
"\u001b[1;38;2;255;0;0m Circuit Summary \u001b[0m\n",
"╭──────────────────────┬──────────────────────────────────────────────────────────╮\n",
"│\u001b[1m \u001b[0m\u001b[1;38;2;59;59;149mInfo\u001b[0m\u001b[1m \u001b[0m\u001b[1m \u001b[0m│\u001b[1m \u001b[0m\u001b[1;38;2;59;59;149mvalue\u001b[0m\u001b[1m \u001b[0m\u001b[1m \u001b[0m│\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ \u001b[1mNumber of qubit\u001b[0m │ 12 │\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ \u001b[1mTotal number of gate\u001b[0m │ 12612 │\n",
"│ Barrier │ 2560 │\n",
"│ Noise Channel │ 0 │\n",
"│ Measurement │ 0 │\n",
"├──────────────────────┼──────────────────────────────────────────────────────────┤\n",
"│ \u001b[1mParameter gate\u001b[0m │ 640 │\n",
"│ 44 ansatz parameters │ \u001b[38;2;72;201;176ms_0, d1_0, s_1, d1_1, s_2, d1_2, s_3, d1_3, s_4, d1_4...\u001b[0m │\n",
"╰──────────────────────┴──────────────────────────────────────────────────────────╯\n"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"total_circuit = hartreefock_wfn_circuit + ansatz_circuit\n",
"total_circuit.summary()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Next, you need to provide a reasonable initial value for the variational parameter. The `PQCNet` built in the preceding text uses 0 as the initial guess by default, which is feasible in most cases. However, using CCSD's computational data as a starting point for UCC may be better. Use the [uccsd_singlet_get_packed_amplitudes](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.uccsd_singlet_get_packed_amplitudes.html) function to extract CCSD parameters from `molecule_of`."
]
},
{
"cell_type": "code",
"execution_count": 19,
"metadata": {},
"outputs": [],
"source": [
"init_amplitudes_ccsd = uccsd_singlet_get_packed_amplitudes(\n",
" molecule_of.ccsd_single_amps, molecule_of.ccsd_double_amps, molecule_of.n_qubits, molecule_of.n_electrons)\n",
"init_amplitudes_ccsd = [init_amplitudes_ccsd[param_i] for param_i in ansatz_parameter_names]"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Just like the previous method, we can get the `grad_ops` with MindSpore Quantum, and optimize it with scipy."
]
},
{
"cell_type": "code",
"execution_count": 20,
"metadata": {},
"outputs": [],
"source": [
"grad_ops = Simulator('mqvector', total_circuit.n_qubits).get_expectation_with_grad(\n",
" Hamiltonian(hamiltonian_QubitOp.real),\n",
" total_circuit)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"`init_amplitudes_ccsd` (coupled-cluster coefficient computed by CCSD) is used as an initial variational parameter:"
]
},
{
"cell_type": "code",
"execution_count": 21,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Step: 5,\tenergy: -7.878223283110822\n",
"Step: 10,\tenergy: -7.880288481802435\n",
"Step: 15,\tenergy: -7.8820356683191415\n",
"Step: 20,\tenergy: -7.882302370884269\n",
"Step: 25,\tenergy: -7.882349803535784\n",
"Step: 30,\tenergy: -7.882352702053716\n",
"Step: 35,\tenergy: -7.88235270792034\n",
"Step: 40,\tenergy: -7.882352708346457\n"
]
}
],
"source": [
"energy_list = []\n",
"res = minimize(fun, init_amplitudes_ccsd, args=(grad_ops, energy_list), method='bfgs', jac=True)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The final optimized result is shown as below."
]
},
{
"cell_type": "code",
"execution_count": 22,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Ground state: \n",
"-7.882352708346457\n",
"\n",
"FCI: \n",
"-7.882362286798721\n",
"\n",
"Optimized amplitudes: \n",
"[-2.38636816e-04 1.89071890e-03 -3.52372349e-02 1.60368129e-02\n",
" 1.34665623e-08 9.09430837e-04 -4.17115460e-10 1.41641317e-02\n",
" -1.28144111e-08 9.08694024e-04 4.20204892e-10 1.41698036e-02\n",
" 5.47710797e-04 4.26820871e-04 -2.87169433e-03 5.38109417e-02\n",
" 2.34627619e-04 4.32417990e-07 -1.26325933e-08 -1.39715964e-08\n",
" -4.27213491e-07 1.24814811e-08 1.37984739e-08 1.32946010e-05\n",
" 7.99032002e-04 -1.04121599e-04 1.25575236e-09 -1.24558942e-09\n",
" -5.50005457e-02 3.09142613e-03 1.53511697e-08 -3.40030804e-08\n",
" 3.29912584e-08 3.67660373e-07 -2.59584687e-08 2.22720149e-10\n",
" 4.47254154e-11 1.06718989e-09 3.09086800e-03 -3.54929776e-07\n",
" 1.95048358e-08 5.27817957e-09 -1.03614664e-09 3.72830373e-04]\n"
]
}
],
"source": [
"print(f\"Ground state: \\n{res.fun}\\n\")\n",
"print(f\"FCI: \\n-7.882362286798721\\n\")\n",
"print(f\"Optimized amplitudes: \\n{res.x}\")"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Summary\n",
"\n",
"In this case, the ground state energy of the LiH molecule is obtained by using scipy in two methods. In the first method, we use the [generate_uccsd](https://www.mindspore.cn/mindquantum/docs/en/master/algorithm/nisq/mindquantum.algorithm.nisq.generate_uccsd.html) function packaged by MindSpore Quantum to generate a quantum neural network that can solve this problem. In the second method, we build a similar gradient operator step by step. The final results are the same."
]
},
{
"cell_type": "code",
"execution_count": 23,
"metadata": {},
"outputs": [
{
"data": {
"text/html": [
"\n",
"\n",
" \n",
" Software | \n",
" Version | \n",
"
\n",
"mindquantum | 0.9.11 |
\n",
"scipy | 1.10.1 |
\n",
"numpy | 1.23.5 |
\n",
"\n",
" System | \n",
" Info | \n",
"
\n",
"Python | 3.9.16 |
OS | Linux x86_64 |
Memory | 8.3 GB |
CPU Max Thread | 8 |
Date | Mon Jan 1 00:25:57 2024 |
\n",
"
\n"
],
"text/plain": [
""
]
},
"execution_count": 23,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"from mindquantum.utils.show_info import InfoTable\n",
"\n",
"InfoTable('mindquantum', 'scipy', 'numpy')"
]
}
],
"metadata": {
"kernelspec": {
"display_name": "base",
"language": "python",
"name": "python3"
},
"language_info": {
"codemirror_mode": {
"name": "ipython",
"version": 3
},
"file_extension": ".py",
"mimetype": "text/x-python",
"name": "python",
"nbconvert_exporter": "python",
"pygments_lexer": "ipython3",
"version": "3.9.16"
}
},
"nbformat": 4,
"nbformat_minor": 2
}