{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Nanostructures in equilibrium with KWANT\n", "\n", "### © Marko Petrović, University of Delaware\n", "[PHYS824: Nanophysics & Nanotechnology](https://wiki.physics.udel.edu/phys824) \n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## What is covered in this notebook\n", "- How to create a square lattice nanoribbon with KWANT.\n", "- How to create a system with a predefined shape.\n", "- How to compute energy-momentum dispersion of a semi-infinite lead.\n", "- How to build a honeycomb lattice (e.g., graphene).\n", "- How to compute LDOS of graphene nanoribbon. \n", "- How to implement a magnetic field. \n", "- How to introduce vacancies." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 1. Creating a closed system on a square lattice " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To create a tight-binding system, KWANT uses Builder object, which is then populated by Site objects forming a lattice. Detailed explanation\n", "of what is a Site object and what is a Builder object is given on the\n", "[KWANT website](https://kwant-project.org/doc/1/tutorial/faq)." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import kwant \n", "import matplotlib.pyplot\n", "import numpy as np\n", "\n", "sys = kwant.Builder() # Create an empty Builder object\n", "lat = kwant.lattice.square(a=1) # Create a square lattice object\n", "t = 1 # Set the hopping value\n", "\n", "# Populate the builder with sites\n", "for i in range(8):\n", " for j in range(4):\n", " sys[lat(i, j)] = 4*t # Adding a site with onsite potential 4*t\n", " \n", "kwant.plot(sys); # Plot the system" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# Populate hoppings\n", "for i in range(8):\n", " for j in range(4):\n", " if i > 0:\n", " sys[lat(i-1, j), lat(i, j)] = -t # Add hopping -t along the x \n", " if j > 0:\n", " sys[lat(i, j-1), lat(i, j)] = -t # Add hopping -t along the y\n", "\n", "\n", "kwant.plot(sys); # Plot the system" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 1.2 Creating a system with a custom shape " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Besides adding sites to a system using lattice indices, a faster way to add sites is to use a shape function. The function should return True or False depending on whether the site is inside of the desired shape or not. \n", "\n", "For example, if we want to create a ring system, we would first define the ring shape analytically. " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.2.1 Defining a shape function" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from math import sqrt\n", "\n", "def ring_shape(position):\n", " x, y = position\n", " outer_radius = 15\n", " inner_radius = 7 \n", " radius = sqrt(x**2 + y**2) \n", " in_outer = radius < outer_radius\n", " in_inner = radius < inner_radius\n", " return (in_outer and not in_inner)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.2.2 Populating a ring shape with sites" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ring_sys = kwant.Builder() # Create an empty Builder object\n", "latt = kwant.lattice.square() # Create a square lattice object\n", "t = 1 # Set the hopping value\n", "onsite = 0.\n", "\n", "# Populate a ring system with sites using a 'shape' function \n", "# of the lattice object. Besides the shape function, you also\n", "# need to supply a point inside of the system (any point)\n", "ring_sys[latt.shape(ring_shape, (10, 0))] = onsite\n", "\n", "kwant.plot(ring_sys);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.2.3 Using 'neighbors' method to add hopping" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A faster way to add hoppings is through the predefined 'neighbors' function of the lattice object." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ring_sys[latt.neighbors(1)] = -t \n", "\n", "# ring_sys[latt.neighbors(2)] = -t # Adding hopping to 2nd neighbors \n", "kwant.plot(ring_sys, dpi=100);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 2. Creating and adding semi-infinte leads" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# Define the lead direction\n", "left_direction = kwant.TranslationalSymmetry((-1, 0)) \n", "\n", "# A lead is a Builder object with translational symmetry\n", "lead = kwant.Builder(left_direction) \n", " \n", "for i in range(4):\n", " lead[lat(0, i)] = 4*t # Add sites\n", " lead[lat(1, i), lat(0, i)] = -t # Add hoppings\n", " if i > 0:\n", " lead[lat(0, i-1), lat(0, i)] = -t \n", "\n", "kwant.plot(lead);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "sys.attach_lead(lead) # Attach a left lead\n", "sys.attach_lead(lead.reversed()) # Attach a right lead\n", "\n", "kwant.plot(sys, num_lead_cells=5);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ring_sys.attach_lead(lead)\n", "kwant.plot(ring_sys, num_lead_cells=5);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 3. Energy-momentum dispersion (or subband structure) of a nanowire" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Before you want to compute any quantity using KWANT, you need to finalize your Builder object. After a Builder object is finalized, you can not add or remove sites or change hoppings (the form of the Hamiltonian is fixed). So before computing the dispersion of a lead, you need to finalize it first." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "lead = lead.finalized() \n", "print(lead)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# You need to pass a finalized lead object\n", "kwant.plotter.bands(lead, dpi=100); " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 3.1 How to plot the dispersion in a limited k-range" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "from math import pi\n", "k_values = np.linspace(-pi/3, pi/3)\n", "kwant.plotter.bands(lead, momenta=k_values, dpi=100);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 4. Graphene flake" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = kwant.Builder()\n", "gr_lat = kwant.lattice.honeycomb(a=1)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "def ribbon_sys(pos):\n", " x, y = pos\n", " in_x = abs(x) < 15 \n", " in_y = abs(y) < 10.\n", " return in_x and in_y\n", "\n", "def ribbon_lead(pos):\n", " x, y = pos\n", " return abs(y) < 10." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys[gr_lat.shape(ribbon_sys, (0, 0))] = 0.0\n", "graphene_sys[gr_lat.neighbors()] = -t\n", "kwant.plot(graphene_sys, dpi=100);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 4.1 Graphene nanoribbon as semi-infinite leads" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "a, b = gr_lat.sublattices\n", "v1, v2 = gr_lat.prim_vecs \n", "print(v1, v2)\n", "left_direction = kwant.TranslationalSymmetry(-v1)\n", "\n", "left_lead = kwant.Builder(left_direction)\n", "\n", "left_lead[gr_lat.shape(ribbon_lead, (0, 0))] = 0\n", "left_lead[gr_lat.neighbors()] = -t\n", "\n", "graphene_sys.attach_lead(left_lead)\n", "graphene_sys.attach_lead(left_lead.reversed())\n", "kwant.plot(graphene_sys);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "left_lead = left_lead.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "k = np.arange(0, 2*np.pi, 0.01)\n", "kwant.plotter.bands(left_lead, momenta=k, dpi=100);\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 4.2 LDOS in zigzag graphene nanoribbon" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = graphene_sys.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ldos = kwant.ldos(graphene_sys, energy=1.e-5)\n", "kwant.plotter.map(graphene_sys, ldos, method='nearest', dpi=200, cmap='viridis', oversampling=20);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 4.3 Disorderd graphene nanoribbon" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "def qpc(pos):\n", " x, y = pos\n", " in_x = abs(x) < 15\n", " in_y = abs(y) < 10\n", " in_cut = abs(y) > 7 and abs(y)-7 > abs(x)\n", " return in_x and in_y and not in_cut\n", "\n", "def ribbon_l(pos):\n", " x, y = pos\n", " return abs(y) < 10 " ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = kwant.Builder()\n", "gr_lat = kwant.lattice.honeycomb(a=1)\n", "graphene_sys[gr_lat.shape(qpc, (0, 0))] = 0.0\n", "graphene_sys[gr_lat.neighbors()] = -t\n", "\n", "a, b = gr_lat.sublattices\n", "v1, v2 = gr_lat.prim_vecs \n", "print(v1, v2)\n", "left_direction = kwant.TranslationalSymmetry(-v1)\n", "\n", "left_lead = kwant.Builder(left_direction)\n", "\n", "left_lead[gr_lat.shape(ribbon_lead, (0, 0))] = 0\n", "left_lead[gr_lat.neighbors()] = -t\n", "\n", "graphene_sys.attach_lead(left_lead)\n", "graphene_sys.attach_lead(left_lead.reversed())\n", "kwant.plot(graphene_sys, dpi=100);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = graphene_sys.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "\n", "ldos = kwant.ldos(graphene_sys, energy=0.00001)\n", "kwant.plotter.map(graphene_sys, ldos, method='nearest',\n", " dpi=200, cmap='viridis', oversampling=20);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 5. Graphene nanoribbon in magnetic field as quantum Hall insulator" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from math import exp\n", "\n", "def bfield(site1, site2, B):\n", " t = 1\n", " x1, y1 = site1.pos\n", " x2, y2 = site2.pos\n", " phase = -1j*B*(x2 - x1)*(y1 + y2)/2\n", " return -t*np.exp(phase)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import kwant\n", "t = 1\n", "graphene_sys = kwant.Builder()\n", "gr_lat = kwant.lattice.honeycomb(a=1)\n", "graphene_sys[gr_lat.shape(ribbon_sys, (0, 0))] = 0.0\n", "graphene_sys[gr_lat.neighbors()] = bfield \n", "\n", "a, b = gr_lat.sublattices\n", "v1, v2 = gr_lat.prim_vecs \n", "print(v1, v2)\n", "left_direction = kwant.TranslationalSymmetry(-v1)\n", "\n", "left_lead = kwant.Builder(left_direction)\n", "\n", "left_lead[gr_lat.shape(ribbon_lead, (0, 0))] = 0\n", "left_lead[gr_lat.neighbors()] = bfield \n", "\n", "graphene_sys.attach_lead(left_lead);\n", "graphene_sys.attach_lead(left_lead.reversed());" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 5.1 Landau levels" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "left_lead = left_lead.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "par = {'B': 0.2}\n", "k = np.arange(0, 2*np.pi, 0.01)\n", "fig = kwant.plotter.bands(left_lead, momenta=k, params=par, show=False);\n", "ax = fig.axes[0]\n", "ax.set_ylim(-1.0, 1.0);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 5.2 LDOS of edge states" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = graphene_sys.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "params = {'B': 0.2}\n", "ldos = kwant.ldos(graphene_sys, params=params, energy=0.5)\n", "kwant.plotter.map(graphene_sys, ldos, method='nearest', dpi=200, cmap='viridis', oversampling=20);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = kwant.Builder()\n", "gr_lat = kwant.lattice.honeycomb(a=1)\n", "graphene_sys[gr_lat.shape(qpc, (0, 0))] = 0.0\n", "graphene_sys[gr_lat.neighbors()] = bfield # Add field in the main system\n", "\n", "a, b = gr_lat.sublattices\n", "v1, v2 = gr_lat.prim_vecs \n", "print(v1, v2)\n", "left_direction = kwant.TranslationalSymmetry(-v1)\n", "\n", "left_lead = kwant.Builder(left_direction)\n", "\n", "left_lead[gr_lat.shape(ribbon_lead, (0, 0))] = 0\n", "left_lead[gr_lat.neighbors()] = bfield # Add field in the lead \n", "\n", "graphene_sys.attach_lead(left_lead)\n", "graphene_sys.attach_lead(left_lead.reversed())\n", "kwant.plot(graphene_sys, dpi=100);" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "graphene_sys = graphene_sys.finalized()" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "params = {'B': 0.30}\n", "ldos = kwant.ldos(graphene_sys, params=params, energy=0.5)\n", "kwant.plotter.map(graphene_sys, ldos, method='nearest', dpi=200, cmap='viridis', oversampling=20);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Please send corrections to petrovic@udel.edu" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "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.8.3" } }, "nbformat": 4, "nbformat_minor": 4 }