{ "cells": [ { "cell_type": "markdown", "id": "0b054e8d-d0bb-4613-bddd-81b37f3d2a8e", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "# Basics\n", "\n", "## Basic Representation\n", "\n", "States and operators in {py:mod}`quimb` are simply dense numpy arrays\n", "or sparse scipy matrices. All functions should directly work with these\n", "but the class {class}`~quimb.core.qarray` is also provided as a very\n", "thin subclass of {class}`numpy.ndarray` with a few helpful methods and\n", "attributes. The {py:func}`~quimb.core.quimbify` function (aliased to\n", "{py:func}`~quimb.core.qu`) can convert between the various representations." ] }, { "cell_type": "code", "execution_count": 1, "id": "a30d59e4-ebf6-417d-9adf-16de0b1bbbe1", "metadata": { "raw_mimetype": "text/restructuredtext" }, "outputs": [], "source": [ "from quimb import *\n", "\n", "data = [1, 2j, -3]" ] }, { "cell_type": "markdown", "id": "1f139c32-2483-406f-8a4c-399cf1945a34", "metadata": {}, "source": [ "Kets are column vectors, i.e. with shape ``(d, 1)``:" ] }, { "cell_type": "code", "execution_count": 2, "id": "69b02c9c-37c3-427d-99de-b480b7ea4d6d", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 1.+0.j]\n", " [ 0.+2.j]\n", " [-3.+0.j]]" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "qu(data, qtype=\"ket\")" ] }, { "cell_type": "markdown", "id": "4a731dcf-13e7-4888-a7d1-e7defb06edd3", "metadata": {}, "source": [ "The ``normalized=True`` option can be used to ensure a normalized output.\n", "\n", "Bras are row vectors, i.e. with shape ``(1, d)``:" ] }, { "cell_type": "code", "execution_count": 3, "id": "0770dcd6-e29b-4b92-9072-2b4ade641268", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 1.-0.j 0.-2.j -3.-0.j]]" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "qu(data, qtype=\"bra\") # also conjugates the data" ] }, { "cell_type": "markdown", "id": "f2ab8914-8431-4dad-bb4b-7a6c400c60df", "metadata": {}, "source": [ "And operators are square matrices, i.e. have shape ``(d, d)``:" ] }, { "cell_type": "code", "execution_count": 4, "id": "6a3a349c-8648-4934-8000-66d3ddf16805", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 1.+0.j 0.-2.j -3.-0.j]\n", " [ 0.+2.j 4.+0.j 0.-6.j]\n", " [-3.+0.j 0.+6.j 9.+0.j]]" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "qu(data, qtype=\"dop\")" ] }, { "cell_type": "markdown", "id": "fb452078-ab1c-4ba5-81e0-799af6ef0f76", "metadata": {}, "source": [ "Which can also be sparse:" ] }, { "cell_type": "code", "execution_count": 5, "id": "24cab7c2-7dd3-46a7-9acf-f1d3e0cbe5d3", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "qu(data, qtype=\"dop\", sparse=True)" ] }, { "cell_type": "markdown", "id": "a10ffe96-9bc3-4e4e-ace8-0098d0e29564", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The sparse format can be specified with the `stype` keyword. The partial\n", "function versions of each of the above are also available:\n", "\n", "- {func}`~quimb.core.ket`\n", "- {func}`~quimb.core.bra`\n", "- {func}`~quimb.core.dop`\n", "- {func}`~quimb.core.sparse`\n", "\n", ":::{note}\n", "If a simple 1d-list is supplied and no `qtype` is given, `'ket'` is\n", "assumed.\n", ":::\n", "\n", "## Basic Operations\n", "\n", "The 'dagger', or hermitian conjugate, operation is performed with the `.H`\n", "attribute:" ] }, { "cell_type": "code", "execution_count": 6, "id": "7bebe19b-bb7e-4d5d-a55c-63340574aba4", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 0.+0.j ]\n", " [ 0.+0.707107j]\n", " [-0.-0.707107j]\n", " [ 0.+0.j ]]" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "psi = 1.0j * bell_state(\"psi-\")\n", "psi" ] }, { "cell_type": "code", "execution_count": 7, "id": "0cb44044-6832-4365-a4b0-5f6ed7b940c4", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 0.-0.j 0.-0.707107j -0.+0.707107j 0.-0.j ]]" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "psi.H" ] }, { "cell_type": "markdown", "id": "9389db86-974f-47a7-94e0-729a34436be2", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "This is just the combination of `.conj()` and `.T`, but only available for\n", "{mod}`scipy.sparse` matrices and {class}`~quimb.core.qarray` s (not\n", "{class}`numpy.ndarray` s).\n", "\n", "The product of two quantum objects is the dot or matrix product, which, since\n", "python 3.5, has been overloaded with the `@` symbol. Using it is recommended:" ] }, { "cell_type": "code", "execution_count": 8, "id": "ccb3c8b1-3673-496b-9efe-d15e13e15c45", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[1.+0.j]\n", " [0.+0.j]]" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "psi = up()\n", "psi" ] }, { "cell_type": "code", "execution_count": 9, "id": "21c148d7-a55c-44cb-bbf3-8de4c9b33f41", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[1.+0.j]]" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "psi.H @ psi # inner product" ] }, { "cell_type": "code", "execution_count": 10, "id": "6916a93e-bbfa-4b37-92f3-a7816d58af28", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[0.+0.j]\n", " [1.+0.j]]" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "X = pauli(\"X\")\n", "X @ psi # act as gate" ] }, { "cell_type": "code", "execution_count": 11, "id": "ad4186dc-309d-44b9-abed-5d6af268746a", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[0.+0.j]]" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "(psi.H @ X @ psi) # operator expectation" ] }, { "cell_type": "markdown", "id": "ba59088e-d798-455f-95b8-179512b0874f", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "Scalar expectation values might best be computed using the\n", "{func}`~quimb.core.expectation` function (aliased to\n", "{func}`~quimb.core.expec`) which dispatches to accelerated\n", "methods:" ] }, { "cell_type": "code", "execution_count": 12, "id": "4012d11b-f4b9-4dd5-af69-96f51299839d", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "np.float64(1.0)" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "expec(psi, psi)" ] }, { "cell_type": "code", "execution_count": 13, "id": "a67c2946-6e4b-40e6-91ff-f47c17e5f0d7", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "np.complex128(0j)" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "expec(psi, X)" ] }, { "cell_type": "markdown", "id": "cd62be41-0bf5-422e-b5b5-8cb9a37174d6", "metadata": {}, "source": [ "Here's an example for a much larger (20 qubit), sparse operator expecation,\n", "which will be automatically parallelized:" ] }, { "cell_type": "code", "execution_count": 14, "id": "425bae86-798a-4d9d-a7e0-2690ea7a219b", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "psi = rand_ket(2**20)\n", "A = rand_herm(2**20, sparse=True) + speye(2**20)\n", "A" ] }, { "cell_type": "code", "execution_count": 15, "id": "4198dbc2-d143-4514-9563-8785ef1395f4", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "np.float64(1.0004858170636874)" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "expec(A, psi) # should be ~ 1" ] }, { "cell_type": "code", "execution_count": 16, "id": "bda59aa7-95c5-4dbf-9b63-95a3635f1ec9", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "55.1 ms ± 3.82 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)\n" ] } ], "source": [ "%%timeit\n", "expec(A, psi)" ] }, { "cell_type": "markdown", "id": "e68195b9-331c-4e07-b3b2-a1276beb808e", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "## Combining Objects - Tensor Products\n", "\n", "There are a number of ways to combine states and operators, i.e. tensoring them\n", "together.\n", "\n", "Functional form using {func}`~quimb.core.kron`:\n", "\n", "```python\n", ">>> kron(psi1, psi2, psi3, ...)\n", "...\n", "```\n", "\n", "This can also be done using the `&` overload on `qarray` and scipy matrices:\n", "\n", "```python\n", ">>> psi1 & psi2 & psi3\n", "...\n", "```\n", "\n", ":::{warning}\n", "When {mod}`quimb` is imported, it monkey patches the otherwise unused\n", "method of `&`/`__and__` of scipy sparse matrices to {func}`~quimb.core.kron`.\n", ":::\n", "\n", "Often one wants to sandwich an operator with many identities,\n", "{func}`~quimb.core.ikron` can be used for this:" ] }, { "cell_type": "code", "execution_count": 17, "id": "bf597230-30d7-4d54-b876-e93967891fc7", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "(1024, 1024)" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dims = [2] * 10 # overall space of 10 qubits\n", "X = pauli(\"X\")\n", "IIIXXIIIII = ikron(X, dims, inds=[3, 4]) # act on 4th and 5th spin only\n", "IIIXXIIIII.shape" ] }, { "cell_type": "markdown", "id": "5b7a4905-ee44-4a82-89f1-6e1a612c021b", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "For more advanced tensor constructions, such as reversing and interleaving\n", "identities within operators {func}`~quimb.core.pkron` can be used:" ] }, { "cell_type": "code", "execution_count": 18, "id": "ed53d5cc-a2e7-4acb-b867-cb1e61951d3f", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 0 1 0 0 0 0 0 0]\n", " [ 1 0 0 0 0 0 0 0]\n", " [ 0 0 0 1 0 0 0 0]\n", " [ 0 0 1 0 0 0 0 0]\n", " [ 0 0 0 0 0 -1 0 0]\n", " [ 0 0 0 0 -1 0 0 0]\n", " [ 0 0 0 0 0 0 0 -1]\n", " [ 0 0 0 0 0 0 -1 0]]" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dims = [2] * 3\n", "XZ = pauli(\"X\") & pauli(\"Z\")\n", "ZIX = pkron(XZ, dims, inds=[2, 0])\n", "ZIX.real.astype(int)" ] }, { "cell_type": "markdown", "id": "753aebb7-8d79-491b-8add-9b2497e6d81a", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "`ZIX` would then act with Z on first spin, and X on 3rd.\n", "\n", "## Removing Objects - Partial Trace\n", "\n", "To remove, or ignore, certain parts of a quantum state the partial trace\n", "function {func}`~quimb.core.partial_trace` (aliased to {func}`~quimb.core.ptr`)\n", "is used. Here, the internal dimensions of a state must be supplied as well as\n", "the indicies of which of these subsystems to *keep*.\n", "\n", "For example, if we have a random system of 10 qubits (hilbert space of dimension\n", "`2**10`), and we want just the reduced density matrix describing the first and\n", "last spins:" ] }, { "cell_type": "code", "execution_count": 19, "id": "5c6f0168-e595-4472-b78e-e7c551089e31", "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[[ 0.279+0.j 0.003-0.005j 0.012+0.009j 0.018-0.009j]\n", " [ 0.003+0.005j 0.239+0.j -0.005+0.007j -0.005-0.002j]\n", " [ 0.012-0.009j -0.005-0.007j 0.24 +0.j -0.015-0.014j]\n", " [ 0.018+0.009j -0.005+0.002j -0.015+0.014j 0.242+0.j ]]" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dims = [2] * 10\n", "D = prod(dims)\n", "psi = rand_ket(D)\n", "rho_ab = ptr(psi, dims, [0, 9])\n", "rho_ab.round(3) # probably pretty close to identity" ] }, { "cell_type": "markdown", "id": "5ea15f65-79d4-40a1-993a-82b011dc77e1", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "{func}`~quimb.core.partial_trace` accepts dense or sparse, operators or vectors." ] } ], "metadata": { "celltoolbar": "Raw Cell Format", "kernelspec": { "display_name": "Python 3", "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" } }, "nbformat": 4, "nbformat_minor": 2 }