{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# The Su-Schrieffer–Heeger (SSH) model" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Saumya Biswas (saumyab@uoregon.edu)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The celebrated SSH model is analyzed with QuTiP's lattice module below." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "![title](images/SSH_E.png)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The above figure shows a SSH model with 6 sites with periodic boundary condition. The same lattice with hardwall/aperiodic boundary condition would be the folloowing.\n", "![title](images/SSH_Ea.png)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In the secod quantized formalism, the periodic lattice Hamiltonian can be written as\n", "\\begin{eqnarray}\n", "H_{per} = -t_{intra} (c_{-2}^{\\dagger} c_{-1} + c_{-1}^{\\dagger} c_{-2} ) -t_{intra} (c_{0}^{\\dagger} c_{1} + c_{1}^{\\dagger} c_{0} ) -t_{intra} (c_{2}^{\\dagger} c_{3} + c_{3}^{\\dagger} c_{2} ) \\nonumber \\\\ \n", " -t_{inter} (c_{-1}^{\\dagger} c_{0} + c_{0}^{\\dagger} c_{-1} ) -t_{inter} (c_{1}^{\\dagger} c_{2} + c_{2}^{\\dagger} c_{1} ) -t_{inter} (c_{3}^{\\dagger} c_{-2} + c_{-2}^{\\dagger} c_{3} ) \\nonumber \n", "\\end{eqnarray}\n", "The aperiodic lattice Hamiltonian can be obtained by discarding the very last term.\n", "\\begin{eqnarray}\n", "H_{aper} = -t_{intra} (c_{-2}^{\\dagger} c_{-1} + c_{-1}^{\\dagger} c_{-2} ) -t_{intra} (c_{0}^{\\dagger} c_{1} + c_{1}^{\\dagger} c_{0} ) -t_{intra} (c_{2}^{\\dagger} c_{3} + c_{3}^{\\dagger} c_{2} ) \\nonumber \\\\ \n", " -t_{inter} (c_{-1}^{\\dagger} c_{0} + c_{0}^{\\dagger} c_{-1} ) -t_{inter} (c_{1}^{\\dagger} c_{2} + c_{2}^{\\dagger} c_{1} ) \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\ \\nonumber \n", "\\end{eqnarray}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The representation in terms of unit cell blocks become obvious once we resolve the terms into unit cell operators.\n", "\\begin{eqnarray}\n", "H_{per}= \\begin{bmatrix}\n", " c_{-2}^{\\dagger} & c_{-1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_1 \\\\\n", " -t_1 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{-2} \\\\\n", " c_{-1} \n", "\\end{bmatrix} + \n", "\\begin{bmatrix}\n", " c_{0}^{\\dagger} & c_{1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_1 \\\\\n", " -t_1 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{0} \\\\\n", " c_{1} \n", "\\end{bmatrix} \n", "\\nonumber \\\\\n", "+ \\begin{bmatrix}\n", " c_{2}^{\\dagger} & c_{3}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_1 \\\\\n", " -t_1 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{2} \\\\\n", " c_{3} \n", "\\end{bmatrix} \\nonumber \\\\\n", "+ \\begin{bmatrix}\n", " c_{-2}^{\\dagger} & c_{-1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & 0 \\\\\n", " -t_2 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{0} \\\\\n", " c_{1} \n", "\\end{bmatrix} + \n", "\\begin{bmatrix}\n", " c_{0}^{\\dagger} & c_{1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_2 \\\\\n", " 0 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{-2} \\\\\n", " c_{-1} \n", "\\end{bmatrix} \n", "\\nonumber \\\\\n", "+ \\begin{bmatrix}\n", " c_{0}^{\\dagger} & c_{1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & 0 \\\\\n", " -t_2 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{2} \\\\\n", " c_{3} \n", "\\end{bmatrix} + \n", "\\begin{bmatrix}\n", " c_{2}^{\\dagger} & c_{3}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_2 \\\\\n", " 0 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{0} \\\\\n", " c_{1} \n", "\\end{bmatrix} \n", "\\nonumber \\\\\n", "+ \\begin{bmatrix}\n", " c_{2}^{\\dagger} & c_{3}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & 0 \\\\\n", " -t_2 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{-2} \\\\\n", " c_{-1} \n", "\\end{bmatrix} + \n", "\\begin{bmatrix}\n", " c_{-2}^{\\dagger} & c_{-1}^{\\dagger} \n", "\\end{bmatrix} \n", "\\begin{bmatrix}\n", " o & -t_2 \\\\\n", " 0 & 0 \n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " c_{2} \\\\\n", " c_{3} \n", "\\end{bmatrix} \\nonumber\n", "\\end{eqnarray}\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence, $H_{TB}$ can be succinctly put in the form:\n", "\\begin{eqnarray}\n", "H_{per/aper} = \\sum_i \\psi_i^{\\dagger} D \\psi_i + \\sum_{i} \\left( \\psi_i^{\\dagger} T \\psi_{i+1} + \\psi_{i+1}^{\\dagger} T^{\\dagger} \\psi_i \\right) \\label{eq:TB_block} \n", "\\end{eqnarray}\n", "where $D = \\begin{bmatrix}\n", " 0 & -t_{intra} \\\\\n", " -t_{intra} & 0 \n", "\\end{bmatrix}$ and $T = \\begin{bmatrix}\n", " 0 & 0 \\\\\n", " -t_{inter} & 0 \n", "\\end{bmatrix}$. The $\\psi_i = \\begin{bmatrix}\n", " c^0_{i} \\\\\n", " c^1_{i} \n", "\\end{bmatrix}$ is associated with $\\vec{R}=\\vec{R}_i$ and $\\psi_{i+1}$ is associated with $\\vec{R}=\\vec{R}_i + \\hat{x}$, with $\\hat{x}$ being the unit vector along the x direction.\n", "\n", "The equation above can be put into the alternate form (also changing the summation variables to m,n to distinguish them from the imaginary i):\n", "\n", "\\begin{eqnarray}\n", "H_{per/aper} = \\sum_{m,n} \\psi_m^{\\dagger} (D \\delta_{m,n} + T \\delta_{m,n-1} + T^{\\dagger} \\delta_{m,n+1} ) \\psi_n \\label{eq:TB_block_1} \n", "\\end{eqnarray}\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "So far, we have see that a SSH lattice can be put into the unit cell representations and the lattice can be completely defied with two matrices D and T. In declaring an instance of Qutip.lattice.Lattice1d we only need to input the two matrices. " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In declaring an instance of Lattice1d class, the two arguments cell_Hamiltonian\n", "and inter_hop are set to D and T respectively. The matrix structure of D and T can be\n", "obtained from the aide function cell_structures()." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from qutip import *\n", "import matplotlib.pyplot as plt\n", "import numpy as np" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "val_s = ['site0','site1']\n", "(H_cell_form,T_inter_cell_form,H_cell,T_inter_cell) = cell_structures( val_s)" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[['', ''],\n", " ['', '']]" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "H_cell_form" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[['', ''],\n", " ['', '']]" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "T_inter_cell_form" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Guided by cell_H_form and inter_cell_T_form, we can set values to cell_H and inter_cell_T which were initialized to all zero elements by cell_structures()." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "t_intra = -0.5\n", "t_inter = -0.5\n", "H_cell[0,1] = t_intra\n", "H_cell[1,0] = t_intra\n", "T_inter_cell[1,0] = t_inter\n", "H_cell = Qobj(H_cell)\n", "T_inter_cell = Qobj(T_inter_cell)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Using cell_structures() is completely optional. The user could equally well have defined cell_H and inter_cell_T directly." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "H_cell = Qobj( np.array( [[ 0, t_intra ],[t_intra,0]] ) )\n", "T_inter_cell = Qobj( np.array( [[ 0, 0 ],[t_inter,0]] ) ) " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For our SSH lattice with 3 unit cells, 2 sites in each unit cell, and [1] degree of freedom per each site, we can initiate an instance of the Lattice1d class at this stage. " ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "boundary_condition = \"periodic\"\n", "cells = 3\n", "cell_sites = 2\n", "site_dof = [1]\n", "SSH_lattice = Lattice1d(num_cell=cells, boundary = boundary_condition, cell_num_site = cell_sites, cell_site_dof = site_dof, Hamiltonian_of_cell = H_cell, inter_hop = T_inter_cell )" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The model can be visualized with the display functions." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAQ4AAACWCAYAAAAv82Y/AAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAADmFJREFUeJzt3XuUldV5x/HvE1FjS+Ms09aCqRq8JVYaOrCQLFKDFVyUBBItlZW0VCuGUjVpo00qq+peOzaixmqqERTSlZjWC/FSAo26lChWV5paLl5qqBXUrrSMaFXipV4wPv1jb3EMDMw7nHP2ec/5fdaaNWfemTnznHnP+Z2933e/e5u7IyJSxXtKFyAi9aPgEJHKFBwiUpmCQ0QqU3CISGUKDhGpTMEhIpUpOESkMgWHiFSm4BCRyoaVLqAKM1sH/AqwoXQtIm3uUOBZd/+tZtx5rYKDFBrDSxchUgNNfZ3ULTg2ALj7pMJ1iLQ1M1vVzPvXMQ4RqaxuLQ5pIYu2B3AEqzmRH3AWB9DHp3mSn3AIyzicibzKMWzl39mb9/ImzwPj2QRsAvpYymHM4kZgPbDGgz9T8vFI4yg4ZBuL9j7geOBjwFhgDPA041jLWjYzmhUM5wE+TB83s5hhfJolfIAXOYiz+Scu54/ZzFNMZwu38Fme5BDgIOATQK9FewlYA6wG7iSFyVtlHq3sDqvTRD5v99t0jKNxLNqBwHRgBvBR4H7gHtILfK0H3wJgZjcBC9x9bf56DXAcMB+4y91XmtlkoNfdL8k/c5e7T8l/x4BRpEA6GvhdoAdYASwH7vbgr7bmUXe+Zr9W1OLoQhZtL+AE4HTgN4DvA9cAMz34SxXvrgd4vt/X79/RD3lwBzbmj+8CZ1u0w0ih9SXgHyzaDcAiD/5IxRqkxRQcXcSi/TowFziNdNzhSuB7HnzrIO9ilpmNyrd7drceD/44cBlwmUU7INd1h0XbCCwEbvXgb+zu35HGU3B0AYu2P3Au8FngeuA4D/7jIdzV0n5dlfl52xZgv3y7B3huKDV68P8BokW7kNRtOhNYYNHOA27w4D8byv1Kcyg4Olg+2PkXwBnAd4APefBnG/xnlgLj8u1RwMrdubPc+rkFuMWiHQNcBHzZos0HbstdHilM4zg6kEUzi3YK8DhwIDDWg39xqKFhZr1AL6mr0pMPgo4C5vZrgUwGtvT7eiYwzszmDvVxePB/BiYC5wGXAPdYtMOHen/SODqr0mHysYLFwEjgVA++rnBJDZHHlJwBnA98FbhC3ZeBNfu1ohZHh+jXylgHPACM75TQAPDgP/PgVwATSGeE7s1nZaQABUcHsGi/ANwA/DkwxYPHCmdKasWDbwAmkU7p/tCinVS2ou6k4Ki5PIDrfuANYIIHf6hwSU3nwd/KrY8pwCUW7QKLpudyC+mfXWMWbSLwI+A64GQP/lrhklrKgz8IjCe1QG61aL9UtqLuoeCoqdxE/0dgjgf/m249TZkvnDsOeAa4P49ZkSZTcNSQRfsj4OvAZA9+e+l6SsujS/+ENP7j3nxmSZpIwVEzOTQWkEZ/Ply6nnbhwd2DfwX4FrDKoo0oXVMnU3DUSO6eXExqaawvXU878uAXk8LjBxbtV0vX06kUHDWRD4R+A5iq0Ng5D34hsAxYZtH2Ll1PJ1Jw1EA+5XoTcEo3nG5tkHOBPmBRngtEGkjB0eby4K5lwOUe/LbS9dRFnlnsZNLEQV8oXE7HUXC0sfxO+XfAo8ClhcupHQ/+MvAp4ByLNrl0PZ1EwdHeZgNHAp/r1nEau8uDPwX8AXCtRdtvFz8ug6TgaFMWbSSpldF1I0IbzYPfDdxKGvsiDaDgaEO5i7IYWJiHVcvuOweYaNGmly6kEyg42tNs4APAhaUL6RQe/BVgDnC1uiy7T8HRZvJ0f18jTcKjiXobyIOvIi3HcH7hUmpPwdF+zgLu9JCm4JOGC8Bsi3Zw4TpqTcHRRvIQ6c+jd8Sm8eCbgauAWLqWOlNwtJe/Aq7z4E+WLqTDXQpMtWhHlS6krhQcbSJfCj4b+OvStXQ6D/4iadkFtTqGSMHRPj4H3KgV3VtmCXBsXt1OKlJwtAGLtidpacZFpWvpFnk4+nWk/7tUpOBoD58CNmix5ZZbBJyWF+GWChQc7eF00iLL0kJ5/dz/AE4sXUvdKDgKywdFP0K6lkJa75uki+CkAgVHeZ8Ebtco0WJuAz6e5z2RQVJwlDeDNAxaCvDgLwCrAc3XUYGCoyCLNhz4beCO0rV0ueWkAJdBUnCUNRn4Vw/+09KFdLkVwHQtIzl4+keV9THgntJFdDsPvhF4DTikdC11oeAoayywpnQRAqT9MLZ0EXWh4CgkN4t7UXC0CwVHBQqOcg4FXvDg/1u6EAHSmRUFxyApOMrpBTRZT/tYA/Rq8abBUXCUczCwoXQRkuSW33uAfUvXUgcKjnJGkpYoHJCZ9ZrZRjO7xsx6zGyymb1gZl/O35+Zt83t9zvbbZNB20TaLwPa1T7JP3PxDn5vu211puAoZwTpiTogd19L6s5c4+5b3H0l8ASw2Mx688+shG1P6O22NbH+TtRH2i8D2tk+AciBPbP/7+xoW90pOAowY2/+bd5oFj8w2YzPmDGUFdVnAVvy7SdIg8l2tE0GwYy9+eHZe/H3t8/bjX2Cuy8m/e93uq3uhpUuoNuYMRZYwfcXjQCOAE4D+syY7l7p1GwP8Hy/r98/wDbZhW375M5L325tzGRo+6RrKDhaKL+LrWD75vAIYIUZH3Tn9R386iwzG5Vv9zSzxm6jfTI06qq01okM3IceAZwwwPeWuvvN7n4z73RFtgBvr0jWAzw3wDbZuUbuk66hFkdr7epaiCrXSiwFxuXbo4CV+faOtsnAGrlPuoZaHK21scr381mRXlKzuMfMJpMCYW4+uk/etsXd1+5oW8MfQedp2D7J358JjPv5U+Q/v63uzN1L1zBoZrYKwN0nla1kaHJ/+kl23DTugwH709IknbpPmv1aUYujhfITcDrbD/zqA6bX8Qlad9onQ6PgaLF8eu+DTDvzYQ6871rgM6R3NZ32K2TbPpn6hTUcfPf1aJ/skoKjAHdeZ/xV6zj1mPvcuVHvauW58zoTrnyLU467Uvtk1xQc5WxiF8ObpeVGsIvrhyRRcJSzCTigdBGSWLQ9gP2Bp0vXUgcKjnLWA6NLFyHbfAj4Lw+uLsogKDjKWQt8JL/TSXma/7UCBUcheSGgzaQL3aQ8BUcFCo6yNEFu+1BwVKDgKGs1cHTpIrqdRduLtPC3hugPkoKjrDuBaZogt7hjgEc9eNdd5TpUCo6yHgb2AI4sXUiXm0FaP1YGScFRkAd3tOBxUbm1p+CoSMFR3nLSRVZSxlHAW8CjpQupEwVHefcCR1i0gwvX0a1mActy608GScFRmAd/A/gOeSIYaZ18NmUOsKR0LXWj4GgPVwNzLNqQpuSXITsBWO/B15cupG4UHG3Agz9GOsPye6Vr6TKnAwtLF1FHCo72sRA4U2M6WsOijQYOBb5XupY6UnC0jxWkBZS0+lprBOBvPfjW0oXUkYKjTXjwN4FzgYssmvZLE1m08cAE4MrStdSVnqDt5RbA6bAFittJ7gpeBEQP/mrpeupKwdFGPPhbwDnAVy3anqXr6VDHAyOBb5UupM4UHG3Gg68kLQL0xdK1dBqL9l7g68D83DWUIVJwtKd5wJcs2odLF9JhIvBjYFnpQupOwdGGPPhTwPnAty2a1vdtAIt2NHAycLqGl+8+BUf7ugZ4GTirdCF1l7so3wb+zINvLlxOR1BwtKl8oHQOqcsyvnQ9NXcZ6erX75YupFMoONpY7rLMAW61aCMLl1NLFm0ecCwwR12UxlFwtDkPvpw0HH2ZRdundD11YtE+TjogOsOD/7R0PZ1EwVEPC4AngMW6lmVw8vwmNwJ/6MEfL1tN51Fw1EBuYp9KWm1sgcJj5yza/sAdwAIPflfpejqRgqMmPPj/AVOBaaQLtGQHLNovAyuB6z34FaXr6VQKjhrx4M8BU4Dft2gXqOXxbrmlcQ/pUvkLCpfT0RQcNZPHIUwiTXB8mdaeTfIxjXuBm4DzdAaluRQcNeTBnwV+BxgDLLdo+xYuqah89uRfgKs8+FcUGs2n4KgpD/486UrPp4AfWbTDylZUhkX7U9LArtkeXPNrtIiCo8Y8+FYPfgZwOXC/RZtWuqZWsWj7WLSrgc8DE/NVxdIiCo4O4MEXkyb/WWjRlli095WuqZks2gRgHbAfMMGDbyhcUtdRcHQID34f8JukVckesWjHFy6p4XIr42uky+LP8+AnefAXS9fVjcxrdBzJzFYBuPukspW0N4s2Bfgm6dTkuR78vwuXtFvyaeeppC7Zw8AZ+QCxDKDZrxW1ODpQHi05GngaeMiiXWLR9itc1pBYtI8Cq0hXuP5lbmUoNApTcHQoD/6iBz+H1H3ZF3jMos23aD2FSxsUizbGoi0DlgLXAqM9uNZAaRPqqnQJi3Y4aaj6NNIgqYUe/MGyVb1bXgJzJmmFtYNIXZOFmo28uma/VhQcXcai/Rppjo95wE9Ix0JWlGr+5+MXY4CTSBfyPQQsyjVpQuEhUnD0o+BonDyX6SeB2aTV4x4BluePx5o5+jK3LI4FZpCGzr9GOlOyxIP/Z7P+bjdRcPSj4GiOPCfnJN55IQ8D1uSP1flz31DCJK8PcyQwtt/HUcCDtCioupGCox8FR/PlrsOBvPuFPhYYDvQBm/LnPuAVYCtp7MgwYE/S+rcjSIsejSAN0trIO0G0BljnwV9q2YPqQgqOfhQc5Vi0X+TdgTAS2IcUGHuQAuRN4AXeCZdNwDM6VtF6zX6taM0OGRQP/gqwIX9Il6tbcBwKDH87TUVkQGNI6/I0Rd2CQyMGRQbnZZr4eqnVMQ4RaQ8aci4ilSk4RKQyBYeIVKbgEJHKFBwiUpmCQ0QqU3CISGUKDhGpTMEhIpUpOESkMgWHiFSm4BCRyhQcIlKZgkNEKlNwiEhlCg4RqUzBISKVKThEpDIFh4hUpuAQkcoUHCJSmYJDRCpTcIhIZQoOEalMwSEilSk4RKQyBYeIVKbgEJHK/h9pHcWbjpLx8wAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiAAAACtCAYAAACJOMQzAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAFtBJREFUeJzt3XmUZkV5x/HvA8MmIOBC1ERlFFGQIAKJcSFGA0gUXDiSiCyyw8yIKyJRw/ViUEEwoGHQQRBB1KMoohEVMCpRTkS2KLiBgkHEoGLCIqtW/qjbMz3D7N1ddd/3fj/n9JlzWuz36ep6b/3eqrp1I6WEJElSSWvULkCSJA2PAUSSJBVnAJEkScUZQCRJUnEGEEmSVJwBRJIkFWcAkSRJxRlAJElScQYQSZJUnAFEkiQVZwCRJEnFzapdgDQkEfFSYEPgxtq1qKrZwJ0ppS/ULkSqxQAilbUhsH7tIlSdfUCDZwCRyroRIKV0We1CVE9E1C5Bqs49IJIkqTgDiCRJKs4AIkmSijOASJKk4gwgkiSpOAOIJEkqzgAiSZKKM4BIkqTiDCCSJKk4A4gkSSrOACJJkoozgEiSpOIMIJIkqTgDiCRJKs4AIkmSijOASJKk4gwgkiSpOAOIJEkqzgAiSZKKM4BIkqTiDCCSJKk4A4gkSSrOACJJkoozgEiSpOIMIJIkqTgDiCRJKs4AIkmSijOASJKk4gwgkiSpuFm1C9D4izYCWAdYC7g3NemByiWpkmhjFrAu8AdyX0iVS1IF3TVhXWBNcj94sHJJqiCS739Nk2hjA+CZwPbADt2/TyRfaB7ovtbt/vO7geuAKyd9/XDcL0QR8RyAlNJltWuZKdHGGsBTyH//ia9tgA3Js673kgeetYH7gF+S//5XdP9elZr0u/KVlzOEfgAQbWwCbMfi14THkT+Q3E8OousCfwTuBL7H4teE61OT/li+cpVgANGURBvbAocBfwM8AbiWRQPJlcANwD2TLyLdp+CHkweliQFqB/KF6XvA54EzU5N+U+wXKWRcB55o4+HAPsCe5AHndhYfSK4Bfgc8ODHr0QWVdYDNWHyA2ha4DfgP4HTgsnGbKRnjfhDAc4FDgOcBm5L/9pOvCTcB901cE7r/zyxgE/LffvI1YRPgKuAzwMdTk+4o+OtohhlAtMqijXWAVwJzyTMcC4ALgB9MZXkl2tiIfNHZB3g58EVgPvCdcRmAxm3giTb+HJgDvAr4GnAW+e+12uEx2lgTeCrwou5n30PuB+emJt011Zr7YAz7wQbk9+1c8ozGacBXgJ+kJv1hCj/3UcCzgAOAFwKfAuanJl075aJVnQFEKy3a+FPgtcCBwH+RB4V/m4llk2jjkcD+5AHoju61zk5Nun+6X6ukcRh4upmLPYF5wJPJAfQjqUm3zNBrvZA8sD0fOBf4QGrSDdP9WiWNQz8AiDaeArwO2Bv4Bvl9+u8zsWzSXX8OAQ4lz6yeCnzGJZrRZQDRCnVTpAcB7yEPAPNTk35S6LXXAHYG3kyezt0/NemaEq89E0Z94Ik2NgfOJC+dnAB8odSm4mjj8eTB53DgROCkUd0zNAb9YBZwJPl9+SFgQWrSzYVeey3gpcBR5P1EB6Ym/bTEa2t6GUC0XNHGE8jr8I8kD/5Vpj67ELQveeA5DThuFGdDRnXg6YLga4FjgOPIsxCrPbU+xVo2Az5C3tS6f2rSD2vUMRWj2g8Aoo2tyEttdwAHpSb9vFIdawKvB94GtMCpzoaMFgOIlmqJWY+TgRP6cPtstPE48pT/48mDz9WVS1olozjwTJr1WBM4oNTs1/J0/fMw4F2M4GzIiPaDybMe7yDPelQfQKKNp5L754M4GzJSDCB6iGhjQ/Jmr8eQB/nvVy5pMUvMhhyfmnRS5ZJW2qgNPNHG3sApwLuBU2rNeixLNxtyBrA+8LLUpP+pWtBKGsF+8CfkjeZ3AQenJt1Ut6LFLTEb8vrUpHMrl6SVYADRYrrNnxeSb4ed24dZj2Xp9gR8lXxhfFsfPo2tyCgNPNHGEeR19r/r810H3fJQA+wF7FxrSWBVjFg/eCJwCXn/V9vn91m0sTXwZeC9qUmn1q5Hy2cA0ULdp5yvkQPIW/t8oZnQ3ab3FeByYF7fax6VgSfaeBv5bqedU5NurF3Pyog2Xg+8CdgpNen62vUszwj1gy3I4eOk1KRTatezMqKN2cDFwBmpSe+pXY+WzWfBCFh4YuFFwGcZkfAB0J038ULyCawndsszmoJo4w3Aa4AdRyV8AHQD5LuAi7vZMU1BtwH9YuDYUQkfAF2f3RE4oAul6ikDiIg21ge+RJ79eOeohI8J3emILybfrvv2yuWMtGjjAOCN5JmPW2vXs6pSkz4CfIAcQjatXc+o6truYuBfujYdKV3f3Ql4U7Sxf+VytAwGEEE+POgG4M2jFj4mdM8O2QU4KNp4ae16RlG0sQPwXmCX1KT/rl3P6kpNej/wOeCT3f4QrYKuzT4FnJeadHLtelZX14d3AY7v+rZ6xjfnwEUbu5Of2TBnVMPHhNSkX5FPTz0t2nhE5XJGSne8/lnAG1KTfly5nOlwDLAB+eAyrZrDgfXIbTjSur78RuCsro+rR9yEOoMi4ijyAV7fJR/m9VZgY+DJKaXDatYGC/d9fB/YJzXpG5XLmTbRxgeBjVKT9qtdC0BEfIncD75F3iR5EfkpoLNTSlvVrG1CtHEcsBWwx6gH0QnRxpbApcBf9OG20RHpB7PJ16vnpSb9qHY906HbF3Y+cF1qUi+WaPs+NpTiDMjM+llK6a0ppfOA21NKC1JKJ5DXVvvgZOD8cQofnaOB53azO33ww5TSX6WUjiQ/dvzYlNJLyHcbVddNTx/MGMyCTdadkHoicEZPNif3vR+sQT5T5fhxCR8AXZ8+HDi4R0sxfR8bijCAzKyrVvH7xUxaejm6di3TLTXpbvItpKd1szy1fXEVv1/MEksvv6pczkw4ibwU04dPlb3tB52JpZf31y5kunV9u09LMb0dG0oygMyglNLPVuX7pUQbE4/LPqgbrMdOatI3ydOu1c8BSCl9c1W+X9iRwPXkTYdjpzuefX/gXbXviulzP+ja5ljyUfu9Ou12Gn2SvNn+yNqF9HVsKM0AMkx7AteO4dLLko4F/sENqUvXPVV0HvD2cVp6WVK3FHMB+dlGWrqDgc+N09LLkro+/nZgbtf3VZkBZJjmkm+9HWupSb8G/o38CVgP9XLgx6lJP6hdSAHzgcO7Z4Zokq5NDmcY14TryLMgL6tdiwwggxNtbAc8jnzw2BDMB+Z4HsRSDSKIAqQmXQXcSj6wTot7CfCL1KRrahdSyHxy31dlXpRnWERsHBGHAk+KiKMi4kmVS5oDfHiM13mX9J/kJ3juVLOIiHhiRJxDfqT9+yLi+VXraWMr4KnA52vWUVj1gadv/aAzmCDaOR/YsrtNu5oejg3FeQ7IgEQbGwM3Ak8blceWT4do4xBgt9Sk6tOufXkIWbTxr8BvU5OamnWU1G2+vhl4dmrSDVVr6U8/eArwbeAJqUn31qylpGjjXcDGqUlH1K5lyJwBGZY9gYuGFD46nwCeH208unYhfRBtzAL2Jh+ANBjdAHs2+UF7yvYDPjak8NFZAOzjnqC6DCDD8lfA12sXUVp3q/GVQF8OIartacBtqUm/qF1IBV8HnlW7iB4Z6jXhZuA35PeCKjGADMv25IF4iAwgiwy+H/TkZNSqujYYfF+oXcSQGUAGItpYD9iC/OyXIbqCfLFVvugOctDpHtN+D7BZ5VL6YDbw+wEuyU7wmlCZAWQ4tiGf+TC0td4JV+LFZsKQP/WCfWGC/cB+UJUBZDiGfrH5GbBB7eO4a+s2oG4DXF27looceLKhXxOuAp7hRtR6DCADEME6fOuoV/KJCzaNYK8I+vAwprLemdbm8rm3cMa33jfUNohgHc785hu55Lh7eWd68WDb4LPnrM+XT37lUPsBdO3wnXkv5sxLNx9sO7wz3ctVB9zJB3/8/sG2QWWeAzLmItie/LTNx0769q3A7ikN49NPn9qg1vkPfWqDWvrUBjXPAelTO9RiG/SDAWSMdYn+RhZ/k024FZidEveVraqsvrVBjYGnb21QQ9/aoGIQ7VU71GAb9Mes2gUUEXEysG3tMko7lTmbzmP+0t5kkN98r2BMH8M+yR4s/UID8Nhz2Od24twHShVzaT6CGyKKHYV/NvustR/nPGwZ/7P9YCD9AOwLneX2hfnMuZr40G0lC+qBa0jpDaVf1D0gY+wGNl9vBf/Jk4sUUtdyf8ef8uSxfw+sxO9oPxhAPwD7Qme5v+MNbL5uqUKGziWYMRbBXuRjyJdlr5TG+9NO39qg0hJMr9qghr61QcUlmF61Qw22QX8YQMaYa539awP3gNTRtzZwD0g9tkF/DGLacai6N9Hu5DfVZBO7vcf+TWYb2AZgG0ywHWyDPjGAjLnulrLZPOd9Z/PMM74D7EVO+IO51SwlrmSjm2bzin3vYMNbTmCobQCzWevufXjBMQ+w3m/2Z6ht8Lz3fIpnnPVtBtgPYFI7vOAdF/L0T1/EANshJa5k90Oezh6vvp8172sYYBv0gUswAxFt7AicmJo0yCeBRhuPJz/74TGpqdfpa57/sLCGNq4AXpeaejXUFG2cB3wuNWl5+wBmtoZ+9IO9gZenJu1Zq4aaoo3nAienJv1F7VqGyhmQ4bga2DraWKt2IZXsAFxZM3z0yNAfwjXYh/EtwX5gP6jKADIQqUl3AT8HtqpdSyXbky+4GvCzUKKNRwKbANfXrqUHrgceFW08onYhlXhNqMwAMiyDHXjwwVuTDb0fXJ2a9MfahdTWtcHVwHa1a6nEa0JlBpBhuRzYsXYRpUUb6wB/iZ92JlwLbBZtPKp2IRXsCHy3dhE9MtRrwqOBJwDX1a5lyAwgw/IZ4OXRxsa1CylsD+Ca1KRbahfSB6lJ9wOfBQ6oXUtJ3f6nA4CP166lR84FDow2hvFYjkUOBM7r3guqxAAyIKlJvwK+AuxXu5bC5gLzaxfRM/OBOdHGkK4BuwM3pSb9V+1C+iI16Rrgv4HdatdSSrSxJnA4XhOqG9LFR9l8YG60EbULKSHa2AZ4EvCF2rX0zHeB24EX1S6kIIPo0s0nt81Q7Ar8JjXJpbjKDCDD8y3gAeAFtQspZA6wIDWp2JNOR0F3O/JgBp5o42nA1uSlJy3uPOAZ0cYWtQspxCDaEwaQgRnSwBNtPBx4FXB67Vp66lPAs6ONzSrXUcLhwBmpSR6zvYSuTc4gt9FYizZmkzek+7C5HjCADNPHgedHG9vWLmSGvQn4SmrSL2sX0kepSb8HzgT+qXYtM6k7BXcfYEHtWnrsw8B+0caf1S5khr0D+Ghq0j21C5EBZJBSk+4E3gKcFW2sXbuemdCFq7nAkbVr6bl/BnaKNsZyL0i31+l08pHbP69dT191bfMB4PRx3R8WbewK/C25z6sHDCDD9THgF8Dbahcy3bpQdRZwlLfeLl9q0h3AweSBZ6Pa9cyAA4FHA8fXLmQEvAf4E8bw9uyuby8ADu76vHrAADJQ3V6Qw8h3xIzbUsw/AreQQ5ZWIDXpYuDLwIm1a5lO3dLLe4H93YS8Yl0b7Q8cP4ZLMScBF6YmXVK7EC1iABmwbnZgrJZiujA1DzjUB8+tkrcAO4/LUsykpZdTUpO+X7ueUZGa9D3GbCmmW3rZidzH1SMGEJ1NXoo5rnYhU9Xd9XI2Lr2ssm5a+hDywPPY2vVMgyNw6WV1vZe8FPPa2oVMVdeXTycvvdxZux4tzgAycN0swQHAbtHGyG7YjDbWIx82dhkuvayWbinmw8BFo/yE1GhjL+Ao4O9dell1XZvtCRwdbbyqdj2rq3vy8cXAh1x66ScDiEhN+jWwMzAv2ji0dj2rqls++jTwS2CeSy9T8m7ycf1f7maURkq0sRvwL8CuqUk/rV3PqOrablfglK5NR0rXdy/svt5duRwtgwFEAKQm/QLYBXhHtPG62vWsrGjjYcAFwIPAa1KT/lC5pJHWhbejyI8pv6T7FDkSoo09yQdqvTQ16dra9Yy6bu/M7sAZXduOhO4pz18jP/36rX4g6S8DiBZKTboe+GvgiGjjmL5vQuturfsq8GtgT6fbp0d3wZ4H/DtwabTxp5VLWqFo42DgFGCX1KTLa9czLrq2fBF5JuSg2vWsSNdXvwlcArzW8NFvkfz7aAnRxmOALwG3ku8m6d1JotHGC8mfdj8PvDk16Y+VS1opEfEcgJTSZbVrWRnRxlHkE2WPSE36TO16lhRtbAKcDDwX+LsuRPfeCPaDLci3av8H8MbUpN9VLukhulmaDwLvT006oXY9WjEDiJaq21fxdvLD3N4CnN2HTxPRxobACeTHhx+WmnRh5ZJWyagNPADRxrOBjwLfI++x+XXlkgCINl5C3jT7eeDo1KS7Kpe00ka0H2xAvqvoZeT33pcqlwRAtLEpcCrw5+QzX/6zcklaSQYQLVe08UzyqaI3U3k2ZNKsx9eBN6Um/W+tWlbXKA48sPAuo2OBfak8GzJp1mNH4KDUpK/XqmV1jWo/AIg2XkB+H15K5dmQSbMeZwONz3gZLQYQrdASsyEnAGemJt1e8PW3Ji8D7MwIznpMNsoDDyw2G/ID8qfhy0vNjHWfwF8NHMMIznpMNgb9YPJsyLHAJ0r9Lbq9aX8JHA1sibMeI8sAopUWbWxDXo7ZDTgfmJ+adMUMvdbawCvID5TbnEUPFBu5WY/JRn3ggYWzIfPIgfR/gfnAJ7un687E623Zvdbe5A2GJ6UmfXsmXquUcegHANHG84A3kzevn0u+Jvxohl7rYcBe5GvCxsBpwKnOeowuA4hWWbTxaPJDvuYAt5EHoC+mJv12ij83yGFjX/KpnD/sfvYF43KHy7gMPADRxhrkW7fnAs8BziHPjlw31duhu3McJn72lsBHgAWpSTdPqeieGKd+AAufu3MY+cGG15HDwUVTffBbtLEm8HTyYYn7kg8anN/97JHYeK5lM4BotXUXh12Bw8mfgG4nnx+x8GtZoaQLG08CdgC27762A+4iT6+flpr0g5n+HUobt4FnQrSxGXAo8PfkY7yvIfeBK7p/f7KsUNKFje1Y1A92AB4HfIf8BNPzU5Pun9nfoKwx7gdrA3uQP0A8i/xQyMnXhKuWFUq668kWLH5N2Bb4H/JBgwtSk26a4V9BBRlANC26T8Obs2gA2R54JvAH4B7gXuABYN3ua0Pgtzw0sNxWvPiCxnXgmazbJDo5UGxPDhR3kvvBPcAscj9YD1iHfIfN5MDyo9SkB4sXX8hA+sEs4Gksfk3YBriPRdeEB8l9YOKa8EseGlh6d8uvpocBRDOmCyWbsGigWYtFA9Dvpzo9O4qGMPAsTbSxPrABuR+sRx547um+/m+cw8bSDLgfzAI2YlE/mMWifnBXatLdFctTYQYQqaChDjxanP1A8ih2SZJUgQFEkiQVZwCRJEnFGUAkSVJxBhBJklScAUSSJBVnAJEkScUZQCRJUnEGEEmSVJwBRJIkFWcAkSRJxRlAJElScQYQSZJUnAFEkiQVZwCRJEnFGUAkSVJxBhBJklScAUSSJBVnAJEkScUZQCRJUnEGEEmSVJwBRJIkFWcAkSRJxRlAJElScQYQSZJUnAFEkiQVZwCRJEnFGUAkSVJxBhBJklTcrNoFSAMzG1g/ImrXobq2Bu4GLqtdiFSLAUQq687aBagX7sa+oIGLlFLtGiRJ0sC4B0SSJBVnAJEkScUZQCRJUnEGEEmSVJwBRJIkFWcAkSRJxRlAJElScQYQSZJUnAFEkiQVZwCRJEnFGUAkSVJxBhBJklScAUSSJBVnAJEkScUZQCRJUnEGEEmSVJwBRJIkFWcAkSRJxRlAJElScQYQSZJUnAFEkiQVZwCRJEnFGUAkSVJxBhBJklScAUSSJBVnAJEkScUZQCRJUnEGEEmSVNz/A+vGnVPGuDqTAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "H = SSH_lattice.display_unit_cell(label_on = True)\n", "T = SSH_lattice.display_lattice()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Hamiltonian of the lattice can be obtained with the method Hamiltonian()" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "Quantum object: dims = [[3, 2], [3, 2]], shape = (6, 6), type = oper, isherm = True\\begin{equation*}\\left(\\begin{array}{*{11}c}0.0 & -0.500 & 0.0 & 0.0 & 0.0 & -0.500\\\\-0.500 & 0.0 & -0.500 & 0.0 & 0.0 & 0.0\\\\0.0 & -0.500 & 0.0 & -0.500 & 0.0 & 0.0\\\\0.0 & 0.0 & -0.500 & 0.0 & -0.500 & 0.0\\\\0.0 & 0.0 & 0.0 & -0.500 & 0.0 & -0.500\\\\-0.500 & 0.0 & 0.0 & 0.0 & -0.500 & 0.0\\\\\\end{array}\\right)\\end{equation*}" ], "text/plain": [ "Quantum object: dims = [[3, 2], [3, 2]], shape = (6, 6), type = oper, isherm = True\n", "Qobj data =\n", "[[ 0. -0.5 0. 0. 0. -0.5]\n", " [-0.5 0. -0.5 0. 0. 0. ]\n", " [ 0. -0.5 0. -0.5 0. 0. ]\n", " [ 0. 0. -0.5 0. -0.5 0. ]\n", " [ 0. 0. 0. -0.5 0. -0.5]\n", " [-0.5 0. 0. 0. -0.5 0. ]]" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "SSH_Haml = SSH_lattice.Hamiltonian()\n", "SSH_Haml" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Sublattice Projectors and Chiral Symmetry of the SSH model:\n", "We explain now that, the Hamiltonian of the SSH model is bipartite. We can set the sites into two partitions. Each partition consists of the set of every other sites starting from site 0(or 1). Since, the Hamiltonian does not enable any transition between sites of the same sublattice, only between them, the Hamiltonian is said to be bipartite." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\\begin{eqnarray}\n", "\\hat{P}_{A} = \\sum\\limits_{m=1}^{N} |m,A\\rangle \\langle n, A|,\\ \\ \\ \\ \\ \\ \\hat{P}_{B} = \\sum\\limits_{m=1}^{N} |m,B\\rangle \\langle n, B|, \\ \\ \\ \\ \\ \\\n", "\\hat{\\Sigma}_z = \\hat{P}_{A} - \\hat{P}_{B}\n", "\\end{eqnarray}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The projectors $\\hat{P}_A$(or $\\hat{P}_B$) project out a ket vector into its projection into sublattice A(or B). Because, the Hamiltonian is bipartite,\n", "\\begin{eqnarray}\n", "\\hat{P}_{A} H \\hat{P}_{A} = \\hat{P}_{B} H \\hat{P}_{B} = 0\n", "\\end{eqnarray}\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The SSH Hamiltonian is said to have chiral symmetry, since\n", "\\begin{eqnarray}\n", "\\hat{\\Sigma}_z H \\hat{\\Sigma}_z = -H\n", "\\end{eqnarray}\n", "which is a generalization of the symmetry concept of a Hamiltonian, since, ordinarily, a Hamiltonian is said to have a symmetry if an unitary operator leaves it invariant.\n", "\\begin{eqnarray}\n", "\\hat{U} H \\hat{U}^{\\dagger} = H\n", "\\end{eqnarray}\n", "Now, we form the operators $\\hat{P}_A$,$\\hat{P}_B$ and $\\hat{\\Sigma}_z$, and verify the properties of the SSH Hamiltonian." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "True\n" ] } ], "source": [ "chiral_op = SSH_lattice.distribute_operator(sigmaz()) \n", "anti_commutator_chi_H = chiral_op * SSH_Haml + SSH_Haml * chiral_op\n", "is_null = (np.abs(anti_commutator_chi_H) < 1E-10).all()\n", "print(is_null)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence, it is verified that $\\hat{\\Sigma}_z$ and H indeed anticommute and the SSH Hamiltonian has chiral symmetry." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The dispersion relationship for the lattice can be obtained with plot_dispersion() method." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAY4AAAEJCAYAAACDscAcAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAIABJREFUeJzt3Xd8VFX+//HXARJCgBDS6CEkgFSFMEizEyxYQA2gImAjqFhW3UW/uv521WW/i7q733VtxIZiVIqIvQUrUpMgAlITeg0JCYGQkGTO7497JwyQMikz907yeT4eeWTm3rnMx3Fm3jnn3HuO0lojhBBCeKqJ1QUIIYTwLxIcQgghakSCQwghRI1IcAghhKgRCQ4hhBA1IsEhhBCiRiQ4hBBC1IgEhxBCiBqR4BBCCFEjzawuwBsiIiJ0TEyM1WUIIYTfSE9PP6y1jvTksQ0yOGJiYkhLS7O6DCGE8BtKqZ2ePla6qoQQQtSIBIcQQogakeAQQghRI5YFh1JqVhX7EpVSCUqppKq2CSGE8D1LgsP88k+sZF88gNY61XW/om0+KlUIIcQZLAkOrXUykFXJ7glAnnk7C0ioZJsQ/iMlBWJioEkT43dKitUVCVFrdjwdNxTIdbsfXsk2IfxDSgokJUFhoXF/507jPsDEidbVJUQt2TE4hPA/zjIozDn1U3QUTh4zfv74p1Oh4VJYCH98APo1g8CWEBQCweGnfpoGWPPfIYQH7BgceUCYeTsUyDFvV7StnDlukgQQHR3t5RJFo3QiD7I3Qc42yM0yfvL3wNF9UHAAdFnFxx04Wsn2XPioonM9FLSKgpCO0KYzhMWaP3EQ1QdaSoNbWMs2waGUCtVa5wHzAIe5ORZINW9XtK2cOW6SDOBwOLR3qxUNXsEB2JsB+3+FfWvg4AY4uvfUftUUQqOhbVeIvcT4km/VDoLDjBZDUBsIbG20Jt48H3bvOfs5unSC+380WiVF+VCYa7RWjh2Cgn2QvxcObYItX0PZyVPHtYyCdn2gwwDoOBA6xUObLqCUt18VIQCLgkMplQg4lFJJ5hc+wBJgkNY6QynlUEolAHla6wzzmLO2CVFvcrfD9h9h53LYvQKO7DC2qyYQ2QtiLoCo3hDZGyJ7Gl/UnnYn/e8/Th/jAAgOhv+dBeFx1R/vLDNC6/BWo8Vz6Hc4sA6WvwTOEuMxrTtC9BCIHgbdLobIcyRIhNcorRveH+cOh0PLXFWiSsXHIOsH2PoNZH0PebuM7S0jIXoodBkKnQdD+/4QGFz350tJgSeegF27IDoaZs6s+8B4abHREtqbDrtWwK7lp1pFrdpD7MXQ43KIu8xoCQlRBaVUutbaUf0jJThEY3LsEGz6DDZ+CjuWGt0/zUOg20XGX+mxF0NET//+S/3IDsj60Wg9Zf1gdH2pJkYQ9r4Gel9rdLEJcQYJDgkO4VKYCxs+gvUfws5lgDYGms8ZDT2vMLp2GuoZTM4yY5xm69ew+Us4uN7Y3nEg9LvR+AnpaG2NwjYkOCQ4GrfSk7DlK/g1BbalgrMUIs6BvtdDn+uMM5P8uVVRWzmZRmvr98XGgD/KGLsZMNF4XQJbWl2hsJAEhwRH45S9BdLnwG8fGF00rTtA/3Fw7nho169xhkVlDm+D9Qvht3nGacWBraHfDTDoNuMsLdHoSHBIcDQeZaXGuMXq12HHz9AkAHqNhoGTjEHhJk2trtDetDYG1de8a3TplRRCx3gYfKfRlRXQwuoKhY9IcEhwNHwnjkDGO7AyGY7ugTbR4LjdCIxWHq1+Kc5UlA9r50HaG8Zpv8Hh4LgTBt8FrdtZXZ3wMgkOCY6GK2+3cf1CxjtQctw4I2rIPcZAt7Qu6ofWxllnK16BzV9Ak2Zw7gQY8aBxDYtokGoSHLa5clyIKmVvgaX/hnXzjfv9x8Gw6cZ1FqJ+KQXdLjR+cjJh5auQMdc42aDX1XDhIzIO0shJi0PYW/Zm+PFZ43TagBYQP8UIjNAuVlfWuBw/DCtnw6pkKMozLiy8+DHoPMjqykQ9ka4qCQ7/l5MJP/wD1i2AgGA4/y4Y/gC0jLC6ssat6KgRHstfNMaZel4Jl/1ZWn4NgASHBIf/OroffpwFa+YaZ0gNSYLhD8qMsHZTXGAEyC//MQbV+95gBIgnc28JW5LgkODwP8UFxpfQsheNC/Yctxt96a3bW12ZqMqJPKP1sfxlKCsGxx1w8aPSMvRDEhwSHP7DWQYZb8N3M6HwsHHtwMj/B21jrK5M1MSxQ0bXYvoco2vxwodg6HQICLK6MuGhmgSHJWuOCwHA9p9g9kXw2UPG5IJTv4PENyU0/FGrKLjmX3DvCuNsrCVPw0vnG1OcNMA/Ths7CQ7he/l7YP4UePtaY7B13Ntw+xfQSc7Q8XuRPeHm92HSYqPlMe9WeGeMcXacaDAkOITvlJ40rsV4cbCxqt2lT8B9q6DvWJlHqqGJuxTuXgqjnzdWUXxlOHz7/4x1UITfkwsAhW/sXAaf/gEOb4Ze18AVfzeWXRUNV9NmcP5U6DMWlvzVOPlh3Ycw+jljPjHhtyxpcSilEpVSCUqppAr2xSultFIq0/yZbW6fZf4+6xhhYyeOwCf3w1tXQckJuGU+3JQiodGYtIqEMS/Bnd8aa7F/cLPRhXV0n9WViVryeXAopeIBtNap7vfdhGmtldY6DhgHzDK3JymlMoEsnxUr6ub3j+HF82FNinHx3vQVxpxSonHqcj5M+xES/gpbU+GlIZD2FjidVlcmasiKFscEIM+8nQUkuO90BYrJobV2BcVUrXXcGfuFHRUcNP6inD8ZQjpA0g9w+TOyUJAwVlu84CG4dzl0HACf/QHeuc6YKUD4DSuCIxTIdbtf4SXBSqkEYL7bpjCze2uGN4sTdaA1rFsILw+BLd/AyL/AXd9Bh3OtrkzYTVg3mPwJXPsC7F8Lr4ww5sKS1odfsPNZVaO01q6WCVrrZLO1EW6GymmUUklKqTSlVFp2drZPCxUYk+DNnwwf3glhccYZNRc+bAyQClERpWDQFOPaj5gR8OUMo/VxZIfVlYlqWBEceUCYeTsUyKnkceVjH2YoJJp3c4DYMx9sBotDa+2IjJSFfHxq85dGf/WWr4z+6zu+lnUbhOfadIKJC+G6/8K+X43WR8ZcuXDQxqwIjnmc+uKPBVyD5KGuByilzgyGLNfjMLq2ZD4ROyg+Bp88AO/fZKzvnfSj0X8trQxRU0pB/GS4dxl0HAif3GeMkx0/bHVlogI+Dw6tdQaUj2Hkue4DS854aJbbMalAgtnqyHE7RlhlbzrMvtBYiW/EgzB1CbTrY3VVwt+FRhtjH6Oega3fwMvDYNuZXw3CajLJoagZpxOWvQDfPQOt2sENyRBzgdVViYbowHr48C7I3gjD7jNOtmgWaHVVDZZMcii8o+AgvHs9pP7FWEL0nl8kNIT3tO8HSd/D4LuMqdvfSJDTdm1CgkN4JvM7eHUE7FppnEI57m1o0dbqqkRDF9ACrv4n3PQe5O2C2Rcbp3wLS0lwiKqVlRpTZM+9AYIjjL8AB02RSQmFb/W62jjFu11f45TvT+43prARlpDgEJUrOAhzx8LP/4SBtxrrZUT1troq0Vi16Qy3fQ4XPGyclPH6KOm6sogEh6jYjqXGWVN70mDsqzDmRQgMtroq0dg1bQYJfzGu+zi6x+i6+v1jq6tqdCQ4xOm0hl9egLevg+atjdNsB9xsdVVCnK7HKJj2s3Gh6fzJ8M2fjW5V4RMSHOKU4gJYMAW+fdLoU576vdGnLIQdhXaB27+CwVNh2X+NbtVjh6yuqlGQ4BCGw1vhtZHGGtGjnoHx70BQiNVVCVG1ZoFw9fNGd+qe1cYa9nvkGi5vk+AQsPkreO0yKDxsrBU94gE5a0r4lwE3GwtFNQ00Fg3LmGt1RQ2aBEdjpjX8+Jwx11TbGGPdjNiLLS5KiFrqcK7xHu463Jjr6vNHoKzE6qoaJAmOxurkcVhwG3z/N+g/zpjRNjTa6qqEqJvgMJj4IQy/H1a/DnOvh+OVTcAtakuCozHK2w1vXmGcxjjqGWO+KTnVVjQUTZvB5X+D65Nh9yp47VI4+LvVVTUoEhyNjeuDdGQn3DJfxjNEw3XeBLj9SygthjdGGevGiHohwdGY/DYf5lwNga3grlToebnVFQnhXZ0HGdPkhHeH9282rlFqgDOC+5oER2PgdMKSZ2DRVOh8vjF1SOQ5VlclhG+EdDRaHn2uM65R+uQ+KD1pdVV+TYKjoSs5AR/eAT8/b6ywNukjYwBRiMYkMBgS58BFM2DNu/DuDXDiiNVV+S1LgkMplaiUSlBKJVWyf5b5O8nTY0QFjmXD29fChsXGIPi1L8hCOKLxatIELnvCHDRfaUySmJtV/XHiLD4PDqVUPJQvB1t+/wxJSqlMzOVjPTxGAKSkQEyM8SGJ7gRfrzauApdBcCEM502AyR8bF7y+NhJeeOrUZyYmxvgMiSpZ0eKYAOSZt7OAhAoeM1VrHecKCg+PESkpkJQEO3caA4BHSuDzk7CmwOrKhLCXrsPhriWwwQl//Oupz8zOncZnSMKjSlYERyiQ63Y/vILHhJndUjNqcIx44gkoLDx924kiY7sQ4nThcfCDE868uLywUD4z1bDl4LjWOtlsbYQrpTxqXSilkpRSaUqptOzsbC9XaENaw66dFe/btcu3tQjhL/bsrXi7fGaqZEVw5AGu03pCgdPmAzADING8mwPEVncMlIeNQ2vtiIyM9ErhtuUsgy9nQEglYxjRMpWIEBWq7LPRpbNv6/AzVgTHPIwwwPztGvAONbdlubZhdEmlVXaMwDjddsEUWJUMd18HwWdMHRIcDDNnWlObEHY3c+bZn5kA4PJQKMyt8BBhQXBorTMAzC6oPNd9YIm5PxVIMFsdOVrrjCqOadxOHDEmcdv4GVz5D/jHYkhOhq5djTOounY17k+caHWlQtjTxIlnf2aeeRCi9xvTs+fvsbpCW1K6AV5+73A4dFpaA1/MJX8vvHsj5GbC9bOh3w1WVyREw7H9J/hgorF88q2LIKqX1RV5nVIqXWvt8OSxthwcF9XI3gxvXG78NTRxoYSGEPWt20Vw2+fgLDVmkt61wuqKbEWCw9/sSYc3r4Syk3D757LwkhDe0uFcuPMbCA6Hd8bClq+trsg2JDj8SeZ3xhQiQSFw59fQ4TyrKxKiYWsbYyxyFtnTmF137TyrK7IFCQ5/sX4RpIyHsFjjjRwWW/0xQoi6axUJUz6DmBHwURKseMXqiiwnweEP0t6ChXdAZwfc9hm0bm91RUI0LkEhxnhi72vhq8fgu5mNel0PCQ67+/lf8NkfoMco4+yOFqHVHyOEqH/NmhtTsw+8FX561rjo1um0uipLNLO6AFEJrSH1L/DLf6D/OBj7CjQNsLoqIRq3ps3guhehRVtY9l84kQdjX250n00JDjtylsHnj0D6W+C4A0b/05jyWQhhPaWM9W2CQuG7Z+DkcUh8EwKCrK7MZ+TbyG7KSuCjaUZojPgDXP0vCQ0h7EYpuOiPcNVzsPlzeG88FB+zuiqfkW8kOykpgvmTYd0CGPkXGPWULL4khJ0NSYKxr8KOpcb0Pyfyqj+mAZDgsIuTx+H9CbD5Cxj9PFz4sNUVCSE8MeBmGP827FsDb18Dxw9bXZHXSXDYQVE+zL3BmB9n7Ctw/lSrKxJC1ETva+GWD+DwNnhrNBzdZ3VFXiXBYbXCXHj7OtibBolvwYBbrK5ICFEb3RPg1g+N0HjrKshruItBSXBY6Vg2zLkGDm2Em96DvmOtrkgIURcxI2Dyx8aSB29eBTmZVlfkFRIcVjm6D+aMhiPbYeJ86HmF1RUJIepD50HGFCWlJ4xuq0ObrK6o3lkSHEqpRKVUglIqqZL9SebPLLdts1z7fFWn1+TtMpqyR/cbTdvYS6yuSAhRnzqcC7d9AWjjD8QD66yuqF75PDiUUvFQvtJf+X23/QlAqtY6GYg17wMkKaUyMZaW9V+5WcZfIYVHYPJi6Drc6oqEEN4Q1Qtu/xKaBRld0nsbzsKlVrQ4JgCuk52zgIQz9se6bcvi1FrjU7XWca7A8UuHt8FbV8PJYzDlE2PSQiFEwxUeB7d/YUyS+M4Y2L3a6orqhRXBEQq4rwIf7r5Ta51stjYA4gHXGrBhZvfWDB/UWP8ObTKarGUnjZXFOg6wuiIhhC+0jTG6rVpGwNyxsHO51RXVmW0Hx80urAytdQaUB0oqEO7WfeX++CSlVJpSKi07O9vX5Vbt4AaYc7Vx+7bPoV1fa+sRQvhWaBfjs9+6A7x7I2z/2eqK6sSK4MgDwszboUBOJY9L0Fo/CuWhkGhuz+FU91U5M1gcWmtHZGRkfddce/t/M/o3mwYYb5xGsOi9EKICIR2N74DQLpAyDjK/t7qiWrMiOOZx6os/FnANkpcvNKGUStJaP2veTsAY63CNbYRzqvvK3vb9aiz1GtDCeMNE9LC6IiGElVq3M07VDYuF92+CbUusrqhWfB4crq4nMxDyXPeBJW7bZymlMpVSR8xjUoEEs9WR43aMfe1Nh3eug+YhxuBYeJzVFQkh7KBVJEz5FMJ7GOuYb/3W6opqTOkGuPyhw+HQaWkWNkr2pBszZbYINd4gbbtaV4sQwp4Kc43B8kMbYfxcOOdKS8tRSqVrrT061dO2g+N+a/dq480Q3NbonpLQEEJUJDjMmJ6kXV+Ydyts/tLqijwmwVGfdq00WhrB4cbpd6FdrK5ICGFnLdrCpMXQvj/MmwSbPre6Io9IcNSXXSvg3RugVZQxptGmk9UVCSH8QYtQYxaJDucZC7n9/onVFVVLgqM+7FxunJvdqh3c9plx2p0QQngqqA1MWgQdB8LC220fHhIcdbVzOaQkQuv2xpiGhIYQojaC2sCt/hEeHgWHUmqetwvxS+6hMeUzCOlgdUVCCH8WFOIX4eFpi2O+UmqgUuoypVSMF+vxH2e1NCQ0hBD1oDw84m0bHp4Gx7da6zVAWyBZKfVHpdRlXqzL3s4Mjdbtra5ICNGQBIUYa/W4Wh4bP7W6otN4GhwZZneV1lpfrrV+Xmv9XaMMj10rJDSEEN7n3m214DZbhYenwfGo1nqC1nqRa4NSagAwyjtl2dSuFcbZUxIaQghfOCs8PrO6IsDz4Gjrfkcp1Qa4CfC/SVZqa9fKU6Ex5TMJDSGEb7i6rToMgAVTbHGRoMfBoZSa5zYw/hgwG2Na9IZv9+pT12nI2VNCCF9zXefR4TyYPwU2fWFpOZ4GR7rWegJu06FrrbcD+d4py0b2pJlXhEeaF/dJaAghLOC6zqN9f+MK881fWVaKp8ExSCl1FxCrlOoGxCmlQoA23ivNBlyz3AaHmy0NubhPCGGhFqEw6SNo3w/mT4ItX1tShqfBkQwc0Vq/DoSaU+9Ow1jNr2Ham2GGRpjR0pC5p4QQduAKj6g+xqy6Fqzn4WlwTMUcCDev50Br/ZzW+jtvFWapfb8aU6O3CDVaGm06W12REEKc0qKtMTFiVG/4YCJsS63+mHrkaXDkaa2Puu6Yp+LWmlIqUSmVoJRK8nR/dcfUm/1r4Z0xRn/ibZ/J1OhCCHtyTcke2RPevwUyffd3vKfBMV4ptdU8s2o+sKC2T6iUiofy5WDL71e1v7pj6s3+34zQaN7aaGmERnvlaYQQol4Eh8HkTyCip7EMbdYPPnlaT4Njtta6B5CktR4P3F2H55zAqbGRLCDBg/3VHVN3JSfgvfEQ0NJoacjKfUIIf+BaSTAszgiPgoNef8pmHj4uTCn1CpCplEoG6rJQeSiQ63Y/3IP91R1TdwEtYMyLEBYLbWPq/Z8XQgivaRkOUz6BHT9D63ZefzpPWxyZWut7gCXuYx12opRKUkqlKaXSsrOza/Vv7A4bTmEr6Z4SQvif0qAwMqN8MwtUTa7juAzoZg6MD6rDc+YBYebtUCDHg/3VHYPWOllr7dBaOyIjI2tcVFFJGTclr+DOOWmcOFlW4+OFEMIqpWVOHpq/lrEv/cKhgiKvP19NruO4HGN+qgSt9XN1eM55uF2BDrgGvEOr2F/hMfUpKKApM648h5Xbc7jrndUUlUh4CCHsr8ypeWTBWj5du4/7Lu1OVOsgrz+nR8Ghtc7XWj+mtR6vtX7evGq8VrTWGQBKqQSM03wzzF1LKttfxTH1asyATvxz/Hksy8xh6jtpEh5CCFsrc2r+uGAtH/+6j0ev7MW0i+N88rxK6+rHuZVSA4HxrrvAQK31Fd4srC4cDodOS0ur9fEL0/fwp4VruaB7BK9NdhAU0LQeqxNCiLorc2r+tHAtizL28qcrzmH6pd3r9O8ppdLNWUGq5elZVQkY3VUuiTWuyo8kDuqM1poZH/7GtLnpzJ40SMJDCGEbZU7NjIW/sShjL4+M6lnn0KipmsyOu931QyNYh2OcowuzbjiXH7dkc/e76RSXSreVEMJ6TqfmsQ9/48OMPTyU0JP7R/bweQ2eBsdjSqnVSqmvlVLfUIcrx/3J+MFd+McN/flhczb3vJsh4SGEsJTTqfmfRetYkL6HB0b24MEE34cGeN5VNUtrvcR1Ryk10kv12M5N50fj1PD4R+u4990MXr41nubNpNtKCOFbTqfm8Y/WMS9tN/df1p2HLAoN8KDFYa76l2bebmNez5Hp3bLs5ZYh0fxtbD+WbDrE9JQMTpY6rS5JCNGIOJ2aJxav44PVu7nv0u48PKonSinL6qk0OJRSOUqpG8A4Hdf125xKvVF0Vbm7dWhXnhnbj9SNh7hXwkMI4SNOp+bPH6/n/VW7ufeSOB653NrQgKpbHP/QWi/SWu9QSt2olPqjK0g4/QyrRmPS0K48PaYvqRsPMv09CQ8hhHc5nZonP17Peyt3cc8lcfzpinMsDw2oOjjKu6O01h8C3bXWi8xNZ0350VhMHhbD02P68u3vEh5CCO9xhUbKyl3cfXEcM2wSGlB1cAxWSg1w/QBH3G4P9lF9tiThIYTwpjND49Er7RMaUPVZVeMw5oVyr/Zx8/dA4H+8VZQ/mDwsBoD/9/EGpr+XwUu3xBPYzNOzm4UQomLuoTHt4ljbhQZUHRzT3E/BddeYTsetint43JuSwcsTJTyEELVn95aGS6XfcpWFRnX7GhtXt1XqxoPcmyJXmAshasd19pTdQwM8v3JcVOFUeBziXrnCXAhRQ67rNN7zg9AACY56M3lYDM+YFwnePTddpmQXQnjENY3I+6t2M/1S+4cGSHDUq0lDuzLz+n58v9mYGFHCQwhRlTKn5rFFv5VPI/LHy+0fGiDBUe8mDunK/5oTI8piUEKIyrjW05ifZkxYaPU0IjVhSXAopRKVUglKqaRK9ieZP7Pcts1y7fNVnbV18/nRPHvjuSzddpi73pY1zIUQpystc/LI/F9ZlLGXhxJ6+lVogAXBoZSKB9Bap7rfd9ufAKRqrZOBWPM+QJJSKhPI8mW9tTV+cBeeSzyPXzIPc8ec1RSeLLW6JCGEDZSWOXlo/loW/7qPP11xjmVTo9eFFS2OCUCeeTsLY3VBd7Fu27LM+wBTtdZxrsDxB4mDOvOv8eexcnsOt721mmPFEh5CNGYlZU4e+GANn67dx2NX9fL5yn31xYrgCAVy3e6Hu+/UWiebrQ2AeMwp3YEws3trhg9qrDfXD+zM/900kPSdR5jy5ioKikqsLkkIYYGTpU7uey+DL9Yd4M9X9+bui+OsLqnWbDs4bnZhZWitM6A8UFKBcLfuK/fHJyml0pRSadnZ2b4ut0rXndeRF28eyNrdeUx6YxX5JyQ8hGhMikrKuOfddL7ecJC/XtuHuy6Mrf4gG/N0BcAaqWQAO1drvRCjmyrM3BZK5TPtJmitH3X791zH53Cq+6qc2UpJBnA4HLpu/wX176r+HXi5iWL6exlMfH0Fc+8YQtuWgVaXJYTwshMny0iam8bPWw/zt7H9uHVoV6tLqjOvBIdbV1NF5gEO83Ys4BokD9Va55m3k7TWz5q3EzDGOlxdVuGuY/zN5X3bkzzZwbS56dz82grm3jmEyNbNrS5LCOElx4tLufPt1azcnsuziecy3tHF6pLqhc+7qlxdT2Yg5LnuA0vcts9SSmUqpY6Yx6QCCUqpRCDH7Ri/c+k5Ubx122B25BznpuTlHDxaZHVJQggvOFpUwpQ3V7F6xxH+b8KABhMaAEpr2/Xq1JnD4dBpaWnVP9BCq7bncvtbq4ho3ZyUu4bQuW2w1SUJIepJXuFJJr+5it/3HeWFmwcyun8Hq0uqllIqXWvtqP6RNh4cb+jO7xbGu3cN4cjxk0yYvYIdh49bXZIQoh4cPlbMTckr2LS/gNmTBvlFaNSUBIeFBka35f2koZwoKWP87OVsPVhgdUlCiDo4kF/EhNnL2ZFznDduczCydzurS/IKCQ6L9e3Yhg+ShgIwIXkF6/fmW1yREKI2ducWMm72Mg4eLeadO4ZwYY9Iq0vyGgkOG+jZrjXzpw2jRUBTbn5tBek7c6s/SAhhG9sOHWPcq8spKCol5a4hnN8trPqD/JgEh03ERLRkwd3DiGjVnElvrGLp1sNWlySE8MCGfflMmL2cUqfmg6ShnNcl1OqSvE6Cw0Y6hrZg3rShRIcFc8ec1Xyz4YDVJQkhqpC+M5ebklfQvFkT5k8bSq/2IVaX5BMSHDYT1TqID5KG0qdjCPekZLAoY4/VJQkhKvDz1mxufX0VEa2as+Ce4cRGtrK6JJ+R4LCh0OBA3r1rCEO6hfHw/LW8vWyH1SUJIdx8uW4/d85Jo2t4MPOnDaNTaAurS/IpCQ6batW8GW/eNphRfdrxl0828J/UrTTEizWF8DfzVu9i+nsZ9O/chnlJwxrltEESHDYWFNCUVybGc2N8Z/6duoWnP/sdp1PCQwirJP+UyaMfruPCHpHMvfN82gQHWF2SJbwyyaGoP82aNuG5xHNp0yKAN3/ZTl5hCc8mnktAU8l8IXxFa80/vtrE7B+zuPrcDvx7/AACmzXez6AEhx9o0kTx5DW9CW8VyHNfbyb/RAkv3RJPi8CmVpcmRINXWubkfxatY0H6Hm4dGs1T1/WjaRP/WR/cGxpvZPoZpRTTL+3O36/vzw+bD3HrGyvJKzxpdVlCNGhFJWXck5LBgvQ9PDiyB89T5TKIAAAVzklEQVSMkdAACQ6/c8uQaF66JZ51e/IZ9+py9uWdsLokIRqk/MISJr2xktSNB3nqur48NKonSklogASHX7qqfwfm3DGYA/lF3PjKMpkcUYh6tj//BONmL2Pt7nz+e/NApgyPsbokW5Hg8FPD4yKYN20YpU5N4qvLWb1D5rcSoj5sOVjAjS8vY19eEXNuH8w153a0uiTbsSQ4lFKJSqmEStYmRyk1y/yd5OkxjVGfjiEsumc44S0Dmfj6Sr5av9/qkoTwayuzckh8ZRklTs28aUMZ3j3C6pJsyefBoZSKh/LlYMvvnyFJKZWJsda4p8c0Sl3Cgll4z3D6mVOU/PiX/4OYGGjSxPidkmJ1iULYW0oKxMSgmzSh88A+TNjyM4vuGU7fjm2srsy2rGhxTADyzNtZQEIFj5mqtY5zBYWHxzRaYS0DeW/qUJ44ksHgvz8GO3eC1sbvpCQJDyEqk5KCTkqCnTtRWtPp6CEeX/xvunz5kdWV2ZoVwREKuHfIh1fwmDCzW2pGDY5p1IICmnLnl68TXFp8+o7CQnjiCWuKEsLm9OOPowoLT9umTshnpjq2HBzXWiebrY1wpZRHrQulVJJSKk0plZadne3lCu1J7d5d8Y5du3xbiBB+oPBkKeySz0xteOXK8UoGsHO11gsxupxcy2OFAjkVHOt6bA4QW90xYIQNkAzgcDga54RO0dFG99QZSjp1pnHOqCNExQ4dLeLOt9N4JSSCzkcr+EMzOtr3RfkRr7Q4zBbDmT8Lzd3zMMIA87drwNu1bFaWaxtGl1RaZceIM8ycCcHBp206EdCcJ8+/heWZZ2WtEI3Sxv1HGfvSL2RmH+PIE3896zNDcLDxWRKV8nlXldY6A8Dsgspz3QeWmPtTgQSlVCKQo7XOqOIY4W7iREhOhq5dQSno2pXCF18hbcRVTH5zJQvSKmmWC9FIfL/pEONeXU6Z1syfNoz+M+496zNDcrLxWRKVUg1xjQeHw6HT0tKsLsM28k+UMD0lg6XbDjPtolhmXNlL5tsRjYrWmjeWbufvX2ykd4cQXp/ioEObxrX4UnWUUulaa4cnj7Xl4LioX21aBPDW7YO5dWg0s3/KYtrcdI4Xl1pdlhA+UVLm5PGP1vG3zzdyeZ/2LLh7mIRGHUlwNBIBTZvwt7H9eeq6vny36SA3vrKM3bmF1R8ohB/LPX6SW19fyfurdjP90jhenhhPcKCsJlFXEhyNzJThMcy5/Xz25Z1gzEu/sDJLBs1Fw7TpwFGue3Epa3bn8e8J5/GnK3rRRLpo64UERyN0Uc9IFk8fQWhwABNfX8m7K84+hVcIf/b1hgPc+PIyTpY6mT9tGNcP7Gx1SQ2KBEcjFRvZisXTR3BBjwj+vHg9/7PoN4pLy6wuS4g6cTo1//pmM9PmptO9XWs+vf8CBnQJrf5AUSMSHI1YSFAAb0wZzPRL43h/1W5uSl7BwaNFVpclRK0cLSph6jtpvPDdNsY7OjMvaSjtQoKsLqtBkuBo5Jo2Ufzpil68PDGezQcKuPqFpTLuIfzOpgNHue6/S/lxSzZPj+nLrBvPJSigqdVlNVgSHAKA0f07sHj6CEKCmnHL6yt5/ecsGuI1PqLh+fjXvVz/0jKOnyzj/aShTB4WI0u8epkEhyjXs11rFt83gpG9ovjb5xu57701FBSVWF2WEBUqLi3jLx+v58EPfqV/pzZ8fv8FDI4Jq/5AUWcSHOI0IUEBzJ40iEev7MWX6/cz5sVf2HTgqNVlCXGaPUcKGT97BW8v38mdF3QjZeoQomQ8w2ckOMRZlFLcc0kc700dSkFxKWNf+oX5abul60rYwnebDnL1C0vJOnSMV2+N58lr+hDQVL7KfElebVGpobHhfP7ABQzs0pYZC3/j4flrOSZTlQiLnCx1MvPz37ljThqdQlvw6f0XcGW/DlaX1SjJtfeiSlGtg3j3riG8+N02/rNkC2t35/HfWwbKeszCp3bnFnLf+2tYuzuPycO68vjo3nLWlIWkxSGq1bSJ4sGEHrw3dSjHT5Zy/UvLeHPpdum6Ej7x8a97Gf2fn8nKNrqmnh7TT0LDYhIcwmNDY8P58sGLuKhnBE9/9ju3z1lNdkFx9QcKUQvHikt5eN6vPPjBr/Rs35ovHrhQuqZsQoJD1EhYy0Bem+zgmTF9WZ6Zw5X/9xPf/n7Q6rJEA7N6Ry5X/ecnFv+6lwdH9mBe0lC6hAVXf6DwCUuCQymVqJRKqGhtcqVUvFJKK6UyzZ/Z5vZZ5u+K1jMXPqSUYtKwGD69/wKiQoKY+k4aj334m6zxIersZKmTZ7/axITZywGYP20YD43qSTM5a8pWfP5/QykVD+VLxJbfdxOmtVZa6zhgHDDL3J6klMrEWJNc2EDPdq1ZPH04d18cx7y03Vz1n59luhJRa661wF/+IZPEQZ358sGLcMgFfbZkRYxPAPLM21lAgvtOV6CYHFprV1BM1VrHnbFfWKx5s6Y8dlUv5iUNA+Cm11bw9Ke/U1QiM+0Kz5SWOXnxu61c9+JSDhUUkTxpEM8mnker5nLSp11ZERyhQK7b/fCKHqSUSgDmu20KM7u3ZlTy+CSlVJpSKi07O7v+qhUeOb9bGF/94UImDe3Km79sl9aH8MjG/Ue54ZVlPP/NFq7o255vHrqYy/u2t7osUQ07dxyO0lq7WiZorZPN1ka4GSqnMfc7tNaOyMhInxYqDMGBzXh6TD/eu2sIpU4nE5JX8OTi9TLflThLcWkZ//pmM9f+dyn78k7w8sR4XrwlnrCWgVaXJjzglbZgJQPYuVrrhRjdVK6Oy1Cgsj9Ly8c+zH/PdXwOEFuP5Yp6Nrx7BF//4SL++c0W3vxlO6kbD/LUdX3lL0kBwKrtuTz+0Tq2HTrGDfGdePLqPrSVwPArXgkOrXVyFbvnAQ7zdizgGiQPdbUwlFJnBkMWkGbeDncdI+wrOLAZT17Th6vP7cDji9aRNDedy/u046kxfenQpoXV5QkL5BWe5H+/2MS8tN10btuCObcP5pJzoqwuS9SCz7uqtNYZUD6Gkee6Dyw546FZbsekAglKqUQgx+0YYXPx0W359P4LePTKXvy0NZuEf/5I8k+ZlJQ5rS5N+IjTqZm/ejeX/fNHFmbsYdrFsXz70MUSGn5MNcRpIxwOh05LS6v+gcKnducW8tSnG0jdeIgeUa14akxfhsdFWF2W8KL1e/N58uP1rNmVx+CYtjw9ph+9O4RYXZaogFIqXWvtqP6REhzCAks2HuSvn25gd+4JrurXnsdH95arghuY7IJinv96M/PTdxPeMpDHR/fm+oGdZGU+G6tJcMiJ0sLnRvZux4juEbz+cxYvfZ/Jkk2HuPOCbtx7SRytgwKsLk/UQVFJGXOW7eDF77ZRVFLGnSO6cf/IHrRpIf9fGxJpcQhLHcgvYtZXm/hozV7CWgby4Mge3DIkWhbm8TNOp+aTtft47uvN7M07wWW9onji6t7ERbayujThIemqkuDwO7/tyePvX2xkRVYu3SJa8tConlzTvwNNmkjXhp1prflxSzbPf7OZ9XuP0rdjCI+P7s2I7jJ25W8kOCQ4/JLWmu83H+LZrzaz6UABvdq35pHLzyGhd5T0jdvQyqwcnv9mM6t3HKFLWAseHtWTMed1krD3UxIcEhx+zenUfPrbPv717RZ25hTSv1MbHhjZQwLEJlZk5fCf1K0sz8ohqnVz7h/ZgwmOLgQ2k+5FfybBIcHRIJSUOfkoYy8vfr+NXbmF9O0Ywr2XdOfKfu1pKn/V+pTWmp+2Hubl77excnsuka2bM+2iWCYO6UqLQFmNryGQ4JDgaFBKypwsXrOXl3/IZPvh48SEB5N0URw3xHeSJUS9rLTMyRfrD/DqD5n8vv8o7UOCSLoolluGRMtr38BIcEhwNEhlTs03Gw7w8g+ZrNubT1jLQCYOiWbS0K5EhQRZXV6Dkl9Ywgerd/H2sh3syy8iNrIld18cx9gBnaRLqoGS4JDgaNC01izPyuHNpTtYsukgzZooruzXgUlDuzI4pq2Mg9TB+r35pKzcyeI1+zhRUsaw2HDuvKAbl/WKkkHvBk4uABQNmlKK4XERDI+LYMfh47y9fAcL0/fw6dp99GzXigmDo7l+YCeZottDBUUlfP7bfj5YvZtfd+cRFNCEMed1YvLwrvTt2Mbq8oQNSYtDNAgnTpbx6dp9pKzcydo9+QQ0VST0bsf1AztxyTlR0r1yhjKnZnlmDh+t2csX6/ZzoqSM7lGtuOX8aG4c1Fmu9G6EpKtKgqNR23yggAVpu/lozV5yjp+kTYsARvdvz9X9OzIkNqzRXpXudGoydh3hi3UH+PS3fWQXFNO6eTOuOa8j4x2dGdAlVLr5GjEJDgkOgXE21i/bDvPxr/v4esMBCk+WERocwKje7RjZux0X9Iho8OtaF5WUsSIrhyUbD/H1hgMcKigmsGkTLu0VydgBnbi0V5ScHSUAGeMQAoCApk245JwoLjkniqKSMn7aks1X6w/w1YYDLEjfQ2DTJgyJDePCHhGM6B5B7/Yhfj8ArLUmM/sYv2zL4eet2fyyLYcTJWW0CGjKJedEcmW/9lzWK0omkxR1YlmLQyk1S2v9aCX7EjGWmI11rSZY0bbKSItDVKWkzEnajiN8v/kQ3206xLZDxwAIaxnI4Ji2DI4JwxETRp8OIbYfGyktc7Ll4DHSduayescRVm3P4eDRYgA6t23BZb2iuLRXFMNiw6VlIapk+64qcw3xR7XWcRXsi8cIh4Xm41wJcNq2qlYBlOAQNXEgv4hlmYf5ZVsOq3fksiu3EIDApk3o1aE1/Tu1oXeHEM5p35qe7VrXbuA4JQWeeAJ27YLoaJg5EyZOrNE/cby4lK2HjrHlQAG/7z/Kur35bNiXT1GJsZpi+5AgHDFtGdE9ghFxEUSHyxonwnO2Dw4ApdS3WutRFWyfBXyrtU41l5eNx1hn/LRtWutnK/u3JThEXRw6WkTaziOs3ZPHuj35rNuTT0Fxafn+iFbNiQkPpmt4Szq3bUGHNkG0bxNEVOsgwloG0rZlAM2buf11n5ICSUlQWHhqW3AwJCeXh8fJUid5hSfJLTzJ4YKT7M8/wf78IvYeOcGOnOPszCnkwNGi8sNbBDSlX6cQ+ncK5dzObXDEtKVTaAsZ3Ba15u9jHKFArtv98Eq2CeEVUSFBjO7fgdH9OwDGuMG+/CK2HChg88ECtmcfZ0fOcZZuy+ZQQTEV/e3VvFkTWjZvRnBgUxb+4xHau4cGQGEhB+57mGu2R3G8uIwTJWUV1hLRKpCY8JaM6B5Bt4hgerYzWj1dwoJlvi5hGTsGR62YXVhJANHR0RZXIxoSpRSdQlvQKbQFl/aKOm1fSZmTQwXF7M87weFjxeQeLyH3eDEFxaUcLy7leHEZ7fIOVfjvtsvP5oq+7WnZvBmtmjcjrGUgYS0DCW8ZSIc2LYgKaS7jEsKWvBIc5pf4mXK11gs9ODwPCDNvhwI55u2KtpUzB8yTweiqqlHBQtRSQNMm5aFSqUejYefOszar6GhmXt/fi9UJ4R1eCY7qznqqiFIqVGudB8wDXP1ssUCqebuibULY38yZFY9xzJxpXU1C1IEl5xqap9Y6zmiZLAFwnS1lDoLnaa0zKtrm65qFqLWJE42B8K5dQSnjt9vAuBD+Rq4cF0IIUaOzqux9dZMQQgjbkeAQQghRIxIcQgghakSCQwghRI1IcAghhKiRBnlWlVIqGzj7iquqRQCHvVBOfbBrbXatC+xbm9RVc3atza51Qe1q66q1jvTkgQ0yOGpDKZXm6alovmbX2uxaF9i3Nqmr5uxam13rAu/XJl1VQgghakSCQwghRI1IcJxS4/m1fMiutdm1LrBvbVJXzdm1NrvWBV6uTcY4hBDCD3lz+e3qNPoWh7niYGX7EpVSCe6TMVa0TQiX6t4fVr6nPKgtyfyZ5bZtlmufhXWdVYMvXrOqnkMpFa+U0kqpTPNndmW1eqm2JCCxkn3xAFrrVLdaz9pWl+dv1MFh9YtfRV22fMPa9QPuYW1e/1Ks7v1h8XuqutoSgFTzL9FY8z5AklIqE8iyoq6KavDFa+bBc4RprZXWOg4YB7jeV159vVzM/0+VPccEjJYF5mMSKtlWa406OKx+8Sti1zesXT/gnjyPD78Uq3t/WPKe8rC2WLdtWeZ9gKla6zjXa2tBXRXV4IvXrMrnOOP1cGitXe8hb79envD68tuNOjiqYdXa53Z9w9r1A+7J8/jqS7G694dV76nKnruc1jrZrd87HnCtSxBmtuRmWFFXJTX44jXz6DnMP0Lmu23y9utlCxIc9mPXN6xdP+DVPo+FX4p+x2ytlS+eZr52qUC4W0vNp+xQQxVGmSuXAraptaLltytbkrtWvLJ0rF1U0nft1bXPfWiU+1/KbmdOjFJKJVjRVD6zBl8/f3Uq+lI0t9fXa1bdh9PK95SnXxwJrjN1zM+P6/OSw6mWms/qqqSGev0SrE1dbsq7RX30elVK+XD57QYdHHZd+7yaQLPsDVuXurz9Aa+n18zbX4oVvmd8+YGuQ20opZK01s+atxMwuvVcrbNwL9VWXV0V1ZBW0TE+rgul1JnvGV+8XpjPXb78ttv33BJgkNY6QynlUGcstV3Rttpq0MFRHate/GoCzbI3bB3r8uoHvB5eM69/KVb2nsGHH+ja1mZun6WUehQjhMdprVPNs9UAcrxRmwevWYU1ePs1q64ut4dmuR3j9dfL7bkWAgvP2DbI7fZZn5e6XrvhTi4AtCHzr+EsTr94J931xjCD41Gt9TS3Y1ynFce6viAtqOusGio6xte1mR/+BRjjIKd9KZ5ZrxCiehIcQgghakTOqhJCCFEjEhxC2EQFY1dC2JIEhxA1YE4T8q2Hj61wOpvK/l2MM8Iq2hdqx9ObReMlwSFEzWQB1Z4tY7Ye8qp7nBtHZWfhmGeGVRgqQlhBgkOImkkAVnvwuMR6vggztZJrWYTwOQkOIWpmMJBhTlVSVZdVnOuGUipWKbVAKZWulPrWvJ3ott99GhT3mX5nuB5ntjrcrx8QwjKN+gJAIWohHog3L8DytEURqrUeZ4ZAXgUtEYfbtSexGNPJTHNNkeL2uDCEsAFpcQhRM3mcPjV7tcyrkOM5fWbeyiQC6ebtwW6zH8PpEzkKYRkJDiE8ZH75f4sx1cMo8wyrWPP3rGrGICaYg9+ndTeZLQz3cDht+pMzTtGVFoewBQkOITwXC6SZrYAcAPN2FsbcRO5TqlR2RlWeUsr9DKkzZ+WdDSS4zaflHhzS4hC2IFOOCFEHrll2MVoDWW6rECZ5MjdXDR4Xy6mxFSEsJcEhhBeYrYqEqr7oPXmM22MTJTSEXUhXlRBeYJ4+e2a31JkS8ODMrArGQYSwlLQ4hBBC1Ii0OIQQQtSIBIcQQogakeAQQghRIxIcQgghakSCQwghRI1IcAghhKgRCQ4hhBA18v8BIBfeNlveOsgAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "SSH_lattice.plot_dispersion()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "plot_dispersion() plots the 3(since number of unit cells in 3) points in k-space \n", "(the first and last one are the same) over the dispersion relation of an infinite\n", "crystal. " ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [], "source": [ "[V,[S0,S1,S2,S3,S4,S5]]=SSH_Haml.eigenstates()" ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([-1. , -0.5, -0.5, 0.5, 0.5, 1. ])" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "V" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First, the eigen-values are the same as the ones obtained from the dispersion calculation. Second, they are symmetric about the value 0.\n", "The second is a consequence of the chiral symmetry of the Hamiltonian, as we explain now.\n", "\\begin{eqnarray}\n", "\\hat{\\bf{1}} = \\hat{P}_A + \\hat{P}_B, \\ \\ \\ \\ \\, \\hat{P}_A = \\frac{1}{2}(\\hat{\\bf{1}}+\\hat{\\Sigma}_z), \\ \\ \\ \\ \\, \\hat{P}_B = \\frac{1}{2}(\\hat{\\bf{1}}-\\hat{\\Sigma}_z) \\nonumber \\\\\n", "H |\\psi_n\\rangle = E_n | \\psi_n \\rangle \\implies H (\\hat{\\Sigma}_z | \\psi_n \\rangle) = - \\hat{\\Sigma}_z H | \\psi_n \\rangle = -E_n (\\hat{\\Sigma}_z|\\psi_n\\rangle)\n", "\\end{eqnarray}\n", "So, if $|\\psi_n\\rangle$ is an eigenstate with eigenenergy $E_n$, $\\hat{\\Sigma}_z | \\psi_n \\rangle$ is also an eigenstate with energy $-E_n$ and the eigen-spectrum is symmetric about 0.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here, S0 is the eigenvector with eigenvalue -1 and S5 is the eigenvector with eigenvalue +1. So, we can verify if S5 is the same eigenvector(withinn a phase factor) as ($\\hat{\\Sigma}_z*$S0)." ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]]\n", "Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]]\n", "Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]]\n" ] } ], "source": [ "print(S0)\n", "print(chiral_op*S0)\n", "print(S5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Clearly, S5 is the same eigenvector as ($\\hat{\\Sigma}_z*$S0)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Since, $\\hat{\\Sigma}_z | \\psi_n \\rangle$ and $| \\psi_n \\rangle$ are eigenvectors of a Hermitian oerator,H with distinct eigenvalues, they must be orthogonal.\n", "\\begin{eqnarray}\n", "E_n \\ne 0 \\implies 0 = \\langle \\phi_n| \\hat{\\Sigma}_z | \\phi_n \\rangle = \\langle \\phi_n \\hat{P}_A |\\phi_n \\rangle - \\langle \\phi_n| \\hat{P}_B |\\phi_n \\rangle\n", "\\end{eqnarray}\n", "i.e., an eigenstate with nonzero eigenvalue has equal support on both sublattices. We, now check this for S5.\n" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "True\n" ] } ], "source": [ "is_null = (np.abs(S0.dag()*S5) < 1E-10).all()\n", "print(is_null) # Are S0 and S5 orthogonal?" ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "True\n" ] } ], "source": [ "dimH = chiral_op.dims\n", "identity_H = Qobj(np.identity(6), dims=dimH)\n", "identity_H\n", "PA = 0.5*( identity_H + chiral_op)\n", "PB = 0.5*( identity_H - chiral_op)\n", "support_S5_A = S5.dag() * PA * S5\n", "support_S5_B = S5.dag() * PB * S5\n", "print(support_S5_A == support_S5_B) # Does S5 have equal support on A and B?" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We discuss the implications of $E_n = 0$ for an eigenvector later in the context of edge states later." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Unsurprisingly, diagonalizing the Hamiltonian gives the same spectrum of eigen-values\n", "as the one obtained from the plot_dispersion() function. We shall soon illustrate, translational symmetry is a very useful computational hack. Here, we see how they can produce the eigen-values and eigen-vectors of the Hamiltonian ($6\\times6$, in our example) from diagonalizing a $2\\times2$ matrix. The reduction in size by a factor of 3 comes from the fact that the lattice of 3 cells repeats itself infinitely on both ends." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Using Translational Symmetry:\n", "Any periodic lattice Hamiltonian can be diagonalized easier exploiting translational symmetry, this feature is not specific to SSH model alone." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\\begin{eqnarray}\n", "|\\psi_n(k) \\rangle = |k \\rangle \\otimes | u_{n}(k) \\rangle \\nonumber \\\\\n", "| u_{n}(k) \\rangle = a_n(k)|a\\rangle + b_n(k)|b\\rangle \\nonumber \\\\\n", "\\end{eqnarray}\n", "The vectors $| u_{nk}(r) \\rangle \\in H_{internal}$ are the eigenstates of the bulk momentum space Hamiltonian $H(k)$ defined as\n", "\\begin{eqnarray}\n", "\\langle k | H_{bulk} | k \\rangle = \\sum\\limits_{\\alpha,\\beta \\in {A,B}} \\langle k, \\alpha | H_{bulk} | k, \\beta \\rangle | \\alpha \\rangle \\langle \\beta | \\nonumber \\\\\n", "H(k)|u_{n}(k) \\rangle = E(k)| u_{n}(k) \\rangle\n", "\\end{eqnarray}\n", "In a lattice with N unitcells, since $|\\psi_n(k) \\rangle$ is required to be periodic over N cells, it needs to be invariant with a translation by N cells, and the valid/good quantum number for k are $0,1,2,...,\\frac{2\\pi}{N}$.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The dispersion, $E(k)$ can be obtained from get_dispersion() method. A list of $H(k)$ at the valid quantum numbers can be produced from the bulk_Hamiltonian_array() method." ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([Quantum object: dims = [[2], [2]], shape = (2, 2), type = oper, isherm = True\n", "Qobj data =\n", "[[ 0. +0.j -0.25-0.4330127j]\n", " [-0.25+0.4330127j 0. +0.j ]],\n", " Quantum object: dims = [[2], [2]], shape = (2, 2), type = oper, isherm = True\n", "Qobj data =\n", "[[ 0. -1.]\n", " [-1. 0.]],\n", " Quantum object: dims = [[2], [2]], shape = (2, 2), type = oper, isherm = True\n", "Qobj data =\n", "[[ 0. +0.j -0.25+0.4330127j]\n", " [-0.25-0.4330127j 0. +0.j ]]], dtype=object)" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "(knxA, qH_ks) = SSH_lattice.bulk_Hamiltonians()\n", "qH_ks" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The array of $|u_n(k) \\rangle$(in terms of its expansion in {a(k),b(k)}) at the good quantum numbers, k can be produced with the method array_of_unk()." ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([[[-0.70710678+0.j , -0.35355339+0.61237244j],\n", " [ 0.70710678+0.j , -0.35355339+0.61237244j]],\n", "\n", " [[-0.70710678+0.j , -0.70710678+0.j ],\n", " [-0.70710678+0.j , 0.70710678+0.j ]],\n", "\n", " [[-0.70710678+0.j , -0.35355339-0.61237244j],\n", " [ 0.70710678+0.j , -0.35355339-0.61237244j]]])" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "(knxA, vec_kns) = SSH_lattice.cell_periodic_parts()\n", "vec_kns" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In both cases, knxA is simply an array containing the valid values of k, in units of $2\\pi/a$, a being the length of the unit cell i.e. 1. a is always 1 in all methods()." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([[-2.0943951],\n", " [ 0. ],\n", " [ 2.0943951]])" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" } ], "source": [ "knxA" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "bloch_wave_functions() yields an ordered array for the eigenvalues and eigenvectors(which are bloch wave functions) for the Hailtonian of the lattice." ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([(-0.5, Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.57735027]\n", " [-0.28867513]\n", " [ 0.28867513]\n", " [ 0.57735027]\n", " [ 0.28867513]\n", " [-0.28867513]]),\n", " ( 0.5, Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[ 0.57735027]\n", " [-0.28867513]\n", " [-0.28867513]\n", " [ 0.57735027]\n", " [-0.28867513]\n", " [-0.28867513]]),\n", " (-1. , Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]\n", " [-0.40824829]]),\n", " ( 1. , Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]\n", " [-0.40824829]\n", " [ 0.40824829]]),\n", " (-0.5, Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[-0.57735027]\n", " [-0.28867513]\n", " [ 0.28867513]\n", " [ 0.57735027]\n", " [ 0.28867513]\n", " [-0.28867513]]),\n", " ( 0.5, Quantum object: dims = [[3, 2], [1, 1]], shape = (6, 1), type = ket\n", "Qobj data =\n", "[[ 0.57735027]\n", " [-0.28867513]\n", " [-0.28867513]\n", " [ 0.57735027]\n", " [-0.28867513]\n", " [-0.28867513]])], dtype=[('eigen_value', '" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "SSH_lattice_TrI.plot_dispersion()" ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [], "source": [ "SSH_H_t = SSH_lattice_TrI.Hamiltonian()\n", "D = SSH_H_t.eigenenergies()" ] }, { "cell_type": "code", "execution_count": 25, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAY4AAAEICAYAAABI7RO5AAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAGe5JREFUeJzt3c9yG9eVx/HfUU244CKBCXErMtDMPqbgF0igzAMMHU2VuJErgqqy0SYl2eFalaKz0mYmoqeUjbCgrcwDxFQeIKaYeQEysrYUaTgLL6wqn1n0bbIFNshuEv0P+H6qWOy/wBGE7sNz7+1uc3cBAJDVlaoDAAA0C4kDAJALiQMAkAuJAwCQC4kDAJALiQMAkAuJAwCQC4kDAJALiQMAkMu/VB1AEa5everLy8tVhwEAjfHy5cs37r6YZdupTBzLy8va2dmpOgwAaAwz+zrrtjRVAQByIXEAAHIhcQAAciFxAAByIXEAAHIhcQBAUwwG0vKyZCZduRL9Hv25ejXarkBTORwXABpnMJDu35cOD6P5K1ekH36IkkHak1rHPb318FD66KNo+vbtQkKl4gCAqgwGUYVgJq2tnSQNKUoa0vgEcZbvv5fW1ycTYwoqDgAo2mg1UYbXrwt7aRIHABShimSRdO1aYS9N4gCASak6WcTm5qRHjwp7efo4AOAyzuqnqEK7LT19WljHuETiAID8ykgWV8Lp2ez0sqUl6dmzqON89OfNm0KThkRTFQCcbzCIRil9/fX44bGT0G5Ljx8XfuK/LBIHAKRJJoukSSeNhiSLJBIHAIwaDKR+X/ruu8m9ZgMTxDiV9HGY2aqZ9cysn7JuxczczPbCz5OwfCP8PrUPAExEfEuPtbXJJI12+6QvooS+h7KUnjjMbEWS3H07OZ+w4O7m7tclfShpIyzvm9mepP3SggUw/ZL3f1pbO900ldeUJoukKiqOW5KGYXpfUi+5Mk4oQdfd40Rx192vj6wHgIuLm6RIFrlUkThako4S8+20jcysJ+nzxKKF0Lz1oMjgAMyASTRJJYfEzkCySKrzdRw33T2uTOTum6HaaIek8g4z65vZjpntHBwclBoogAaYRJNUMlm8ejVTySKpisQxlLQQpluSxl05c9z3EZLCapg9lNQZ3Tgklq67dxcXFycZL4Cmu2yT1Px8lDBmOFkkVZE4tnRy4u9IijvJW/EGZjaaGPbj7RQ1be0UHCOAaTCpJqnNTRJGQumJw913peM+jGE8L+nFyKb7iX22JfVC1XGY2AcA3kWTVOEquQDQ3TdTlt1ITO9Lujey/nkJoQFossteuDc/T3WRQZ07xwEgG5qkSsUtRwA022WqjKWl6LkVJItcSBwAmukyD02iSepSaKoC0DyDgXTnzsWSBk1Sl0biANAc8QOU1takt2+z78coqYmiqQpAM8RVRp6EQZNUIag4ANTbZaoMkkYhqDgA1FfeKmNuTnr6lGRRMCoOAPVzkSqj3SZplISKA0C9UGXUHhUHgHqgymgMKg4A1aPKaBQqDgDVocpoJCoOANWgymgsKg4A1Vhfp8poKBIHgHLFzVNZHrA0NxfdKuTNG5JGjdBUBaA8eZqn2m3p8WMSRg2ROAAUL88t0OnLqD0SB4BiUWVMHRIHgMkbDKLO7yz9GLGlpeiW56i9SjrHzWzVzHpm1h+zfiP87mfdB0BNxI9yzZM05uaiR7iiEUpPHGa2Iknuvp2cH9E3sz1J+zn2AVAH6+v5nv/NUNvGqaLiuCVpGKb3JfVStrnr7tfjRJFxHwBVGgyk5eXslQZDbRurisTRknSUmG+nbLMQmqUe5NgHQFXyNk9RZTRaLS8AdPfNUG20zSxTdWFmfTPbMbOdg4ODgiMEIOmkylhby9Y8NT9PlTEFqkgcQ0kLYbol6Z2B3SEBrIbZQ0md8/aRjpNN1927i4uLhQQOICFvlcGjXKdGFcNxtyR1w3RHUtzh3XL3oaI+jJ2wvh3W76TtA6BCeTrBGWo7VUqvONx9V5JCE9Qwnpf0IqzfltQLVcehu++esQ+AsuXtBJ+fZ6jtlKnkAkB330xZdiMx/TzLPgBKFjdP5ak0Hj2ieWrKcOU4gPPlvRJ8fp7+jClWy1FVAGqETnCMoOIAcDY6wTGCigNAOjrBMQaJA8BpNE/hDDRVAThBJzgyoOIAEKHKQEZUHAAidIIjIyoOYNbRCY6cSBzALKN5ChdAUxUwi+gExyVQcQCzhioDl0TFAcwaOsFxSVQcwKwYDKSrV+kEx6WROIBZMBhId+5Ih6cenpmO5imcgaYqYJoNBtL9+9kTBp3gyICKA5hWVBkoCBUHMK3W16W3b8/fjg5w5ETFAUybPJ3gc3N0gCM3Kg5gmsTNU1kqjXZbevyYpinkVknFYWarZtYzs/6Y9f3ws5FYthGvKytOoDHiKmNt7fykMTcnPXsmvXlD0sCFlJ44zGxFktx9OzmfWN+TtO3um5I6YV6S+ma2J2m/zHiB2svTCd5uS0+fkjBwKVU0Vd2S9GWY3pfUk7SbWN8JP5thfScsv+vuz8sKEmgMOsFRsiqaqlqSjhLz7eRKd98M1YYkrUjaCdMLoXnrQQkxAvVHJzgqUttRVaEJa9fdd6XjhLItqZ1ovkpu3zezHTPbOTg4KDtcoFw0T6FCVSSOoaSFMN2SNO6b33P3h9JxUlgNyw910nx1LCSWrrt3FxcXJx0zUA90gqMGqkgcWzo58XckxZ3krXgDM+u7+6dhuqeor2M7rG7rpPkKmB1UGaiJTInDzH5uZstm9jMz+62ZLV/0DeOmp5AQhvG8pBeJ5Rtmtmdm34R9tiX1QtVxmNgHmH55qgwp6gSnykCBso6qarn7KzP7yt0/MLP/kPTqom+a6PxOLrsRfm9Lei9lPSOqMHvyXNAn0QmOUmRtqvrWzH6uUBVI+qageABI+asMieYplCZr4jhSNDT2v8zsrqSbxYUEzKDBQFpelsykK1eihJH1rrZ0gqNkWZuqupL+VdKvJH0W5gFMQvwM8Phxru7Z9+V+U6hA1sSx5+6fmdn77v6tmeX4ZgM4NhhEV3p//XVUXeRJEklzczRLoTJZE8cNM5OkVkgaNyT9tbCogGn0m99If/zjSbK4aNKgykDFsiaOTUmfKLru4m/u/ofiQgKm0GDwbtK4CKoM1ESmxOHu30r6OJ43s2V3f1VUUMDUWV+/XNKgykCNZEocZvb75KykX0j6oJCIgGn0+nX+fZaWomsySBaomaxNVSbpSZjuSPqqmHCAKXXtWra72ErS/Ly0uUnCQG1luo7D3T9293+EnxfiAkAgn0ePooQwzpVwKC4tkTRQe1mbqv4iKW6gHSqqOBhVBWQVJ4L19ajZ6to1mqHQWFmbqjZCpQHgom7fJlFgKmRtqnonaZjZz4oJBwBQd2MrjtA8lezLMEXNVSbpfUn/VmxoAIA6OqupamzzlJm9X1A8AICaG9tUdUbSWJb004LiAQDUXNZRVXcl3VP0vG+TtCfpfwuMCwBQU1lHVcndu2b2C3d/YWa/KDIoAEB9ZX6Qk5n9VpKb2a8VdY4DAGZQ1uG4f5b0wt3/qqip6u+FRgUAqK2sfRxb7n5Lktz9s8u+qZmtKroCvePum1nWn7cPAKAcWZuqPjezn5nZzy978Z+ZrUiSu28n589af94+AIDyZE0cX7r7/0l6T9LvzOy/L/GetxRVDpK0L6mXYf15+wAASpJ1VNWume1J+kLS3fBgp4tqSTpKzLczrD9vHwBASbImjoehg7y2zKwvqS9J165dqzgaAJheeUZVTcpQ0kKYbim6qPC89eftI3ffdPeuu3cXFxcnGC4AICnzBYATtCWpG6Y7kuIO75a7D8etH7MMAFCyrJ3jE+Puu5JkZj1Jw3he0otx68/YBwBQsgtVHGb2Y3f/50XfNO06DHe/cc56rt0AgBrIegHg+4qGxCafx/HvBcYFAKiprBVHT9KTxPxqAbEAABoga+J46e7/iGfM7MuC4gEA1FzWxPGxmW0ougjPFD3IiUfHAsAMypo4NsJzOH7i7t/yPA4AmF1Zh+N2wv2p7prZjxV1kgMAZlDWimPP3T8zs/fd/Z9mVmhQAID6ypo4boRk0TIzl3RD0l8LiwoAUFtZE8empE8U3e7jb+7+h+JCAgDUWdbE8aG7fxzPhH6O30n6S3icLABgRmTtHH/PzLbMbDnMf6LogsD3iggKaLzBQFpelsykK1ei32bS1avROqDBsiaOl+GZ450w3wkXBA7P2AeYLclksbYmff11tNwTgxAPD6WPPiJ5oNGyJo4bZvZrRcNyfyrpemiu+klxoQENMC5ZnOX776X19cJDA4qSNXFsSvrG3f9HUsvdu5LuiYoDs+giyWLU69cTDwsoS6bO8fCM8T+H6b+b2TIjqzCTBgOp35e+++5yr8PjjdFgYxOHmW1JuivpuqQNSd/EqxTdVp17VWF2DAZR89JFqotRc3PSo0eXfx2gImdVHB+Hq8SHku65+z8S96p6v6wAgcpNqsqQpHZbevxYun378q8FVGRsH0d8G/XwuzdyryqG4WL6DQbR8Nm1tYsnjaUl6dmzaGSVu/TmDUkDjce9qoA0g4F054709m3+fZeWoqYoEgSmFPeqAtKsr+dLGiQLzJA8w3F/Kek/JfUYUYWpFTdPZe0En5+PmqJevSJpYGbkGY778bkbZmRmq4quAem4+2bK+n6YvO7uD8OyDXd/aGb9tH2AS8vbPEWVgRmVteKYGDNbkSR3307OJ9b3JG2H5NAJ85LUN7M9SftlxosZkOwEPy9pzM2ddHZTZWBGlZ44JN3SyRXn+5J6I+s7iWX7Ork/1l13vx4nHGAi4irj8PD8bdtt6elTkgVmXtbO8UlqSTpKzLeTK0eaoVYkbYXphVB9rLj7p6MvGpq3+pJ0jatycZ7BQLp/P1vCkKJmqVevCg0JaIoqKo5MQhPWrrvvSlFCCdVGO9F8dSys77p7d3Fxsexw0SR5qgyJK72BEYVUHInO7aQjd3+uqJlqISxrSRp39PYSHeP9xP6HOmm+AvLLM9SWK72BUwpJHOeMetqS1A3THUlxJ3nL3Ydhuh83R4XqYl/STtinHe8D5JKneWpujv4MYIzSm6ripqeQEIbxvKQXieUbZrZnZt+EfbYV3fZkVdJhYh8gGzrBgYmponM8tSJx9xvh97ZS7oUVmqmAfKgygImrJHEApchzQR99GUBmJA5Mr6yd4Ay1BXKp7XBc4MLy3G+KobZAblQcmC40TwGFI3FgOtAJDpSGxIHmo8oASkXiQPPRCQ6Uis5xNNNgIC0vS2Z0ggMlo+JA8wwGUr8vffddtu1pngImisSB5llfz5Y06AQHCkFTFZojbp7K0jTF/aaAwlBxoBnyNE/RCQ4UiooD9RZXGWtr2ZLG/Dyd4EDBSByor7jKyNI0JUWVxuYmzVNAwWiqQn1l7QSXaJ4CSkTFgfrJ0wku0TwFlIzEgXqheQqoPZqqUA+DQdQ0lafKIGEAlaDiQPWoMoBGoeJAdfJWGRKd4EANVFJxmNmqmfXMrD9m/Ub43c+6Dxomb5Uh0QkO1ETpicPMViTJ3beT8yP6ZrYnaT/HPmiSPENtJZqngBqpouK4JWkYpvcl9VK2uevu1+NEkXEfNMFFhto+exY1T5E0gFqoInG0JB0l5tsp2yyEZqkHWfcxs76Z7ZjZzsHBweSixeTQCQ5MhVqOqnL3zVBttM0sU3UR9um6e3dxcbHgCJHLRe43RZUB1FYho6rGdGAfuftzRU1OC2FZS9Jhyr7xtoeSOuftgxrL+9ClpaWoA5yEAdRWIYnD3TfPWL0lqRumO5LiDu+Wuw8V9WHshPXtsH4nbR/UGENtgalVelOVu+9KUmiCGsbzkl6E9duSema2KunQ3XfP2Ad1xFBbYKqZu1cdw8R1u13f2dk5f0MUI8+oKYnmKaAGzOylu3fP37KmneNoqMFAunqVobbAlCNxYDIGA+nOHekw47gFhtoCjcW9qnA5g4F0/372hMFdbYHGo+LAxVFlADOJigP55a0yGGYLTBUSB/KJq4y3b7NtPzfHMFtgytBUhWziEVNra9mTRrstPX1K0xQwZag4cL6LVBkkDGBqUXFgPKoMACmoOJCOKgPAGFQceBdVBoBzUHHgBFUGgAyoOECVASAXKo5ZR5UBICcqjllFlQHggqg4ZkXyiXxmUp7nsFBlAEggccyC0ed+50ka7bb0+DFJA8Axmqqm2WAQPY1vbe0kaWQ1Nxc9ZOnNG5IGgHdQcUyr0SojD6oMAGeg4pg2VBkAClZJ4jCzVTPrmVk/Zd2KmbmZ7YWfJ2H5Rvh9ap+ZFycLsyhhZH3mdxIjpgBkVHriMLMVSXL37eR8woK7m7tfl/ShpI2wvG9me5L2Swu2ziaRLJaWogrDnSoDQGZVVBy3JA3D9L6kXnJlnFCCrrvHieKuu18fWT+b4v6LiyQLKXru97Nn0VP5SBYAcqoicbQkHSXm22kbmVlP0ueJRQuheevBmO37ZrZjZjsHBweTi7ZOLtN/EeO53wAuqc6d4zfdPa5M5O6bodpoh6TyjrC+6+7dxcXFUgMt1CSapCSqDAATU8hw3DEd2Efu/lxRM9VCWNaSdDjmZY77PsLrxfsfSupMMNz6SV7lPQlLS9Fzv0kYACagkMTh7ptnrN6S1A3THUlxJ3krrjDMbDQx7EvaCdPteJ+pQrIA0BClN1W5+6503IcxjOclvRjZdD+xz7aknpmtSjpM7NNs8Y0GL9sMFUuOkqJJCkBBKrlyPK0icfcbiel9SfdG1j8vIbTiDQbS/fvS4bgWuguYn6fDG0Bp6tw5Pj1GK4tJJg1GSQEoGfeqmqQiqok09F8AqBCJI6uyksI4JAsANUHiiFWdGJKuXJF++IFkAaCWSBxS/uduF4XbmQNoABKHFF0/UVXSIFkAaBgShyS9fl3u+5EsADQYw3El6dq14t+j3eYW5gCmAolDijqgf/SjybxWMkEkf0gWAKYEiUOKTuh/+lN00k9zJXxMyVt6jPshQQCYcvRxxG7f5oQPABlQcQAAciFxAAByIXEAAHIhcQAAciFxAAByMXevOoaJM7MDSRd9nN5VSW8mGM6kEFd+dY2NuPKpa1xSfWO7SFxL7r6YZcOpTByXYWY77t49f8tyEVd+dY2NuPKpa1xSfWMrOi6aqgAAuZA4AAC5kDhO26w6gDGIK7+6xkZc+dQ1Lqm+sRUaF4ljhLvX8otAXPlVHZuZbYzMr5pZL22ZmfUrjq0/ujyeLjO2lLg23H0zGUMVn9nI57JiZm5me5IemtmT5DZl/1+mKfq7T+IIqjqAxzGzfvip9EBOietUDHX47JIHc/ip9GAO77eajE+S3H07Ee+pZRXF1pO0HU42nURy64eT434VcaXFUMVnlhLXgrubu1+X9KGk+Bgt/fNKOUecOhaLOD5JHKrmy3hOPLU4kMeo/EAeoxYHcyz83yXf85akYZjel9Qbs6yK2DqJ994P85J0192vx/+3FcSVFkPpn9loXCOfR9fd43WlfV5p54gy/zghcUQqOYDPUIsDeYzKD+Q0dTiYz9GSdJSYb49ZVjp330w0baxI2gnTC+GE9KCKuMbEUIvPTDo+eX+eWFTm55V2jijtjxMSR6Q2X0aJA/kyKj6YGy38Nbrr7rvS8fdwW1J7tG+mLHWI4Qw33T0+KZca65hzRGl/nJA4aowD+UIqO5jPMZS0EKZbkg7HLKtSz90fSsft53G7/qFOqt7SjImhTp/ZcbNPVZ/X6DmiLDzIKVKnL2PSOweypCN3f64KD+SUGOr22b1zMKvizyxhS1J8JW9HUtx0lrasdGbWd/dPw3RPUbNGXOm2K4otLYYd1eAzM7PR71JVn9fxOULjj8WJH59UHJEtnZxUKj2AY2MO5Diutk6+pGVKi6E2n92Yg7mSzyz89dmNR7LEfxGG/8uhu++mLasitvD+G2E02jch3m1JvbDtYRmxpXxmp2Ko4jMbjSthtMO87M9r9ByRdiwWcnxyr6ogfCn2JXVqMP6/J+kLRW2TC5I+dPftRCncib8wFcR2Koa6fHYhcTx093uJZZV/ZsCknXGOOHUsFnF8kjgAALnQVAUAyIXEAQDIhcQBAMiFxAEAyIXEgcYys46ZfZFj+41JX0EebiCXvPgrbZtccU5SuGr+yyreG9OLxIHGCvejuptjl60CwrgVhjiOHR9/gTgnJlxfMDx3QyAHEgcaK1y3UeXV4MeStzkZVac4gUngliNoug1JN8MFUQ/D/IqiW07HVxk/kLSrd29H0pLUD8s7ii6k+kzST8Pvr0YvGBy9kCrcJ6hjZqvhtibjXns7EWdy3YqiauDzke33x/1bRuJZDdvc0MmdULcUbjExerFX/Bm5+834GQ7u/jAl3p3wGi1FV2dXficF1AsVBxorNAENw/S2omdybEt6rugkGp/sd8Py5Anwk8TyG+HEf1fSryQ9SUkaDyTth+13zOxBOJnvJ5PGmNc+jlPRCT4+GX8QTu6j26f+W1L+/c/D+w8VJZvf6+Tq4Hsp2yf//U/GxRtiPH79tPfGbCNxYJocpSy7MWZ5R1IrVA17ZtYKJ8p7Sr+n1Qc6OYkOw/w4p147XhEnmUSFNG77tJjTfBEqj4WQyHbDa2fd/9T7K6qAbprZS0VVB/AOEgem3Uud3B006StFf63vStp092E44d7VydMDk/Z1chLthP3HOfXa8Ypwct529+3Ew6bGbp/B54oqkqNQXXUST3tL61eJXzu57p33l/RLd3/o7sfVB5BEHwcaK5yEVxKdz52wrBeWt0JfxAMzk6L+gptmtunun4blC+G1OopuFHczDJ/dSNyuOu4LiLdfCfvH799LNgOlvPZRHKe775rZSzPbV2heStleY/4tpxJKSHhH4XUX4njCa6+E6mXFzFZCYvgqrG8puptrZ/T9JV1PbDPaDAdwk0OgTKGvJK5wOpLuJRMU0AQkDqBEib/kh2LUEhqKxAEAyIXOcQBALiQOAEAuJA4AQC4kDgBALiQOAEAuJA4AQC7/DyCKl264iluWAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.plot(D,'ro')\n", "plt.xlabel('index of eigen values')\n", "plt.ylabel('eigen values')\n", "plt.show()\n", "plt.close()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Unsurprisingly, we see a gap in the dispersion relationship and also in the spectrum of\n", "eigen values of the Hamiltonian." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "However, a calculation of the winding number yields 0, since it has trivial topology." ] }, { "cell_type": "code", "execution_count": 26, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAANMAAADDCAYAAAAY0SRCAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAEgFJREFUeJzt3Xl0VGWax/Hvm53se4BASMIeYYQQEBqxg0RbW7tFZJo5h9PuRHum9ZxpW1ub1pnpTZ2e7unTfUZtenDUFhEbhZ5xxoXFQlFp9kWRLYGELWQjZCXrO39UFZY5WSrJW7m3qp7PP6lb703uU5X88t576733VVprhBBDF2J1AUIECgmTEIZImIQwRMIkhCESJiEMkTAJYYiESQhDJExCGCJhEsKQMKsLGKrU1FSdnZ1tdRkANDU1ERMTY3UZ/ZI6B2bPnj3VWuu0/tbz+zBlZ2eze/duq8sAwOFwUFhYaHUZ/ZI6B0YpVebNerKbJ4QhEiYhDJEwCWGIhEkIQyRMQhgiYRLCEAmTEIZImIQwRMIkhCESJiEMkTAJYYiESQhDJExCGCJhEsIQCZMQhkiYhDBEwiSEIRImIQyRMAlhiIRJCEMkTEIYImESwhAJkxCGSJiEMETCJIQhEiYhDJEwCWGIhEkIQywLk1JqqVKqSClV3Ev7s66vPbYLYTeWhEkplQ+gtd7sudxNsVKqBCgdztqEGCyreqZlQJ3rcSlQ1MM6K7TW492BE8LurApTIlDrsZzSwzrJrt3Ax4apJiGGxLaTnWmtVwEopW5QShV59lCu46higIyMDBwOhzVFdtPY2GibWvoidfqGVWGqA5JdjxOBGs9GV1hqtdbrXW25nu2uoK0CKCgo0HaYXQ7sM9Ndf6RO37BqN28dXwYkF3CfiEh0PVfqfg7nLqA95tkUog+WhElrvRdAKVUE1LmXgS2u9s1AkVJqKVDj0S6EbVl2zOQ+Jur23CyPx+uHtyIhhkZGQAhhiIRJCEMkTEIYImESwhAJkxCGSJiEMMS2w4nsqL2zi6MVDRw8c4mSqkbKapqpamylvqWd1vZOWltbSdjtIH5EOGlxkWSnRDMhPZbpmYlMyoglLFT+dwUyCVM/6prbePezCjZ/cYFPSmpobusEYER4KONSokmLiyQrOZrIsBAqKipITInnUks75TXNfHisitaOLgBiI8OYPyGFoqkZfGPaSOKjwq18WcIHJEw90Fqz82Qtr3xaxvuHK2jv1IxJGsGS/EyuyUlhxthEMhNHEBKivvJ9DsdFCgu/vDSrq0tTXtvMgTN17CitYdvRKt77/AIrN37GzdNGcue8bGaNSxrulyd8RMLUzfbj1fxm01H2lteRMCKcO+dlc/vMTK4aHY9Sqv8f4CEkRJGdGkN2agy3zchEa83+03Vs3HeWt/ae5S/7zzE3N5kf3DCZOTnJ/f9AYWsSJpfymmb+5X8+Z8uRSkYnRPGzxdNYmj+GERGhxrahlGJmVhIzs5J47KYpvL7rNC9sK+E7f/iUb04fyU9uyWN04ghj2xPDK+jDpLXmlU/LeOadI4SGKJ64eQp3z88mMsxciHoSExnGfdfmsPyaLFZ9WMpzjhN8dKyaJ2/N428Lxgy4FxTWC+owXWpp54d/PsCmwxconJzG00umMypheHuGqPBQHl40kcUzMnnszQM89uZBtp+o5ukl04mJDOpfj98J2t9WeU0z97y0k7KaZp68NY9752db2htkpUTz2v1zec5xgt9sOkZJVSOr75rNyIQoy2oSAxOUH3wcv9DAHS98QnVjG6/efw33XZtji92qkBDF96+fyOq7Z3Oquok7nv+Espomq8sSXgq6MJVUNfJ3q3YAsP7BeczN7eleLtZaODmddQ/Mo6mtg2V/2MHp2marSxJeCKowVdZf5s7VO1EK1hXPZWJGnNUl9WpaZgJrV8ylua2Du17cycWmNqtLEv0ImjC1dXTx4Kt7uNjcxkv3zCE3Ldbqkvo1dVQ8q++ezZmLLTy0dh8dnV1WlyT6EDRheuadI+wtr+NXS69mWmaC1eV4bXZ2Mj9fPI3tJ6r53dYTVpcj+hAUYfr4RDUvfnySu+aN45a/GWV1OQP2ndljWZKfyX98cIJ95RetLkf0IuDDdLm9k5UbDpGTGsPjN0+1upxB++dvX0V6XCRPvHVIdvdsKuDD9OLHJzlV08xPb7vK6NCg4RYfFc5Tt+ZxpKKBtTvLrS5H9CCgw3SppZ0XHCUsnJzGgolpVpczZDdNG8mc7GR+t/UELa5LQYR9BHSY1vy1jPrLHTxy42SrSzFCKcUjN06iqqGV9XtOW12O6CZgw9Te2cXLn5xiwcRUvzp71585OclcPTaRFz8+hdba6nKEh4AN00fHq7hQ38p3546zuhSjlFLcOXccJ6ub2F0mZ/bsJGDD9Jf950iKDqdwcrrVpRh38/SRREeEsnHfWatLER4CMkwdnV18cKSSoqkZRIQF3kuMjgjjuolpbPmiUnb1bCTw/tKAg2cvUX+5IyB7JbeFU9KoqL/MicpGq0sRLgEZpj2nnMcSs3MC92Yls7Od94yQ4yb7CMgwfX7uEqMSokiPC9wL63JSY4iLDOPzc5esLkW4DPhKW6XUozin0Nyktd5qvqShO1HVyIR0+48KHwqlFBMyYimplIsH7WLAPZPW+lda6yeA8UqpJYPdsFJqqWs29eLBtPflfN1lxiZHD7Y0vzEmKZrzl1q8/4Y1ayA7m69ffz1kZzuXhTEDDpNS6hml1NPA+MFuVCmVD1em27yy7G17Xzq7NLXNbaTGRg62PL+RGhtBVUOrdyuvWQPFxVBWhtIaysqcyxIoYwZzzLROa/2E1vpx4OQgt7sM54zr4JwMumiA7b2qb2lHa6gYyH9sPxUTEUZLe6d3p8dXroTmbpe/Nzc7nxdGDObuRO4/crTW+wa53USg1mO5+40Y+mx37foVA2RkZOBwOK60VTQ5L09wHD6LwzG8Z7oaGxu/UouvnT3dRpeGDxwOQvq5IczXy8vpaQ1dXs62Yax5IIb7/Rwqr8KklFoB5AObgL1KqSVa67d8WlkfXJNLrwIoKCjQhYWFV9oaWzt4/KP3uGNOLoWFU4a1LofDgWctvrar9QihpaVcv3Bh/ytnZTl37bpRWVnDWvNADPf7OVTe7uaVAo8DyvU1d4jbrcN5RhCcvVDNANt7FR0eSmiIIsT6O3f5XMPlDmK8vUbrF7+A6G4nZaKjnc8LI3rtmZRS7wMlOHujZGCX1vpN4E0D210HFLge5wLuEw2JWuu63tq9ERKiSI2NoLLeywNzP1ZZ30p6vJefpS1f7vy6ciW6vByVleUMkvt5MWR99Uw/4sveaAKwXim1Tin1Q6XUjKFsVGu9F0ApVQTUuZeBLf20e2VMUjRlQXCvubLaZsYkDeB2zsuXw6lTbNu6FU6dkiAZ1mvP5HFy4Su9kVJqJs7T4vuHsmHXcU/352b11e6tiemxvPd5BVprW9yp1Rc6OrsorWpk/nj73UQzWA3mQ9t9rt0925o+JoGLze2U1QRu7/TF+QZaO7qYPiZwLnz0dwE5Ns89CPSvJ70+b+F33K/N/VqF9QIyTBPTYxmdEMWWLyqtLsVnPjhayaSMWJkczUYCMkxKKW7Iy2DbsSoaLrdbXY5xVQ2t7CitpWhqhtWlCA8BGSaAxTMzae3o4u2D560uxbiN+87S2aVZkp9pdSnCQ8CGacbYRKaMjOPlTwLrLj6dXZqXPz1FwbgkJqTbdxaPYBSwYVJKcd+1ORypaOCDo4Fz7PT2wXOcudjC/QtyrC5FdBOwYQLnrl5WcjT/9t4xurr8v3dq7+zit5uPM2VkHDfmjbS6HNFNQIcpPDSER26cxOHz9by+y//vgPpfH5/kZHUTj35jMiHBMPjQzwR0mAC+ffVo5uQk8+y7R7hQf9nqcgatvKaZ324+zqIp6SySs3i2FPBhUkrx9JLptHZ08uj6g365u9fR2cUP3thPqFL8dPE0q8sRvQj4MAGMT4vlJ7fk8eGxKn7vh7PvPfPOEXaXXeTnt08jUz6kta2gCBPA8muyWJKfyb9vPsZ/HzhndTleW7uznP/cfpK7v5bNbTPkcyU7G8xl635JKcUvb5/OmYstPPLGfuIiw1g4xd53fP3fg+dZueEQhZPTWHmL/856GCyCpmcCiAoP5Y93FjApI44H/rSHdz+rsLqkXm3Yd4aH1u5lZlYSzy3PJzw0qH5VfinofkMJI8J57f655I2O53tr9rB6+0lbjZDQWvP7Lcf5x3UHmJOTzCv3ziE6Imh2IPxa0IUJICE6nLUr5nJjXgY/e/swD7++n8bWDqvL4lJzO8V/2sOvNx1j8YzRvHzvHGIiJUj+Imh/UyMiQnl++Sye31bCr98/yt6yizxzx3TL5r7dfPgCP95wiNqmNp66NY975mcH7FXCgSooeya3kBDFPyycwJ8fnEdkeAjfXb2T7726h7Ka4bt/94nKBu59aRf3v7Kb5JgINvz9fO69NkeC5IeCtmfyNGtcMv/38AL++GEpzzlKeP/wBW6fmcl91+YwdVS8T7Z56MwlVn1UytsHzxETEcaPvzmFu7+WE5CTswULCZNLVHgoDy2ayLLZY3lhWymv7Sxj/Z4zFIxL4raZmdx01UjS4oZ2//IL9Zd559B5Nuw/x4HTdcREhPLAdeNZsSCHlCC4N3qgkzB1kx4fxVPfyuPhRRN4Y/dp3th9hic3fsaTGz9jWmY81+SkcPXYRCZlxJKdEkNUeM83gWxp6+RkdRPHKxvYV17HjtIajlQ0ADB1VDz/9K087pg1hvio8OF8ecKHJEy9SIyOoPi68axYkMsX5xvYeuQCHx6r5tUdZaze/uV8BXFRYSRFRxAZFkJTUzOhO7dS19ROg8fZwRHhocwYm8iPbprCoqnpTMqQi/oCkYSpH0op8kbHkzc6nu9fP5H2zi6OX2ikpKqRspomqhvbqGtuo7Wji+rqFsaMSiZhRDhpcZFkp8QwPj2GCWmxhMmHrgFPwjRA4aEhV8LVnfNG80O62a3wY/LvUghDJExCGCJhEsIQCZMQhkiYhDBEwiSEIRImIQyxLExKqaVKqSLXzOk9tT/r+tpjuxB2Y0mYlFL5AFrrzZ7L3RQrpUpwTk4thO1Z1TMtwzmjOjjDUtTDOiu01uPdgRPC7qwaTpQI1Hos9zQxa7Jrguh8rfW/eja4dv2KATIyMnA4HL6qc0AaGxttU0tfpE7fsO3YPPcE0UqpG5RSRZ49lKttFUBBQYEuLCy0pshunGPzCq0uo19Sp2/4LEy9nDio1Vqvx7mL556MNRH4yuSzru91r1sD5PqqTiFM8VmY3D1LL9YBBa7HuYD7RESi1roO53HUbld7irtdCDuz5ASE1novgOuYqM69DGxxtW8GipRSS4Eaj3YhbMuyY6aeei6t9SyPx+uHtyIhhkZGQAhhiIRJCEMkTEIYImESwhAJkxCGSJiEMETCJIQhEiYhDJEwCWGIhEkIQyRMQhgiYRLCEAmTEIZImIQwRMIkhCESJiEMkTAJYYiESQhDJExCGCJhEsIQCZMQhkiYhDBEwiSEIRImIQyRMAlhiIRJCEMkTEIYImESwhCltba6hiFRSlUBZVbX4ZIKVFtdhBekzoEZp7VO628lvw+TnSildmutC/pf01pSp2/Ibp4QhkiYhDBEwmRWX1OP2onU6QNyzCSEIdIzCWGIhMkQpdRSpVSRUqq4l/ZnXV97bLeopj7bh4sd37vBkDAZoJTKhyuzxF9Z7qZYKVUClNqhJi9r9jk7vneDJWEyYxlQ53pcChT1sM4KrfV49x+NDWrypubhYMf3blAkTGYkArUeyyk9rJPs2pV5zCY1eVPzcLDjezcoEqZhorVe5frPmqKUsqoX8Ev+8t6FWV2Av+jl4LdWa70e525Ksuu5RKCmh+91r1sD5PqyVpc+a/KifbjY8b0bFAmTl7TWfX2AuA5wjyHLBdwH04la6zqcxwK7Xe0p7nYf66+mHtstYMf3blBkN88ArfVeANcuSJ17Gdjiat8MFCmllgI1Hu1W1tRb+7Cy43s3WDICQghDpGcSwhAJkxCGSJiEMETCFASUUvlKqU1W1xHoJEzBoRSbj2sLBBKm4FAESM/kYxKm4DCbL8e32foyBn8mYQoO+cAbOEcSzLK4loAlYQoOda6hObK750MyNi/AuS622+VavAF4VinlHiyaiPN6opJ+xh4KL0iYAp/nINYSIN81AhulVCLO8W4SJANkbF6Qcl/agPPyh1K7X8XqDyRMQhgiJyCEMETCJIQhEiYhDJEwCWGIhEkIQyRMQhgiYRLCkP8HHLoXEYH0IusAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" }, { "data": { "text/plain": [ "0.0" ] }, "execution_count": 26, "metadata": {}, "output_type": "execute_result" } ], "source": [ "SSH_lattice_TrI.winding_number()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The Topologically Nontrivial Insulator\n", "The other limit ($|t_{inter}| < |t_{intra}| $) is the interesting one, although the same\n", "calculations on a lattice with periodic boundary condition reveal nothing interesting." ] }, { "cell_type": "code", "execution_count": 27, "metadata": {}, "outputs": [], "source": [ "t_intra = -0.5\n", "t_inter = -0.65\n", "H_cell = Qobj( np.array( [[ 0, t_intra ],[t_intra,0]] ) )\n", "T_inter_cell = Qobj( np.array( [[ 0, 0 ],[t_inter,0]] ) ) \n", "pSSH_lattice_nTrI = Lattice1d(num_cell=100, boundary = \"periodic\", cell_num_site = 2, cell_site_dof = [1], Hamiltonian_of_cell = H_cell, inter_hop = T_inter_cell )" ] }, { "cell_type": "code", "execution_count": 28, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYkAAAEJCAYAAABhbdtlAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAGRVJREFUeJzt3c1yHOd1xvHnQLIWWMgwPrYEDMV7ixzdgAM5+wQSUwWyVHKVoGTvSE64ZiVMsk8JlZJLRWFBW/ENCPQNZMhcQEqwqC04NJwFF1aJJ4vuFpvD7um3Z/q7/7+qKc5HD/DWcKYfnPf022PuLgAAsqy1PQAAQHcREgCAXIQEACAXIQEAyEVIAAByERIAgFyEBAAgFyEBAMhFSAAAcr3a9gBWtb297Xt7e20PAwB648GDB4/dfSdk296HxN7enqbTadvDAIDeMLNHodsy3QQAyEVIAAByERIAgFyEBAAgFyEBAMhFSADLOj2V9vaktTVpezu6hF7f24ueD3QcIQFkKQoAM+nmTenRI8ldms2iS+j1R4+i55sRJOg0QgJIJMEQEgBSdHsVyfNDgoTAQEsICYxPVpWQDgZp9QCoSjKOvMqD8EDNCAmMQ1GVIHUnGPJkVR5UG6gZIYHhygoGqfthUFa62nj/faoMVIqQwLA0HQxm0b9bW9HFrPh6+nlV+/ZbqgxUipBA/9UVDEUBsLsr3b0b/Z7Hj6PLs2fF192j5+3u1h8kWT0NAgMlEBLop7qDISQAvv5aOjpa7vccHUXPDwmS9LhWQWBgCYQE+uf0VDo+Xi0YsqqEdDCsEgBVSYIkr/KQlg+PdGAcHxMUyEVIoB/Sh62+95709Gn5n1FUJXQhGPJkVR5Z4fHaa+V/9tOn0o0bVBXI1FpImNmdBY8dmtmBmR03OSZ0TN5hq999F/4zsoKhy2FQ1nx4fPrp8lNUTEMhQyshEe/8D3MeuypJ7n6Wvo2RqKLXMPRgWCRrikoKDwz6FpjTSki4+4mk85yHr0u6jK+fSzpoZFBo3yq9hjEHQ54qA4O+xWh1sSexIelJ6vZWWwNBQ5Lq4caNcr2GV17pZsO5i/ICIxR9i9HqYkhgDPKmlUKtr0uffdb9hnMXJYHx+efR61gG01Cj08WQuJS0GV/fkDSb38DMjs1sambTi4uLRgeHCiw7rZSeUjo5IRhWdXQUvY5MQ2GBzoSEmW3EV+9J2o+v70s6m9/W3U/cfeLuk52dnaaGiFUtM61Er6Feq/YtmIYavLaObjqUNJk7xPW+JLn7w3ibA0mXyW303Hz1EIJgaNYqfQuqisEy7/kZMSeTiU+n07aHgTynp9KtW+XCYX2d6aSuSMK9zAEFu7vS7dv8/3WYmT1w90nItp2ZbsIAlake6Dd00zJ9C6qKQSEkUL2yvQemlbptmWkoehWDQUigWmWqh/X16DBMgqE/yh4+S1XRe4QEqrFM9cC0Un/NT0MtQlXRa4QElrfMgjiqh+FYpqpgIV7vEBJYzjIL4qgehqlMVcFCvN4hJLCcW7fCD4ukehi+ZU718fRp9D5CpxESKCeZYgpd90D1MC5lqgopeh8x9dRphATCceQSQnAE1KAQEigWeuQSC+KQVmYhHkdAdRYhgcVCqwcWxCFL2YV4VBWdQ0ggW5l1D7u7BAOKJYFRFBRUFZ1CSOBlZXsPt2/XPyYMx+3b9Cp6hJDAc6yaRhNYrd0rhAQiHLmEJnEEVG8QEmNH9YA2UVV0HiExZlQP6AKqik4jJMYs9NQaVA9oQtmqglN6NIKQGKPQU2tQPaBpZaoKTunRCEJibMosjqN6QFtCqwqmnmpHSIxFaIOa6gFdEVpV0NCuFSExBlQP6DOqilaZh3xZTIdNJhOfTqdtD6PbQvoPyak1gC7jvVwJM3vg7pOQbakkhqxMg5pTa6APQk7pQUO7UoTEUDHFhCFi6qlxhMTQ0KDG0NHQbhQhMSRUDxgTqopG0LgeEpp6GCve+6XQuB4bGtQYu9CG9toa008lERJ9xxQTED715M70U0mERF/RoAZeVOa8TzS1gxESfUT1AORLVxVmi7elqihE47qPaNIB4fi8vITG9VDRoAbKY5X2SgiJvmCKCVgO6ylW0kpImNmhmR2Y2XHO43fifzMfHxUa1MDqWKW9tMZDwsyuSpK7n6Vvzzk2s68knTc5ts6hegCqRVVRWhuVxHVJl/H1c0kHGdt84O5vJEEyWiHfQZ003AgIIExSVRQFBd+jLamdkNiQ9CR1eytjm814OuqjhsbUHcn00toaDWqgTjS0g3Syce3uJ3EVsWVmL1UaZnZsZlMzm15cXLQwwpqkp5eKDk1miglYDVNPQdoIiUtJm/H1DUmz9INxABzGN2eS9ud/QBwiE3ef7Ozs1DrYRoVML9GgBqpTpqE90qmnNkLinp7v+PclJQ3sjfi+8+Q+RVNRw18pF7L+wYzqAahLSFUx0qmnV5v+he7+0Mwm8TTSpbs/jB+6L+mau5/Fh8hK0iz1+DAlU0yLKoiRrQYFWnF0FF0W/cGWTD0l248Ap+VoW1EFsb5O9QA0aQR/uHFajj4ImWJiegloHlNPL2h8ugkaxV8qQK8x9fQ9KokmlTnFBusfgPYVraUYwWk8CImmcIoNoH9YS0HjujGc0x7otwF9hmlcdwnfAQEMw0hP40FI1IkpJmA4Rjr1REjUqeg0G5xiA+iXEZ7Gg5CoA2sggGEb0VoK1klUjTUQwDiMZC0FlURVWAMBjNPA11IQElWgQQ2M18Ab2qyTqMKAjp8GsIKe7AtYJ9EU1kAASBvgWgpCYllMMQGYN8CpJ0KirDINatZAAONTZi1FDxrahEQZVA8AQg2kqiAkQoRWD9LzphQBASCpKoqCosNVBSFRJLR6kGhQA8gW0tCWOllVEBJ5ylQPElNMAPKFTj1JnasqCIksZasHGtQAioQ2tBMdqSoIiTSqBwB161lVMc6QSMJgbU3a3o4uZtLNm1QPAOq3TFVx82a0n0r2WWtrjYRHUEiY2b1aR9Gk9FSSuzSbRRcpul2E6gFAVcpUFcn+KdlnuTcyJRVaSfzGzN40s5+Z2V5to2lC0RcB5aF6AFCHslXFvJq/4Cg0JL509/+R9CNJJ2b2SzP7WW2jqtM335R/DtUDgLqVqSrmLbNfCxQaEg/jKSd395+7+7+7++97GRRXroRvS/UAoEnLVhVl9mslhYbEx+5+3d1/l9xhZj+V9HY9w6pR0aIWs+hfqgcAbZmvKpL9UpaaF/GGhsSP0jfM7IeS/lbSl5WPqG7pF99M2tqKLmbRfXfvRg0hqgcAbUqqCvdov5S3z6r5j9mgLx0ys3+QNFFUUXxtZv8s6UTSm+nqog2d+NIhAOiROr506IG7X5e0H9/ed/c/SPrTMgMEAPTDq4HbXTOzfUkysz9IesPMXpf0w9pGBgBoXWglcSLpj+7+n5I24jLlQ0mXtY0MANC60EriA0VBoXi9hNz93+oaFACgG0IriUt3/7/kRnz469LM7NDMDszseJnHAQDNCA2Jd83sf83snpn9RtJvl/2FZnZVktz9LH079HEAQHNCQ+ITd/+JpGN3f1fS363wO6/reS/jXNJByccBAA0JDYlNM/sPSR/ERzUFnC4114akJ6nbWyUfBwA0JDQkvnL3v5d0P92baIuZHZvZ1MymFxcXbQ8HAAYrNCSuxSfz+3HctL62wu+8lLQZX9+QNCv5uNz9xN0n7j7Z2dlZYSgAgEXKrJP4uaLzNR2sePjrPaVWbktKGtQbix4HADQvaJ2Eu/9J0q+S22b2+rLTTu7+0MwmZnag6NDah/FD9yVdW/A4AKBhQSFhZm9Keje5KelNSX+17C9195OM+64tehwA0LzQFdcHildcxw5rGAsAoGNCQ+JBfNZXSZKZ9e97JAAApYWGxK/M7I6i9Qsm6ceSflLbqAAAnRAaEnfc/X5yw8z+sqbxAAA6pDAkzGxP0jS+/kNFayS+qnVUAIBOyF0nYWYzM/tr6ftDYOXuf3L332uFE/x1wumptLcnra1J29vRZW0tuu/0tO3RAcCLWtxnLVpM9y/u/rv4O63/xsx+mYSGXjzSqV9OT6XjY+nRo+gLxmez6OIe3XfzZvQF4wQGgDYlwWAW7Zfy9lnHx7XuqxaFxPdTSu7+X5L+wt1/F9/10qkyeuPWLenp0/zHPT53YQMvPgBkSv8xKz3fL2V5+jTar9VkUUi8ZWY/TS6S/pi6/lZtI6rbN9+Eb/v0qXTjBlUFgGYk1cONG4v/mJ1XZr9W0qLG9TuKzp1kqfv+Kf73TUn/WNeganXlyvN0DpVUFZJ0dFT9mAAgqR7KhEPiypXqxxNbVEl86O7X3f3d+YtW+9Khdt2+La2vl38eVQWAOixbPSTW16P9Wk1yQyK9LqLMY513dCSdnEi7u1FDaGsrukjR7SL0KgBUZb73sEiyf0r2WWbRfuzkpNYZjtBThQ/L0ZH09dfSs2fS48fRxV26ezd60YtQVQBYRdnqYXc32j+5P99nPXsW7cdqngIfZ0jkScLj88/DpqSoKgCUVaZ6WF+P9kcNhEEeQiJLekqqSM2HnwEYmKLD8BMNTCWFICTylKkqHj1i6gnAYskUU1EF0YHqIY2QKBJaVTD1BCBP6BRTR6qHNEIiRGhVQUMbQFpog7pj1UMaIVEGVQWAUD2uHtLMF50TpAcmk4lPp9Pmf3HI3OLubvSXAYDx6fA+wsweuPskZFsqiWWFrNymoQ2MT5kGdY0rpatCSCyLqScA8wYyxZRGSKyChjYAaRAN6jyERBWoKoDxGmD1kEbjumodblYBqEEPP/M0rttEQxsYh4E1qPMQElVj6gkYvoFPMaUREnUo09Dm5IBA/xSdpK+HDeo8hESdQqoKpp6A/giZYhpA9ZC26DuuUYWjo+iy6I3Fd2gD3RfyHdQda1BXgUqiKUUNbdZSAN1UZg1EjxvUeQiJptDQBvpnRA3qPKyTaEMPj6sGRmmgn1XWSXQdaymAbhvJGogQhEQbmHoCuospphe0EhJmdmhmB2Z2nPP4nfjfzMcHgbUUQDeNaA1EiMZDwsyuSpK7n6Vvzzk2s68knTc5tlaErqVYW2P6CajTCNdAhGijkrgu6TK+fi7pIGObD9z9jSRIBi+pKhYFhTvTT0BdQqaYkgb1iAJCaickNiQ9Sd3eythmM56O+ijrB5jZsZlNzWx6cXFRyyBbEdLQZj0FUJ2Rr4EI0cnGtbufxFXElpm9VGnEj0/cfbKzs9PCCGuSnnoyW7wtVQWwGhrUQWo5LUdOw/mJu3+haKppM75vQ9Is47nJtjNJ+3WMsbOS03hIxfOjSVN7pG9eYCVFDWqpl2sgqlZLJRH/pT9/+SJ++J6e7/j3JSUN7I34vvPkPkVTUT1bKVch1lMA1WMNRCmNTze5+0NJiqeRLpPbku7Hj59JOjCzQ0mz1OPjw3oKoFpMMZXGaTn6YqRnoAQqVVRBrK+PIhw4LccQ8d0UwPJYA7E0vk+iT/huCqA8qvCVUEn0Ed9NASyWVA5ra9J777EGYgWERB/R0AbypZvT7tJ33+VvyxRTIRrXfTfQ890DSwv5TEij/lzQuB4T1lIAkdD1DxJTTCUQEn3H1BMQtv7hlVei090wxVQKITEEZb6bgoY2hqTMCfo++0x69myUZ3JdBSExJFQVGBNWTzeCxvVQ0dDG0PEeXxqNa9DQxnBxgr5GERJDxdQThogppsYREkNGQxtDUaZB/fnnNKcrREiMAVUF+ozqoVWExFgkVUVRUFBVoCtCqwfpeYOagKgcITE2IQ1tiaoC7QqtHiQa1DUjJMYmdOpJoqpA88pUDxJTTA0gJMYotKGdoKpAE8pWDzSoG0FIjFnZquLWrfrHhPG6dYvqoYMIibErU1Ww+A51KLM4juqhcYQEIhwmizZweGvnERJ4jsV3aAqL43qDkMDLqCpQJ6qHXiEkkI3Fd6gai+N6iZDAYiy+QxVYHNdbhAQWY/EdVsHiuN4jJFCMxXdYBovjBoGQQDiqCoSgehgUQgLlLFNV3LwpmREYQ5YEg1n0/031MBiEBJZTpqpIvkedaahhmp9WSv6/F6F66A1CAssrW1VITEMNSdlpJYnqoYcICayuTFWRoKrotzJN6QTVQy8REqgGVcU4UD2MTmshYWZ3Fjx2aGYHZnbc5JhQgfmqwqz4OTS3u22ZpnTy/0710HuthES88z/MeeyqJLn7Wfo2eiSpKtylu3dpbvfZsk3pu3ejbakeeq+VkHD3E0nnOQ9fl3QZXz+XdNDIoFAPpqH6iWklxLrYk9iQ9CR1e6utgaBCyza3mYZqzjLTSgmmlQariyFRyMyOzWxqZtOLi4u2h4NQy1QVTEM1Y5lpJYnqYQRqCYl4Jz5/yexBZLiUtBlf35A0m9/A3U/cfeLuk52dnaqGjaYs09yWmIaqwzLTSjSlR6WWkIh34vOXLxY9x8w24qv3JO3H1/clndUxRrRsmeZ24tEj6f33pe1taW2N0Chr1WklmtKj0tbRTYeSJnOHuN6XJHd/GG9zIOkyuY0BW2Ya6ttvpdks2lnRuyiWFwxMK6GAeeibpKMmk4lPp9O2h4GqnJ5Kt25FOzGz8J1YWvK83d3oy2vGulNb9bXkdRwsM3vg7pOQbXvZuMaArTINlUg3u8dWYaxaMSSYVkKMkEB3LTMNNW/ogZGEwtpa1KP5xS+WDwaJaSW8hJBA96WPhjKTtrak114r/3OGEhhZ1YJ71KP585/L/zyOVsIChAT6Iakqnj2THj+WPv20/CG0aVmBsb3drSOm5quE7e3Vp5ES6WBgWgkLEBLop7zexSqBMZtlHzFVd3gUhUFSJcxmL463LIIBS+DoJgxLFUdHFUl+7lZ8xpgnT6TNzeWuz2b1jTM9Vo5QQkqZo5sICQxXE4HRNT/4gfT661EQXblCMCATh8ACUrVTUl2Wnkb69a+jns2zZ0wloRKEBMZhaIFBfwENISQwPlmBkRxam/QZuhIeyTiSsZkRDGgUIYFxmz+09vHjdsKjKAySsTGNhIYREkCWMuGx6nXCAB32atsDAHrl6IgdN0aFSgIAkIuQAADkIiQAALkICQBALkICAJCr9+duMrMLSSW+yf0F25IeVzicqjCuchhXOYyrnCGOa9fdd0I27H1IrMLMpqEnuWoS4yqHcZXDuMoZ+7iYbgIA5CIkAAC5xh4SJ20PIAfjKodxlcO4yhn1uEbdkwCALjOzO+7+cc5jh5IuJe27+0nefasaXSVhZncWPHZoZgdmdrzoPkAqfm+09X4KGNdxfLmTuu9O8liL43ppDG2/XmZ21czczL6KL5/kjbWGcR1LOsx57KokuftZapwv3VfFOEYVEl150TN+d1ffpL37UCfjqXsnWPTeaOv9FDCuA0ln8V+Z+/FtSTo2s68knVc9ppBxZY2hC6+XpE13N3d/Q9I7kpL3VK2vVzymkwU//7qiikHxNgc5961sVCHRlRc9ratv0r5+qBvcCRa9N1p5PwX8jv3UfefxbUn6wN3fSF7XFsaVNYbWX6+512Pi7sn7p+7Xq8iGpCep21s5961sVCFRoLEXfU5X36S9/FCruZ1g0XujrffTwt/h7iepueqrkqbx9c24OvuohjEVjitnDK2/Xon4j43fpO6q+/XqDEKifV19k/byQ93iTrBX4grsobs/lL5/3c4kbaWqr0Z1YQwLvO3uyR8nXRjrpaTN+PqGpFnOfSsb1JcO5cw3P3H3LwKenvcCV/6iL+nt9F/BqaMZ3jazgzbK3vkxNP37F8naCcb3V/V6FX0g23o/he4oDpKjZuLPTfI5mel59dXYuHLGUMtOr8y4Ur6f1mzo9cpkZhtxWN2TlKy23peUvJ+z7lvJoEJimUO+mnjRC8KrtTfpKuOq80Nd0etV904w8/3Sxoe45LhkZsfu/q/x9QNF03JJxbXV0riyxjDNek7D45KZzb9fmni9ksNZJ/H/V7Jvuy/pmrs/NLNJ/P93mfwxlHXfqgYVEkXaetELwqu1N+mK46rtQ13B61X7TjDv/aIWPsRlxhXff8fMPlYUtu+4+1l8xJgkzdoYV94Y2n69Upuep55T++sV/54vJH0xd9+11PWXPidVrY1IYzFdB8R/6Z7rxUUxD5I3RBwSH7v7h6nnJIfy7ic7xBbG9dIYsp7T5LjiD/tvFfUtXtgJzo8VQDFCAgCQi6ObAAC5CAmgBRl9JqCTCAkgR3wqjS8Dt8083Uvez1V0ZFbWYxtdO5wY40ZIAPnOJRUeuRJXBZdF26VM8o6IiY/QygwQoA2EBJDvQNJ/B2x3WPFixrOctSJA4wgJIN9bkh7Gp/NYNO30RnLFzPbN7Ldm9sDMvoyvH6YeT58qJH1G24+S7eJqIn2MPtCaUS2mA0q6KulqvKgptFLYcPd34h3+ZUaFMUmt7dhXdLqVD5PTiKS22xTQAVQSQL5LvXi68ULxCt6revEMtHkOJT2Ir7+VOsOv9OJJDIHWEBJAhnhH/6Wi0yK8HR/ptB//e6egZ3A9bky/MGUUVw7pIHjhFCFzh8VSSaATCAkg276kafzX/UyS4uvnis7Xkz7lSN6RTZdmlj5Saf7ss59IOkidXyodElQS6AROywEESs4mq+iv/PPUt+Mdh5ynqsR2+3reCwFaRUgAK4qrhYNFO/WQbVLbHhIQ6Aqmm4AVxYeszk8tzTtQwBFSGX0LoFVUEgCAXFQSAIBchAQAIBchAQDIRUgAAHIREgCAXIQEACAXIQEAyPX/PtouIWVNEjEAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "pSSH_lattice_nTrI.plot_dispersion()" ] }, { "cell_type": "code", "execution_count": 29, "metadata": {}, "outputs": [], "source": [ "pSSH_H_nt = pSSH_lattice_nTrI.Hamiltonian()\n", "nD = pSSH_H_nt.eigenenergies()" ] }, { "cell_type": "code", "execution_count": 30, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYkAAAD7CAYAAACfQGjDAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4wLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvqOYd8AAAFFZJREFUeJzt3U9z3eZ1x/Hf0TRccJHQvORWZOh2H4nsG7Bv0jdAV5kRN/JU10ttOrYarjUdNt1ok46vO8pGWChR+gJiKm+gtNI3IFrWVqJMZ6GFNaPTBQAJurwgH9x78f/7mdEIuADFMxCBw3OeB4C5uwAAmOZS3QEAAJqLJAEAyEWSAADkIkkAAHKRJAAAuUgSAIBcJAkAQC6SBAAgF0kCAJDr7+oOYF5ra2u+ublZdxgA0BrffPPNC3dfD9m39Ulic3NTR0dHdYcBAK1hZt+F7ku7CQCQiyQBAMhFkgAA5CJJAABykSQAALlIEgDQRFEkbW5KZtKlS/Hfk3/W1uL9StT6KbAA0DpRJN26JZ2cxOuXLklv3sQX/mlvC817g+jJifTpp/Hy9eulhEolAQBViKL4N38zaW/vXYKQ4gQh5SeD8/z4o7S/v5gYp6CSAIBFmqwSqvDsWWn/NEkCAOZVR2LIuny5tH+aJAEAs6g7MaSWlqQ7d0r75xmTAIBQ540r1GEwkO7dK23QWiJJAMD5qkgMl5JLsdnZzzY2pPv340HtyT8vXpSaICTaTQDwviiKZwt9913+lNRFGAyku3dLv8jPiyQBANnEkLXoBNGSxJBFkgDQb1EkjUbSq1eL+zdbmAzy1DYmYWYH52zbNbOhmY2qjAlAj6SPvdjbW0yCGAzejR1UMFZQlVqSRHLx383ZdlWS3P0wuw4Ac8s+D2lv72x7qaiOJoasWpKEu48lHedsvibpNFk+ljSsJCgA3Za2lUgMhTRxCuyKpJeZ9UFdgQDogEW0lbLTUHuQGLIYuAbQPXmzlYrY2IjvZO5RQpimiUniVNJqsrwi6cydK8mYxkiSLpf4zBIALTTvbKXlZWk87n1ySDWm3WRmK8niA0lbyfKWpMPJfd197O477r6zvr5eVYgAmmxRbSUSxHvqmt20K2lnYorrI0ly98fJPkNJp+k6AOSaZ1A6O97w9CkJYkIt7SZ3fyjp4cRn25nlceVBAWifecYeaCsFaUy7CQAKmbd6IEEEIUkAaJdZxx5oK82kibObAOB987zgh7bSXKgkADRbFEk3bsyWIGgrzY1KAkAzUT00ApUEgOahemgMKgkAzUH10DhUEgCaYZbqIfseaBJEKagkADTD/r70+nXYvh1681vTUUkAqFcUSWtrYTfFLS3F9zr07HHddaKSAFCftMUUUkFQPdSCJAGgekUGqJeWpHv3SA41IUkAqBbVQ6uQJABUo+j01o2N+BlLqBVJAkD5ilQPUtxiunOn3JgQhNlNAMqTzlza2ys2vZUxiMagkgBQjlmqB5JD41BJAChH0ZvjSBCNRJIAsFjcHNcptJsALA7TWzuHSgLAfNLXiZqFDVBTPbQKlQSA2UWRNBqFv2ua6qF1SBIAZre/H54guDmulWg3ASgubTGFDE5L3BzXYlQSAIqhxdQrVBIAwqTVw95eWIJYXmaAugNIEgAullYPoe0lXifaGbSbAFyMAereopIAkK/oAPXyMgPUHUOSADAdLSaIdhOASVEUt5eKVA8kh86ikgDwDtUDJlBJAHiHAWpMoJIAwAA1cpEkgL6jxYRz1JIkzGzXzIZmNsrZfpD8PXU7gAWY9Q7qp09JED1SeZIws6uS5O6H2fUJIzN7Ium4ytiA3qB6QKA6Bq6vSfo6WT6WNJT0eGKfm+7+sNKogD5hgBqB6mg3rUh6mVkfTNlnNWlHfV5RTEA/MECNgho5cO3u46QdNTCz4eR2MxuZ2ZGZHT1//ryGCIEWosWEGdSRJE4lrSbLK5JOshuTBLCbrJ5I2pr8B5IksuPuO+vr66UGC3RGaIuJAWpk1JEkHujdhX9LUjqAvZJ8dpx+prgVdVRpdEDXFGkxUT1gQuUD1+7+2Mx2kjbSqbung9aPJG27+2EyRVaSTjLbARRV5C1yDFBjiloey+Hu4ymfbWeWmdkEzCOKpFu3pJOTi/eVGKBGrkYOXAOYQxRJN26EJwhaTDgHD/gDumZ/X3r9OmxfWky4AJUE0BVRJK2tcQ8EFookAXQBLSaUhHYT0AUhLaalJenePRIDCqGSANostMU0GJAgMBMqCaCt0hbTRRUEg9OYA5UE0FahLSYGpzEHkgTQNrSYUCHaTUCb0GJCxagkgDZIq4e9PVpMqBSVBNB0odWDFLeY7t6lxYSFCUoSZvaR4kd4ryh+3ehDd39aYlwAUqGP2aDFhBKEtptWkqTwlbv/p6TtC/YHMK8ij9mgxYSShLabfkiqiUfJ+vclxQNAosWExgitJF5Kuirpd2Z2U9IvywsJ6LGiA9T370svXpAgUJrQSmJH0t9L+mdJXyXrABaJ6gENFJoknrj7V2Z2xd1/MDMvNSqgT4q+RY4BalQoNElsJ++cXkkSxLakv5QWFdAXRaoHiQFqVC50TGIs6VeSfi1p6O6/LS8koOOiSNrclMzCxh5SPGYDNQiqJNz9B0m303Uz2+Q+CWAGUSSNRtKrV+Ffw3sgUKPQm+n+Pbsq6WNJ/1hKREAXRVF8U1zoq0VTDFCjZqFjEibpy2R5S9L/lhMO0BFFB6MnUT2gIULbTbczq98mN9YBmKboYPQkqgc0SGi76c+S0mmvp4orCWY3AdOEPmtp0vKyNB6THNAooe2mA3d/dPFuAPTsWfGv2diIp7aSINAwQVNgJxOEmf2inHCADrh8OXzf5eX40RpPn5Ig0Ei5lUTSYso+yM8Ut5xM0hVJ/1BuaEBL3bkT/vY4qgc03HntptwWk5ldKSkeoP3Si352dtOlS9KbNyQGtE5ukjgnQWxK+rmkv5YTEtAB16+TCNAJobObbkr6TNKJ4nbTE0n/U2JcAIAGCH7HtbvvmNnH7v7IzD4uMygAQDMEv3TIzP5VkpvZvygeuAYAdFzoFNg/SXrk7n9R3G5iPAIAeiB0TOKBu1+TJHf/at5vama7iu/c3nL3cdHtAIBqhLab/mBmvzCzj+a9kc7MrkqSux9m10O3AwCqE5okvnb3/5P0gaTfmNl/zfE9rymuEiTpWNKw4HYAQEVCZzc9NrMnkv4o6WbyEqJZrUh6mVkfFNwOAKhIaJL4Ihm8bgQzG0kaSdLlIs/JAQAUUmR206KcSlpNllcU36BXZLvcfezuO+6+s76+vsDQAABZoWMSi/RA8dvtlPydDlCvnLcdAFC9ypOEuz+WJDMbSjpN1yU9umA7AKBiwY/lyDKzn7r732b9ptPufXD37fO2AwCqF3oz3RXFU1Oz75P4pxLjAgA0QGglMZT0ZWZ9t4RYAAANE5okvnH3b9MVM/u6pHgAAA0SmiRum9mB4pvcTPFLh3h9KQB0XGiSOEjeI/Ezd/+B90kAQD+EToHdSp7XdNPMfqp4ABuAJEWRtLkpmcXvsjaL/6ytxduAFgutJJ64+1dmdsXd/2ZmpQYFNF4USfv70nffvf+5Z35/OjmRPv00XuZ912ip0Epi28w+kvTz5FHh2xd9AdBZUSSNRmcTxDQ//hgnE6ClQpPEWNKvJP1a0tDdf1teSEBDpW2lvT3p1avwr3v2rLSQgLKFtps+cffb6UoyLvEbSX9OXmkKdFtaPRRJDimeVIwWC60kPjCzB2a2maz/m+Kb6z4oIyigMWatHlJLS9KdOwsPC6hKaJL4JnnH9dunsyY3152e8zVAuxUZe5hmMJDu3WPQGq0W2m7aNrMtSTKzbyV9mLScflZaZEBd8mYuXWRjI64aSArokCID19+7+39LWnH3HUmfiUoCXTNL9bC8LN2/Lz19SoJA5wRVEsk7rf+ULP/VzDaZ4YRO2t8vNvZA9YCOy00SZvZA0k1JH0o6kPR9uknxo8J5dhO6I4qkW7fiG+BCLC9L4zHJAZ13XiVxO7m7+lTSZ+7+bebZTVeqChAoXRRJN25Ir1+H7U/1gB7JTRLpo8GT5HDTzK5KemJmYzH1FV1A9QBciGc3oZ+oHoAgRabAStKKmbniZzdxpzXaa38/LEFsbMSzloCe4tlN6Jcoih/hHTLFlbulgUJTYG9fuCPQZEVaTIOBdPcu7SX0Xmi7CWivIgPUS0s8SgPIIEmg26gegLmQJNBtDFADcwkduAbahQFqYCGoJNA9tJiAhSFJoHtCWkwMUANBaDehO0JbTLwMCAhGJYFuCG0xMUANFEIlgXZLq4e9vbAWEwPUQCFUEmgvBqiB0pEk0F7cAwGUjnYT2od7IIDKUEmgXWgxAZWqpZIws10zG5rZKGf7QfL31O3osdB7IO7fl168IEEAc6o8SSSvQZW7H2bXJ4zM7Imk4ypjQ0NFkbS5KZlxDwRQsTraTdckfZ0sH0saSno8sc9Nd39YaVRopiiSRiPp1auL92WAGli4OtpNK5JeZtYHU/ZZTdpRn0/7B8xsZGZHZnb0/PnzUoJEzdLqYW8vLEEwQA2UopGzm9x9nLSjBmY2zNm+4+476+vrNUSIUqXVQ8jsJYkWE1CiUtpNOQPOL5MW0qmk1eSzFUnvvS4s+dp03xNJW2XEiAbb3w+rHiRaTEDJSkkS7j4+Z/MDSTvJ8pakdAB7xd1PFY9THCXbB+l29EAUxQkitIJYXqbFBJSs8naTuz+WpKSNdJquS3qUbD+UNDSzXUknme3osqItpo0NaTymxQSUrJab6aZVGu6+nVlmZlPfhLaYlpdJDkCFGjlwjR5JZzGFVBBUD0DleCwH6sM9EEDjUUmgekXvgWCAGqgNSQLVYoAaaBXaTagW90AArUIlgWoUGaCWaDEBDUGSQPloMQGtRbsJ5ZnlDmqSA9AoVBIoB9UD0AlUEigHA9RAJ1BJYLEYoAY6hSSBxaHFBHQO7SbMjwFqoLOoJDAfqgeg06gkMJsokm7dkk5OLt43xQA10DokCRQXRdKNG9Lr1+FfwwA10Eq0mxAuiqS1tfjprUUSBC0moLWoJBBm1uqB5AC0GpUEwuzvUz0APUSSwPnSFlPI7KWlJen+fck9HqAmQQCtR7sJ+Yq0mAYD6e5dEgPQMSQJnFVkeuvSknTvHskB6CiSBN5H9QAggySBWNGb47gxDugFkgSKT29dWuLGOKAnmN3UZ7PcHDcYMAYB9AiVRF/NUj2QHIDeoZLoG6oHAAVQSfTBLE9slageAJAkOm+WZy5JTG8FIIl2Uzel75k2K/7E1vTRGi9ekCAAUEl0TvqmuFevin8t1QOACVQSXZFWD3t7xRME1QOAHLUlCTM7OGfbrpkNzWxUZUytM9lWCn3PdBYzlwCco5YkkVz8d3O2XZUkdz/MrmNC2laaNTGkj/SmegBwjlqShLuPJR3nbL4m6TRZPpY0rCSotpinrbS8TFsJQCFNHJNYkfQysz6oK5DGWERbiTfFAZhBK2c3Je2qkSRdvny55mhKNs9sJYn3TAOYSylJImfA+aW7Pwz48lNJq8nyiqQztwkn7aqxJO3s7PiscTZaFMXvlZ6lakhtbMRPayVBAJhRKUkiuYgXYmYr7n4q6YGkneTjLUmHi4yt0UgMABqmrtlNu5J2JiqOR5Lk7o+TfYaSTtP1zlrEeIP0blD66VMSBICFqWVMImk7PZz4bDuzXLgSaZVFVAxZVA8AStLKgetWmvVJrHlIDAAqQJIo06ITg8RsJQCVauJ9Eu2WvtQnHWNYZILgXgcAFaOSmFUZVcI0tJUA1IgkMU1VCSAPiQFAQ/QzSdSdBLIuXZLevCExAGik/iWJWV/nuWi84AdAC/QvSezv15cgSAwAWqZ/SeLZs2q/H4kBQIv1bwpsFU+N5aU+ADqif0nizh3pJz9ZzL+VTQbZPyQGAB3RvyRx/br0+9/HF/hpLiWHZGNjegIgGQDokf6NSUjxhZ2LOwBcqH+VBAAgGEkCAJCLJAEAyEWSAADkIkkAAHKZu9cdw1zM7LmkWd8DuibpxQLDWRTiKq6psRFXMU2NS2pubLPEteHu6yE7tj5JzMPMjtx9p+44JhFXcU2NjbiKaWpcUnNjKzsu2k0AgFwkCQBArr4niXHdAeQgruKaGhtxFdPUuKTmxlZqXL1OEu7eyP904iqu7tjM7GBifdfMhtM+M7NRzbGNJj9Pl6uMbUpcB+4+zsZQxzGbOC5XzczN7ImkL8zsy+w+Vf9fTlP2z34vk0RdJ2seMxslf2o9aafEdSaGJhy77Imb/Kn1xE2+3242Pkly98NMvGc+qym2oaTD5MKylUlko+RCeFxHXNNiqOOYTYlr1d3N3T+U9Imk9Byt/HhNuUacORfLOD97lyTqOlnPiacRJ22O2k/aHI04cVPJ/132e16TdJosH0sa5nxWR2xbme99nKxL0k13/zD9v60hrmkxVH7MJuOaOB477p5uq+x4TbtGVPmLSO+ShGo6Wc/RiJM2R+0n7TRNOHEvsCLpZWZ9kPNZ5dx9nGlPXJV0lCyvJhefz+uIKyeGRhwz6e2F+g+Zj6o8XtOuEZX9ItLHJNGYHzyJk3YeNZ+4rZb8lvnY3R9Lb38ODyUNJsdSqtKEGM7xS3dPL8CVxppzjajsF5E+JolG4qSdSW0n7gVOJa0myyuSTnI+q9PQ3b+Q3va70z78id5Vs5XJiaFJx+xt66au4zV5jahKH1861KQfvKz3TlpJL939oWo8aafE0LRj996Jq5qPWcYDSekdsFuS0vbXtM8qZ2Yjd/+PZHmouDWRVrCDmmKbFsORGnDMzGzyZ6mu4/X2GqH8c3Hh52cfK4kHencBqfVkTeWctGlcA737gazStBgac+xyTtxajlnyW+VOOqMk/U0v+b88dffH0z6rI7bk+x8ks8K+T+I9lDRM9j2pIrYpx+xMDHUcs8m4MiYHs6s+XpPXiGnnYinnZy+f3ZT8ABxL2mrA/PqhpD8q7iWuSvrE3Q8z5exW+sNRQ2xnYmjKsUuSxBfu/lnms9qPGbBo51wjzpyLZZyfvUwSAIAwfWw3AQACkSQAALlIEgCAXCQJAEAukgQAIBdJAgCQiyQBAMhFkgAA5Pp/TZxLkv+MKUQAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.plot(nD,'ro')\n", "plt.ylabel('eigen values')\n", "plt.show()\n", "plt.close()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Topologically Nontrivial Insulator: Hardwall Boundary Condition\n", "To reveal the topological nontriviality, we have to put the insulator next to a \n", "boundary, i.e. use a hardwall boundary condition and the spectrum of eigenvalues of\n", "the Hamiltonian shows the conspicuous edge states which plotted as a eigenvector\n", "are dense towards the edge considerably. However, the dispersion looks completely the\n", "same." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "t_intra = -0.5\n", "t_inter = -0.65\n", "cell_H = Qobj( np.array( [[ 0, t_intra ],[t_intra,0]] ) )\n", "inter_cell_T = Qobj( np.array( [[ 0, 0 ],[t_inter,0]] ) ) \n", "apSSH_lattice_nTrI = Lattice1d(num_cell=100, boundary = \"aperiodic\", cell_num_site = 2, cell_site_dof = [1], cell_Hamiltonian = cell_H, inter_hop = inter_cell_T )" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "apSSH_lattice_nTrI.winding_number()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "apSSH_H_nt = apSSH_lattice_nTrI.Hamiltonian()\n", "[n_D,Vx] = apSSH_H_nt.eigenstates()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "plt.plot(n_D,'ro')\n", "plt.ylabel('eigen values')\n", "plt.show()\n", "plt.close()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We observe the much celebrated mid-gap edge states, who appear at mid-gaps. The nontriviality of the topological insulator is revealed with their placement next to a boundary. The edge states are the consequences of the bulk-bundary correspondence. Plotting the absolute values of the eigen-functions reveal that they are localized at the edges. We plot the absolute values of the eigenfunctions corresponding to the two eigenvalues found at mid-gaps." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "xA = [i for i in range(200)]\n", "Es0 = np.abs(Vx[99])\n", "Es1 = np.abs(Vx[100])\n", "fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(10, 4), sharey=True)\n", "ax1.plot(xA, Es0, label=\"State0\")\n", "ax2.plot(xA, Es1, label=\"State1\")\n", "ax1.legend()\n", "ax2.legend()\n", "plt.show()\n", "fig.suptitle('Mid-gap Edge states')\n", "plt.close()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "It is quite obvious, the edge states are localized towards the edges.\n", "The eigenstates are not unique, they are only unique within a phase factor. The routine used in QuTiP only yields eigenstates that are one of the possibilities within the phase factor." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us check out the two edge state eigenvalues as well as the three preceding and three succeeding ones." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "n_D[96:104]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "It is obvious that, n_D[99] and n_D[100] are values 0 within numerical precision. There is an interesting property of the eigenstates of a chiral Hamiltonian with eigenvalue 0." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\\begin{eqnarray}\n", "H|\\psi_n \\rangle = 0 \\implies H(\\hat{P}_{A/B}| \\psi_n\\rangle) = \\frac{1}{2}H(\\hat{\\bf{1}} \\pm \\Sigma_z )| \\psi_n \\rangle = \\frac{1}{2}(\\hat{\\bf{1}} \\mp \\Sigma_z )H| \\psi_n \\rangle = 0\n", "\\end{eqnarray}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Therefore, $\\hat{P}_{A/B} |\\psi_n\\rangle$ is an eigenstate with eigenvalue 0 as well. As a result, the edge-states can be chosen to have support in only one sublattice." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "chiral_op_nTrI = apSSH_lattice_nTrI.distribute_operator(sigmaz())\n", "dimH_nTrI = chiral_op_nTrI.dims\n", "identity_H_nTrI = Qobj(np.identity(200), dims=dimH_nTrI)\n", "identity_H_nTrI\n", "PA_200 = 0.5*( identity_H_nTrI + chiral_op_nTrI)\n", "PB_200 = 0.5*( identity_H_nTrI - chiral_op_nTrI)\n", "\n", "xA = [i for i in range(200)]\n", "Es0 = np.abs(PA_200*Vx[99])\n", "Es1 = np.abs(PB_200*Vx[99])\n", "fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(10, 4), sharey=True)\n", "ax1.plot(xA, Es0, label=\"State0\")\n", "ax2.plot(xA, Es1, label=\"State1\")\n", "ax1.legend()\n", "ax2.legend()\n", "plt.show()\n", "fig.suptitle('Mid-gap Edge states')\n", "plt.close()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "$\\hat{P}_{A} |\\psi_n\\rangle$ and $\\hat{P}_{B} |\\psi_n\\rangle$ are orthogonal eigenstates with eigenvalue 0, and we can choose the two edge states to be localized at two edges of the lattice. Just to be sure, we also check if $H\\hat{P}_{A/B} |\\psi_n\\rangle $ is indeed 0." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "is_null = (np.abs(apSSH_H_nt*PA_200*Vx[99]) < 1E-10).all()\n", "print(is_null)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "qutip.about()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "qutip.cite()" ] } ], "metadata": { "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.5.4" } }, "nbformat": 4, "nbformat_minor": 2 }