{ "cells": [ { "cell_type": "markdown", "source": [ "# Monitoring self-consistent field calculations\n", "\n", "The `self_consistent_field` function takes as the `callback`\n", "keyword argument one function to be called after each iteration.\n", "This function gets passed the complete internal state of the SCF\n", "solver and can thus be used both to monitor and debug the iterations\n", "as well as to quickly patch it with additional functionality.\n", "\n", "This example discusses a few aspects of the `callback` function\n", "taking again our favourite silicon example." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "We setup silicon in an LDA model using the ASE interface\n", "to build the silicon lattice,\n", "see Creating slabs with ASE for more details." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using PyCall\n", "\n", "silicon = pyimport(\"ase.build\").bulk(\"Si\")\n", "atoms = load_atoms(silicon)\n", "atoms = [ElementPsp(el.symbol, psp=load_psp(el.symbol, functional=\"lda\")) => position\n", " for (el, position) in atoms]\n", "lattice = load_lattice(silicon);\n", "\n", "model = model_LDA(lattice, atoms)\n", "kgrid = [3, 3, 3] # k-point grid\n", "Ecut = 5 # kinetic energy cutoff in Hartree\n", "basis = PlaneWaveBasis(model, Ecut; kgrid=kgrid);" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "DFTK already defines a few callback functions for standard\n", "tasks. One example is the usual convergence table,\n", "which is defined in the callback `ScfDefaultCallback`.\n", "Another example is `ScfPlotTrace`, which records the total\n", "energy at each iteration and uses it to plot the convergence\n", "of the SCF graphically once it is converged.\n", "For details and other callbacks\n", "see [`src/scf/scf_callbacks.jl`](https://dftk.org/blob/master/src/scf/scf_callbacks.jl).\n", "\n", "!!! note \"Callbacks are not exported\"\n", " Callbacks are not exported from the DFTK namespace as of now,\n", " so you will need to use them, e.g., as `DFTK.ScfDefaultCallback`\n", " and `DFTK.ScfPlotTrace`." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In this example we define a custom callback, which plots\n", "the change in density at each SCF iteration after the SCF\n", "has finished. For this we first define the empty plot canvas\n", "and an empty container for all the density differences:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots\n", "p = plot(yaxis=:log)\n", "density_differences = Float64[];" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "The callback function itself gets passed a named tuple\n", "similar to the one returned by `self_consistent_field`,\n", "which contains the input and output density of the SCF step\n", "as `ρin` and `ρout`. Since the callback gets called\n", "both during the SCF iterations as well as after convergence\n", "just before `self_consistent_field` finishes we can both\n", "collect the data and initiate the plotting in one function." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using LinearAlgebra\n", "\n", "function plot_callback(info)\n", " if info.stage == :finalize\n", " plot!(p, density_differences, label=\"|ρout - ρin|\", markershape=:x)\n", " else\n", " push!(density_differences, norm(info.ρout.real - info.ρin.real))\n", " end\n", " info\n", "end\n", "callback = DFTK.ScfDefaultCallback() ∘ plot_callback;" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Notice that for constructing the `callback` function we chained the `plot_callback`\n", "(which does the plotting) with the `ScfDefaultCallback`, such that when using\n", "the `plot_callback` function with `self_consistent_field` we still get the usual\n", "convergence table printed. We run the SCF with this callback ..." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy Eₙ-Eₙ₋₁ ρout-ρin Diag\n", "--- --------------- --------- -------- ----\n", " 1 -7.844855073505 NaN 1.97e-01 4.3 \n", " 2 -7.850323461935 -5.47e-03 2.92e-02 1.0 \n", " 3 -7.850618768967 -2.95e-04 2.81e-03 1.5 \n", " 4 -7.850646704398 -2.79e-05 1.26e-03 2.3 \n", " 5 -7.850647338276 -6.34e-07 5.86e-04 1.0 \n", " 6 -7.850647502876 -1.65e-07 4.63e-05 1.0 \n", " 7 -7.850647511449 -8.57e-09 6.37e-06 2.3 \n" ] } ], "cell_type": "code", "source": [ "scfres = self_consistent_field(basis, tol=1e-8, callback=callback);" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "... and show the plot" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=1}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nO3deVxU9f4/8PfnnM+wiiIisoloUpo7iIoLuywKIuZaaZppdVusbt/s0Xatbt2sfpXVrVu2mGVmmhsIKCq4peW+ZOaCCqIoouzbnOX3xxipoQLOzJnl9fwLDjPnvD3RvPh8zmdhqqoSAACAvRK0LgAAAEBLCEIAALBrCEIAALBrCEIAALBrCEIAALBrCEIAALBrCEIAALBrCEIAALBrCEIAALBrCEIAALBrNwzCurq6f/3rXxMnTly5cuVN3n/8+PGqqqqmX09RlGZUZ9NwKxooioKl/gxUVcWtMMCtaKCqKj4uGpjiVtwwCJ955hmdTvf222+/884727dvv9HLpk+fvmvXrqZfr7q6Gr/cBrW1tfjlNqirq5NlWesqLEJ9fb0kSVpXYREkSaqvr9e6Cosgy3JdXZ3WVVgERVFqa2uNflp+ox8sW7bs9OnTTk5Ojz766KJFi8LCwox+bQAAAM013iIsLS11cnJycnIiosDAwNOnT5u3KgAAADNpPAg55w29drIs63Q6M5YEAABgPn91jR49enT58uUuLi4TJ0708vIiorKysjZt2hw9evSOO+7QrkIAAAATutIi3L17d2hoaFlZ2f79+4ODg4uLi2fMmPHPf/5z/fr1//3vf6dNm6ZtlQAAACZypUU4d+7cJ5988vXXXyeilJSU+fPnv/zyy4sXL966devXX39999133+j9iqIUFBT88ccfVx/08PBo3759w7e559TFJ5T/DRXZVa85VqY+vUNOi+dXHwQAADAzZpjM4OHhkZGRMWjQICL67LPPlixZsnHjxqa8v1+/fsXFxc7OzlcfnDJlylNPPXX1kTcO8XM17MP++prqKhcXl7xK4YHtDp8PrL+7jf1OpaipqXF0dBQErGlANTU1Op2O8xuOYbYfdXV1giDgqTwR6fV6RVEcHR21LkR7kiTp9frrPmbtk6IodXV1zboVLi4uoije/DWciGpqai5fvuzt7W045O3tffbs2SZeo3Xr1h988EFERMTNX/ZWGL2yW/7nft17fVmR7DLtF2VRtNizrUMTr2KTRFF0cnJCEBIR5xxBaODg4IAgNEAQNkAQNlAURafTubi4GPe0AhEJgsAYa5jRrCjKLfOzBV4LEb2d6cGf+fiNyneRYs+26BMFAADtCUTk6Ojo4eFRVFRkOHTu3DkfHx9TXOyBIOHni+KpShXrqQAAgIW40hmVmJi4cuXKIUOGENGqVasSExONfqVjZer4jfLKyPo3fnMcsFL6T39xVi8BrUIAsAHZ2dknTpww0ckVRcF8bk9Pz7Fjx5ro5Hz48OHffPPN7NmzIyIiKioqLly4kJeX9+CDDxr3MoYU/DZSDNSpK4eLs3Yobx2QM84oCyO5N/q9AcDKffbZZ+Xl5V26dNG6ENtUUlJy/PhxEwbh7Nmz3d3dfX19Dx48mJ6e7uLismDBgtatWxvxGsfK1Akb5e8ixR5tWWUlEdGHYWIrTusK1eAV+q/CeYI/WoYAYN1mzpxpuk9qO7d3716jN8+uxmNjYw1f+fr6zpw50xTXcNXRd1Hi3e7XpN2boWKUryoyemCTPCaQvTNQdMDwSQAAMDtzhI+vC7suBQ2G+7FoX7Y3lZ+qpCGrpWNl9junEAAAtKJ9K8zTiVYNF5/oIQxLlz4/gvGkAABgVtoHocGUIGHjSP7JYWXcBvkydqAEAABzsZQgJKK73dmOFO7rQv1WSFuL0E0KAADmYEFBSEROIs0LE+eFCWM3SHP2yDLSEAAATMyygtAgpZOwf4zulwvq0DTpZAXCEABs1pFStbT+Fq8pq6ffS/FJaEKWGIRE1MGZMhL4xC7CwFXSDycwggYAbFNeBSWvk26ShWX1lLROyqswY032x0KDkIgY0ayeQkYCf2WPMiVXrpK0LggAwNhGdGQv9hVHZEmXGhskWFZPyeukZ3sJIzu2cNWRDz/8sLi4+LZKvG1HjhxZsmRJoz86ceLEggULzFtOIyw3CA36e7Ldo7kzp/4rpb0l6BwAAFuT4M9eCRaT1l6fhYYU/GcvIaVTyz+o58+ff/Hixdspb+LEiXv27LmdM5w4cWLNmjWN/ujkyZM//PDD7ZzcKCw9CInITUefDRVfDxFGZElz9ysK0hAAbMvfs9AoKXid6urq8vLy6w6Wl5cbtme/kcOHD1dUNK9nVlXVkpISSbrSjzdy5MiFCxc2/LSioqKuzrImyVlBEBqM7SzsSOFp+Up8lnSuWutqAACM6uosNHoKqqr63HPPdenSpU+fPsOGDTNsuvfLL7907969f//+Pj4+77//PhEtXbp0woQJDe/y9PR88cUXjx8/PnPmzP79+y9evLgp1zp06FCPHj2GDh3q4+Pz6aefEtGKFStGjx594MCBPn36PP744/369fPy8vrPf/5jlH+aUVjTnuCdWrHckfzf++TgFfovwnmLO80BAMwg55xa08zBDXH+wpA0iYjGd2E6gWUUNKMHzJlTlE/jn4qLFy/OycnJy8tzdnZ++OGHn3vuufnz50+YMOHVV1994IEHTp06FRISEhoaWl9fX2nYGIGIiC5fvvzMM8+kpaV99NFHERERTalBVdVx48a98sorkyZNOn78eEhISFRUlF6vr66uliTpwIEDs2fP/vjjj48fP967d28TrW7dAtYUhETEBZoTLEb6CFNy5dRA9vYA0VHUuiYAgMYsPKacr2nesxxJobNVqkr083naWSw3671eTizKp/EPxJUrV06fPt3FxYWInnzyycjIyEOHDtXU1EyZMoWIAgMDx44du2bNmp49ezbxWtXV1UePHiUinU7Xo0ePhuNHjhy5dOnSxIkTiahr164JCQlZWVm+vr6Gn3p4eNx7772GH/n5+eXl5TXrH2g6VhaEBpE+7OA9/JFtcv+V0uJosWdbNA0BwOJ8Hd68v9MNPaILI0VHkb22R06P5x6Oxqnk/Pnzbdu2NXzdtm3b0tLSkpKStm3bMnblw7Ndu3YlJSVXv0VV1Zs8Ozx9+vQTTzxBRJ6enitWrGg4fuHCBXd394bTenh4XL58uSEIr97gz9HRsb7+VjMozcVqnhFep40DLY4S/6+3ELVGmncIEw0BwLpd/VzwRuNIW0xV1YMHDxq+PnDgwB133BEUFJSfn19aWmo4uG/fvqCgoLZt2zbE4bFjxwxBqNPpZPn6tmn37t23bNmyZcuWq1OQiBRFKSgoKCsrM3y7f//+rl27GuffYEpW2SJsMCVIGNCeTcqRNxep84eJxvrrCQDAnP4+OibBnxGJSWslY7ULFy1aFBQU5OXl9fTTTz/11FOdO3ceOXLk5MmTn3rqqW3btu3evfvbb78lot9///2jjz7q3Lnz559/bmjY9erV69NPPy0qKgoNDQ0KCrrlhRwcHB566KHHHntsw4YNhYWFY8eOTUtLM8I/wJSstUXYoJs72zGKB7Wh4BXSZizVDQDW5kZjRI3bLnz11VcPHz68aNGil19++dFHHyWi77//PiYm5quvvqqqqtq1a1e7du3atWuXlZW1d+/eVatWzZ07d/bs2c7OzvPmzYuKijp48OD58+ebcqHAwMCHHnpo/vz5FRUV27Ztc3Z27t69+8SJEzt06PDII480vGz69On+/v5G+IcZg3W3CA0cRXorVIzxVe/Lke/ryl4PEXVWn+8AYBduPlPCiO3CNm3azJ079+ojOp3uqaeeuu5lYWFhYWFhhq8bZjj84x//aNa14uPj4+PjG77t1atXr169iGj27NkNB59++mkiOnbsWLPObCK2kxjD/dju0fzgJXVYunSiHE1DALACy04q/9f7ZvMFE/zZC33F5acwEsKEbKFF2MDLmdLj+YeHlLDV0vuDxPu62k7MA4BNmn7XrT+mkgIYUQvHxr/33nv+/v6bN29u2dubKyoqat++fU1/fe/evV966SXT1dNENhWE9OdS3cO82aQcee0Z9ZMhYiud1jUBAGhk+PDhWpdwM15eXl5eXlpXYUNdo1cL9mT7UnlbR+q9XNp+Ad2kAABwQ7YZhETkzGlemPjuQGF0tjRnj4ylugEAoFE2G4QGYwKFnSk856walymdrUYYAgDA9Ww8CIkooBXbOJIP9xP6r5TS8jHyCgAArmH7QUhEIqPZfYQVsfzpHcrDW+VqbHYPAAB/srVRozcx0IvtSeWPbpVDV0qLo8XeHliqGwCM5syZM7/99pvWVdgmU+9TYUdBSEStdbQoSlx4TInJkF7qKz7ZU0AYAsDtmzx58gsvvDB//nwTnV9V1YYtHexT9+7dTXdy+wpCgylBwkAvdm+OnHtO/SJcbIelugHg9qSkpKSkpJjo5JIk6fV6Z2dnE50f7OIZ4d/d1Yb9Mor3aUe9ftKvPYPRpAAA9stOg5D+3Oz+20g+fYs8a7tcj/GkAAB2yX6D0CDGl+1N5XkV6tA06TiW6gYAsD/2HoRE1N6JVsfxh+4ShqZJ3x5HwxAAwL4gCImIGNHMbsKGEfydA8r4DXJpvdYFAQCAuSAI/9KjLfs1hfu4UL8V0rbz6CYFALALCMJrOIk0L0x8f5Bwz3ppzh5ZRhoCANg6BGEjRncSdo3muefU4RlSYRXCEADAliEIG+fvynJG8pROQr8V0pI8jKABALBZCMIbMmx2vyaev7RLmZIrV2GpbgAAW4QgvIXQ9mx3KleJQldK+0rQTQoAYGsQhLfWWkffRoqvhgiJWdK8Q9jrHgDApiAIm2pcZ2H7KP7jSSUhUyqq0boaAAAwEgRhMwS6sU0jeVgHFrxCn1mAliEAgC1AEDaPYanuRVF85lYs1Q0AYAsQhC0R5cP2pvJTlTR4tXS0DE1DAAArhiBsIU8nWjVcfLKHMCRNmncIDUMAAGuFILwtU4KEzUn866PK2A3y5TqtqwEAgOZDEN6u7u5sRwr3c6F+K6QtRegmBQCwMghCIzAs1f35MHFSjoylugEArAuC0Gji/Nju0fyXC+rQNOlkBcIQAMA6IAiNqYMzZSTwiV2EgaukxScwggYAwAogCI3MsFR3RgL/1x4s1Q0AYAUQhCbR35PtHs2dOfVfKe3FUt0AABYMQWgqbjr6bKj4eoiQkCXN2SNjrW4AAMuEIDStsZ2FX1P4hrNqfJZ0rlrragAA4G8QhCbXqRXLGcGHdGDBK/Sf/K789/D1g2hqZXr2F1nC2BoAAC1wrQuwC4aluhP8hftyZC5QaR292O/KnyD1Ck3YKCf6M46/SQAAtIBPX/MZ5MX2pPJ+nuydg8rT22Uiqldo3AY50Z890h3/IQAAtIEWoVm1caAfosQFR5VHtslHLnMdV0Z0FJCCAAAaQhBqYOqdQk8PNjRN8ndRl8ciBQEAtIRPYQ3UK/T6XuW1vnKNTN2XSeV6rQsCALBjNwvCsrKyPXv2VFdj1L8xNTwXfPwu5cR4kRF1+1EqrMI0QwAAbdwwCBcsWBAVFZWQkHDkyBFzFmTbrhsd4yDQ7+O4jyv1Xi4fKUUWAgBo4IZBOHXq1D179vTq1cuc1di8E+XqPYHXjBEVGe0YxSO8WWQ6tjMEANAAnhGaVXd3NiXo+nuuE2j5cHFRNB+7QVqSh3n1AABmxYlo5cqVGRkZVx9NSkoaNWqURiXZqRhftj6RJ62TT1fSc73xBwoAgJlwIho4cGBAQMDVR318fDSqx6718mA/J4sj18qFVer7g0SBaV0QAIAd4Bs2bIiJiUHyWQg/V5Yzkqeul8ZtkL+LFJ0xzxMAwMSEqVOnSlIju8du2bKlf//+ZWVlM2fOfOSRR8xfmd1q60hrE7ijSDEZ0sVarasBALB1zNfXd+HChTExMS17f3BwcH19fevWra8+OHbs2BkzZjT6+qqqKhcXF8bQ60c1NTWOjo6C0PjjQJXord/4inzxp4j6ji42Ppq0pqZGp9NxjvYv1dXVCYKg0+m0LkR7er1eURRHR0etC9GeJEl6vd7Z2VnrQrSnKEpdXV2zboWLi4soijd/DQ8KCsrLy2txEDo7O0+bNq1v375XH/T393dzc2v09YwxV1dXBCERiaLo5OR0oyAkojcHkbebErdRWD1cDPa05TvGOUcQGjg4OCAIDRCEDRCEDRRF0el0Li4uxj0td3V1LSsra/n7Oe/du/ewYcOMWBM0eLKH0NGVErKkbyJ4YkdbzkIAAK0Ily9f7tChg9ZlwA2lBgqr4/iDm6X5RzDFEADA+ISDBw/27t1b6zLgZgZ5sS3J/J2Dypw9sta1AADYGqFPnz59+vTRugy4ha6t2eYknp6vTt8sS2gZAgAYj7Bs2TKta4Am8XamTUn8fI06cq1UgZ2bAACMRPD29ta6BmgqV06r4nigG4vJkC7UaF0NAIBNwJqWVkZk9NlQ8Z5AIWy19EeZjc8vBAAwAwShVZrdR3i5nxCRLm3Fzk0AALcHQWitpt4pfBvJx6yXlp7E4BkAgJZDEFqx4X5s/Qj+7C/Kfw8jCwEAWghBaN16e7DckeLHh5VZ22UFvaQAAM2HILR6nd3YtmS+t0SdsFGuxYR7AIBmQhDaAg9Hyk7kIqMRWVJZvdbVAABYFQShjXAUaXG02L89G5Im5VeikxQAoKkQhLaDEb09QJxxlzAkTd5XgiwEAGgSBKGtmdVT+GCQEJcprT2DLAQAuDUEoQ26p7OwcjiftlladBzTKgAAbgFBaJsGd2DrR/AXd2HnJgCAW0AQ2qy73dn2UXz1aXXGFuzcBABwQwhCW+bjQpuT+JkqdewGuVrSuhoAAIuEILRxrXSUFsfbO1F0hlRcq3U1AACWB0Fo+7hAnw8TE/xZ2GrpGHZuAgC4FoLQLjCiOcHiC32FYenSz+eRhQAAf0EQ2pEH7xQWRvJ71kvp+chCAIArEIT2Jc6PrY7jM7dK//sdA0kBAIgQhHYotD3bmsw/OKTM2i6jYQgAgCC0R13c2M+j+O6L6tRNsh4tQwCwbwhCO+XhSOtH8GqJErFzEwDYNwSh/XIS6Ydo8S53NixdOlOFXlIAsFMIQrsmMvrvYHH6nULYavnAJWQhANgjBCHQrJ7Ce4OE+Exp0zlkIQDYHQQhEBGN6ywsiuLjN0qLT2DwDADYFwQhXBHtyzaM4M/vxM5NAGBfEITwl55t2fZR4qrT6qztsoJeUgCwDwhCuIavC9ucxI+UqmM3yDXYuQkA7ACCEK7npqO0eO4sUnSGdBE7NwGArUMQQiMcBPouSoz3Z2GrpePl6CQFAFuGIITGGXZumt1HiFoj776ILAQAm4UghJt56C7ho8FCYpaUUYAsBADbhCCEWxjdSUiL4w9tkT4/gimGAGCDEIRwawO92JYk/u5B5fmdmGIIALYGQQhNckdrtiWJbyhUp23Gzk0AYFMQhNBUHZwpN4lfrFVHrpXK9VpXAwBgJAhCaAZXTiuH8ztas2FpUiF2bgIAm4AghOYRGX06RLz3DmFYunykFFkIAFYPQQgtMbuPMCdYiFwjbSlCFgKAdUMQQgtNCRK+i+T3rJeW5GHwDABYMQQhtFysH9swgv/fL8o7B5CFAGCtEIRwW3p5sM1J4tdHFezcBABWCkEItyvQjW1L5vtK1PEb5VpMuAcAa4MgBCNo60jrErlOoJgMqaRO62oAAJoDQQjG4SjS91HiMG8WniblV6KTFACsBoIQjIYRvRUqzuwmDE6T95YgCwHAOiAIwchm9RQ+ChPiM6WsM8hCALACCEIwvtRAYXUcn7ZJ+uIPTKsAAEuHIASTGOTFNo7kb+xT5uzBQFIAsGgIQjCV7u5s+yielq9O3yxLaBkCgKVCEIIJeTvTppG8qEYduVaqwM5NAGCREIRgWq10tGo47+TGYjKkCzVaVwMA8DcIQjA5LtBnQ8URHVnYaumPMgwlBQDLgiAEc2BEc4LFl/oJEenStvPIQgCwIAhCMJ9pdwrfRvLUbGnpSQyeAQBLgSAEsxruxzIS+FPblU8OK40OJcX4UgAwMwQhmFt/T7Y1WXz3oNJ5ib6s/pofnapQw9OlGkmjygDALiEIQQOd3diu0dxNR0E/6otrrxzMr1Tv2SB/OFh05poWBwB2BkEI2vBwpL2puqA27M6l+rM1LL9STV0vfzZU7O/JtC4NAOwL/vYGzTiKtDmJJ62Veq8W/FxoSYzYrx1SEADMDS1C0JLI6NMhoiungkoatU7+xzY564xah9VJAcCMbhiEX375ZWhoaHBwcGpq6oULF8xZE9gPw3PBFVHKwgjyc6XAVuz/HZS9vtMnr5M+P6Kcq9a6PgCwAzcMwh49emzatGnPnj0hISGvvvqqOWsCO9HwXDDYQ03txJ7rLWScUZbH8ryJuilBwtYitcdP+v4rpTl75N0XMQcfAEzlhs8IBw0aZPiiS5cueXl55qoH7MWZKnV0tvxluNivHaupISIaEyioKo1aJ2XE83GdhXGdSVLEHRfU9ALl3hxZr9BwP5YUwOL9BQf06AOA8XAiKioqOn369NVHO3ToEBgYSESXLl164403lixZoklxYMM8ndj3UWI392tGx9zTWejahjVMn+ACDfVmQ73Ft0Lpt8tqer46d79yX44c4cPGdRaSA4S2jhpUDgA2hhPR4cOHly1bdvXRiIiIwMDAsrKy5OTkN998s2fPnhqVBzbLSaTrUtCgj0fjA0d7tGU92rLZfYTiWsosUNLz1Se363u2ZckBwqhOrHtjpwIAaApeXFwcHR0dHR193Q+qqqpSU1Mff/zxlJQUTSoDaFR7J5oSJEwJohpJ3HZeTctXhmcqjgIlBbDkACHSh3F0nAJAcwhJSUmN/uAf//hHYWFhbm7uww8//N5775m5LIBbcuYU68fmhYlnJvHVcWJbR3p+p+zzvX78BnnhMaUc+wADQNOwdu3arV69evDgwdf94Pjx42VlZYav3dzc7rzzzkbfHxYWFh4eftddd119sGfPniEhIY2+vrKy0tXVlTF0ZFF1dbWTk5MgoP1CNTU1Op2OcyMs73C6Ul1XSOkFau456uNBYwNZaiB1dLWa37e6ujpBEHQ6ndaFaE+v1yuK4uiI58AkSZJer3d2dta6EO0pilJbW+vi4tL0twiCcMvE4QMHDty1a9ffg7Br165NuYZer9+7d29hYeE1J+X8Ro8V6+rqOOcIQiKqq6tjjCEIiaiurk5RFFk2wkR6bx1NCaQpgVQtUe55llkozD3A2jrQCH8l0VcN81It/DfPEISKgj04rgSh1lVYBEMQ4rOCiBRFqaurE0Wx6W9xcnK6dRC2a9euqKioxWW5urq++OKLERERTXy9oiguLi4IQgO0CA0YY8ZqETZwIRrbmsYGkazSvhI1LV95fp+aX6km+AvJnViiv9DKIhtdoiiiRWiAFmEDtAgbKIoiCEKzWoRNwaurq1u1amXckwJYDpFRiCcL8RTnBFNehZp2Wv38iPLQZnmAF0vqKIzrwnxd8GcZgF0TTp061bFjR63LADCHLm5sVk8hO5GfnKib2U3YfVHt9ZPUsHgNVq8BsE/8yJEjycnJWpcBYFYejmRYvEZWxe3n1fQC5b4cuVameH8sXgNgd/gHH3zg7u6udRkA2hDZX4vXGDpOP/xNuT9HDvdhyQFCSiehA57LANg6/tBDD2ldA4BF6OLGZvVks3oKF2spo0BJz1dn/6q/ozVLCmDjuwh3Y/EaABuFjXkBruf55+I1tbK4tUhNy1fiM2WHPxevifBhOnScAtgQBCHADTmJFOvHYv3EeWFXVv2es0c+UqpG+wpJASylk9DGQesSAeC2IQgBmqRh1e8LNZR1Rll6UvnHNrlvO5YcIKQGsjvboOMUwFohCAGax8vZ0HEqVEu04aySnq9GrpGdRZYUwMZ1FgZ3YAIyEcCqIAgBWsiFU3KAkBxAn6ri3hI1LV95aofcsHhNgr/ghvVhAKwBghDgdglXLV5zskLNLlQXHvtr8ZqxnZmf9az6DWCHEIQAxtTZjc3sxmZ2E6ok2nhWWZqnvrpX9nFmyZ1YUkdhiHfjy+zWyeT4t2WEGz0IAEaHYeAAJuHKKTlAWBgpFt+v+2yoSEQztsqBP0gPb5XT8pW6q3baOFdNEWukCzXXvL2gSh2WLpXWm7doALuEIAQwLcPiNW+Fir+P5Tkjxbvd2Ye/KV7f6ZPXSZ8fUYpqyMeF3h0gjlwrnau+8paCKnV0tvzxYNEd0zMATA9dowDm07B4TUkdbTyrpJ1Wn9+p7+LGkgLYI92FUeukZZFMYTQ2R/50iDigPZ4sApgDghBAA+3+XPW7VhZzz6mrTytf/aGoRKFprJUDLYpECgKYD7pGAbTkJFKCP/tkiJg/ic8fJqiMimooMUsau0H+8g/lbDX2hgIwObQIASxCQZX64i5lRZQiqcLTv7I4P/ZrsfrSLtnDkSV3YrG+WOMUwFQQhADaM4yO+XSI2Ke1LAj00WDx6R3y6jj+yRBxX4m6vlCds0c+cEmN8GHJAUJSAPN1QccpgNEgCAE0dr6GRmfL/xsihrZndXVEREO92TsDxdHZ0rpEHuLJQjzZ7D5CcS3lnlPSTquzf5V9Xa40EyN9GEczEeD2IAgBNObpRD9EiUHXLtsd6cMWRoqtr5o+0d7pyvgaWRX3lahp+crzO+XTlWqUjxDrx0Z1EryxhzBAiyAIATQmMgpqbPOKu26wo4V41Ypu52to7RklPf/KNIxYv5utXwMAjUIQAlixDle2wiBJEXdcUNMLlKd2yAVVaoS3kBTARnUSMCUf4JYQhAC2gAs01JsN9RYplPIq1PWFanq++uR2fdfWLCmAJQcIwZ5oJgI0DkEIYGu6XFn4m2okcdt5df1Z5f5c+XKdGucnJHdi8f5Ca+wPBXAVBCGAzXLmFOvHYv3Et/5sJi48pkzbJPdtx5IDhFg/FuKJViIAghDAPvzZTBSqJfr5vJqWr6RmKzrBkJQs0V9ohWYi2CsEIYB9cfmzmTgvjPIq1LTT6udHlBlb5ND2LNZXSOnEurmjmQj2BUEIYBTnaNIAABkFSURBVL8adsMwbCOcnq/GZiqOAsX6saQAFucnYGdgsAcIQgC4so1wcgAR0W+X1fR89cPflMm5cmh7ltRRGNOZdXRFMxFsFoIQAK7Roy3r0ZbN7nNl08T1hepbK2UXfmUaRrgPc8CibmBbEIQA0LiGTRM/HSLu/dva3yM7Mj80E8EmIAgB4BaEK4u6Ye1vsE0IQgBoBqz9DbYHQQgALYG1v8FmIAgB4HY1uvZ3fqUa6YO1v8EKIAgBwGiuXvv7ZIWaXaim56uztuvvwNrfYMEQhABgEp3/tvb35Fz50k3X/lZU2n9J7dfu+qw8cEnt7s50GI8DpoHfLAAwLcPa32+FiofH8p9H8aHebOExxXeRfmiaNHe/svui2vDKSome+FnOKFCvfvvmIvWhLXJZvdnrBruBFiEAmM/f1/4es14RGQ3/c+3vNfE8aZ0kKUKiLxHR1iL1nzvk1XHc00nr0sF2IQgBQAM3X/v7tb2KLJO7jj23W14dx31ctC4XbBqCEAA01rD2d2k9ZRcqmQVqQaU6aRN5OLK9qbwDJiaCieEZIQBYCncHGtdZ+Cpc/DGWd2pFl2rpxV2K1kWB7UOLEAAsy9Yi9dkdcnYcVejVAelyaZ2yLBafVGBCaBECgAXZWqQ+vePKc8Gubur+MbqsQnXMeknrusCWIQgBwFKU1tPsnXJ6/F+jY7q2pkNjdOsL1Qc2yYp60zcDtBSCEAAshbsDbU66fnRMoBsVTNKdqVKnbJIlPDEEE0AQAoAFERtbga2NA6XF8eIa9f5cWY8sBGNDEAKAFXDhlBbP62Qas16qlbWuBmwLghAArIODQEtjRHcHNiJLqtRrXQ3YEAQhAFgNLtCCCDHQjY1YK5UjC8FIEIQAYE1ERl+Gi/3asZg1Ukmd1tWATUAQAoCVYUQfhIkRPiw2Qyqu1boasH4IQgCwPozo3YHixC5CeLpUWIUJhnBbEIQAYK1m9xGm3SkMS5fzKpCF0HJYwQ8ArNhzvQU3HcVkyOsSxKA2jU1CBLgVBCEAWLdHuwtuOorKkDPjxV4eyEJoNgQhAFi9+7sKnFFsppQex0PbIwuhefCMEABswcQ7hC+H8eR10vYLeF4IzYMgBAAbkRTAFkTwUeukDWeRhdAMCEIAsB0J/uynWH5fjrT2DLIQmgpBCAA2JdybpcXxqZuk5aewUQU0CQbLAICtCW3PMhK4YW3uKUH4cx9uAUEIADaoXzu2cSSPz5QllR68E1kIN4MgBADb1N2dbUoSYzPkinqa1RNZCDeEXw4AsFmd3VjOSPHjw8rre/G8EG4IQQgAtiygFduSzH/MU57fiY3toXE37BrdsWPH+vXr6+rqhgwZkpCQYM6aAACMyNuZNozgcZmSospzB4hYeAauc8MWYUFBQY8ePQYNGvTqq68uWbLEnDUBABiXlzPlJvEtReqjW2UFMwzhWjdsEY4bN87wRV5e3oEDByZMmGCukgAAjM/dgdaP4KPWSVM2yQvCRY7nQvCnm40a/f7777du3Xr06NFvvvnGbAUBAJiIK6e0OJ6aLd2XK38XKeqQhUBEhq7RV155JeRaK1euJKJBgwalpKQQ0caNGzUuEwDAGFw4pcVzvUKp2VItRs8AERlahK+99tprr73295916dKlS5cuTk5O77777uTJk81eGwCA8TkI9GO0OG2zPCJLWh3HW+m0Lgi0JixbtqzRH6xYseLkyZPHjx//3//+FxYWZuayAABMhwu0IELs7MZGrJXK9VpXA1oTnnvuOUVpZKppfn7+s88++8ILLwwbNuy5554zf2UAAKYjMvoiXOzXjkWvkUrqtK4GNMXc3Nw2bdrUr1+/lr2/T58+Bw8eZOyamTlPPPHEv//970ZfX1VV5eLict3r7VN1dbWTk5Mg4Hk91dTU6HQ6zrHgH9XV1QmCoNOht470er2iKI6Ojia9ikr00j6ee15YFan3dLTQeRWSJOn1emdnZ60L0Z6iKLW1tS4uLk1/i4uLyy0/Zrm3t/eZM2daHITu7u4bNmwIDw+/+qAgCDeJOldXVwQhEQmCgCA0EEURQWig0+kQhAbmCUIimjeU5u5XEnOE9Ymin6slfjQhCBsoisI5b1YQNgXX6XT19fW3cwpBEERRNFZBAABmNruPwBgNS5fXjxC7uFliFoJJ8eLiYl9fX63LAADQ0nO9hdY6il4jZyeKQW2QhfaF19fX9+7dW+syAAA09kh3oZWOojLkzHixlwey0I7wJ554wtXVVesyAAC0d39XQSdQbKaUFscHtEcW2guh0an0AAD2aUIX4cthfNQ6afsFCx1ECkZ3s+GdAAB2KCmALYjgo9ZJG84iC+0Cxu4DAFwvwZ/9FMvvy5HWnkEW2j4EIQBAI8K9WXo8n7pJWn6qkbW3wJZgCjMAQOP6e7KMBD4iS6rU05QgNBtsFoIQAOCG+rVjOSN5XKasV2j6XchC24QgBAC4mW7ubFOSGJshV+jpqZ7IQhuE/6gAALfQ2Y3ljBT/e1h5fS+eF9ogBCEAwK0FtGJbkvmPecrzO7Gxva1BEAIANIm3M20YwdeeUZ/4WcakCluCIAQAaCovZ8oZyXdfVB/dKisIQ1uBIAQAaAZ3B8oewY+Vq5NzZQlPDG0CghAAoHlcOaXF8Yu16n25sh5ZaP0QhAAAzebCKS2e6xVKzZZqMXrGyiEIAQBawkGgpTGihyNLzJIq9VpXA7cBQQgA0EIio68jxC5uLHGtVI4stFoIQgCAlhMZfREuhniy6DVSSZ3W1UCLIAgBAG4LI/pgkJjYkYWnSeeqta4Gmg9BCABgBK+HiFOChOgMqbAKEwytDIIQAMA4ZvcRHrxTGJYu51UgC60Jdp8AADCa/+stuOkoeo2cnSgGtWFalwNNgiAEADCmR7oLrXQUlSFnxou9PJCFVgBBCABgZPd3FXQCxWZKaXF8QHtkoaXDM0IAAOOb0EX4chgftU76+TyeF1o6BCEAgEkkBbDF0XzMeml9IbLQoiEIAQBMJcqHLY3h9+ZIq09jcW7LhWeEAAAmNMybZSTwUeskvUL3dEbbwxIhCAEATKu/J8uI54lZUqVEDwQhCy0OghAAwOT6tmM5I3lcpiwpNP0uZKFlQRACAJhDN3e2KUkcnilX6OmpnshCC4L/GAAAZtLZjW0cIX7yu/LaXoydsSAIQgAA8wloxTYn8WUnled3YmN7S4EgBAAwK29nyh3JN55VH/9ZxgRDS4AgBAAwNw9HWpfI91xUH9kqKwhDrSEIAQA04O5A2SP4iXJ1cq4s4YmhphCEAADacOW0Oo6X1Kn35cp6ZKF2EIQAAJpx4bQ6jksKpWZLNZLW1dgrBCEAgJYcBPoxRvRwZCPWSpV6rauxSwhCAACNiYy+jhC7uLHEtVI5stDsEIQAANoTGX0RLoZ4sug1Ukmd1tXYGQQhAIBFYEQfDBJHdGThadK5aq2rsScIQgAAC/JaiDglSIhaI52ubGSCIWbgmwKCEADAsszuI0wJEu5aKm0pumZSRZVE8ZnS8XKEoZFh9wkAAIvzQl+hSlJjM+TMBAr3IiKqlmh0tjTtTqFra6Z1dbYGQQgAYIne6C+6cErIkpfHsDAPmpgrTQ0S7uuKbjzjQxACAFioF/uKAtGYDUq3NuJL/YTxXZCCJoHbCgBguZ7sKQa60uFSdqICjwZNBS1CAAALZXgu+FJfcmLy1G206Zz6QzR3d9C6LJuDFiEAgCWqliglW5oaJNx7B0vuqOYk8T0X1eAV0r4SNA2NDEEIAGBxamVKXifNuOuv0TED27PMBM4YxWVK8w5hrwpjQtcoAIDFcRTpnQFisOc1MyVCPNnq4aJINHWLvKlI/SpcRDepUaBFCABgcRjRdSlo0KMt69aWbRrJO7pSvxXSLxfQTWoECEIAACvjKNK8MPG9gcLobHSTGgGCEADAKqUGCjtS+OI8JTVbvowNK24DghAAwFp1asU2jeQBrSh4pbQD3aQthSAEALBihm7S9wcJqdnSvEMKwrAFEIQAAFZvdCdhRwr/IU8Zg27S5kMQAgDYgk6tWO5IHtCK+q1AN2nzIAgBAGyEoZt0XpiQki3N3Y9u0qZCEAIA2JSUTsLOFL7yNEaTNhWCEADA1gS0YluSeDd36rdC2o5u0ltBEAIA2CAu0Fuh4rwwYTS6SW8FQQgAYLMM3aSrTiujs+VL6Ca9gVsEYW1tbUVFhXlKAQAAowtoxTYn8e7uFIxu0hu4WRDKshwdHZ2SkmK2agAAwOgM3aQfhgmjs6U5e2SE4XVuFoTvvvvu0KFDzVYKAACYzqhOwq7RPLtQRTfpdW4YhEePHv31118nTJhgzmoAAMB0OrqyTSN5v3YUvEL6+TxahldwIvr6669LS0uvPjp58uRZs2Z98sknly5d0qgwAAAwPi7QnGAxxFNJXS892l14pZ8oNLLvoX3hROTh4aHT6a4+um7duvLy8h9//LGwsDA/P3/x4sWTJk3SqEIAADCy5ABh12g2aaO856K8IEL0cNS6IE3x+vr6vw+HOX36dEM0Ojk5+fr6mr0wAAAwoY6uLHck//c+OXiFtChKHNLBfhuGwttvv/33o506dRo3bty4ceNiY2O9vLwiIiLMXxkAAJiUoZv048HiPeulOXtku511zwIDA0+ePNni9w8YMKB9+/be3t5XH4yNjU1NTW309ZWVla6urozZ758eDaqrq52cnAQBaxpQTU2NTqfjnGtdiPbq6uoEQbjuUYV90uv1iqI4Otp3nx0REUmSpNfrnZ2dTXeJM9U0ZYvQ1pHmD1Y9HCw3DxVFqa2tdXFxafpbHBwcbvkxy0+dOlVTU9PiWyyKYrdu3YKCgq4+GBgYeKP/k3U6nU6nQxDSn7cCQUhEkiQhCA0URUEQNlAUBbeCiAwfmCa9FZ3bUM4IevOAGpZB30WwwV6mu9RtURRFluVm3YqmxA0nokuXLvn5+bWsLAcHh1GjRjW971QURVEUEYT0561AENKft0IURa0L0Z7hVwK3gogURWGM4VYQkaqqiqKY+laIIr0aQsO81Yk58oxuzDJHkxp+JYx+KwQi8vDwMO5JAQDAGsX6sV9SxA1n1fgs6XyN1tWYixAQEGDSrmcAALAi/q4sZwQf0oENWCVtLbLc54VGJEybNk3rGgAAwIIYRpN+OUycmCPbw2hS4fnnn9e6BgAAsDixfmzHKHHDWTUu08a7SQUnJyetawAAAEtk6CYd6s0GrJK22G43KYYsAgDADRm6Sb8KFyfZbjcpghAAAG4hxpf9MkrceFaNy5SKbK6bFEEIAAC35ufKckbyod4seIV+faFNNQwRhAAA0CQioznB4neRfOpmm+omRRACAEAzRPuyPaP5z+fV4bbSTYogBACA5vFypswEPsybBa/QZ1t/NymCEAAAms3QTbooij9o/d2kCEIAAGihKB+2ezT/+bwamyGdq9a6mpZCEAIAQMsZuknDfVjISv066+wmRRACAMBtaegmnb5ZnrNHlq0tDRGEAABgBIZu0u3n1eHW1k2KIAQAAOPwcqasRB7vLwSvsKZuUgQhAAAYDSOa3UdYHM0f3Cw/v9M6ukkRhAAAYGSRPmxvKt9XosZmSGerLT0MEYQAAGB87Z0oM4En+AshK6S1Zyw6CxGEAABgEg3dpNO3WHQ3KYIQAABMqKGbNMZSu0kRhAAAYFqGbtLUTkL/lVKW5XWTIggBAMDkGNGsnsLiKP7QFnnWdlmvaF3QVRCEAABgJhE+bG8qP1qmDs+0oG5SBCEAAJhPeyfKsLBuUgQhAACYlaGb9IdoPsMyukkRhAAAoIFwb7Y3lR8rV2MzpMIqLZuGCEIAANCGpxOtiedjAoXQVVJmgWZZiCAEAADNNHSTztyqWTcpghAAADRm6CY9Xq4OS5dOVZi7aYggBAAA7Xk6UXo8n9RFGJwmZZi3mxRBCAAAFsHQTbpyOH/8Z7N2kyIIAQDAggxoz35N4cfL1aFpZuomRRACAIBlMXST3nuHMDhNWlOg7ixuJA7/KFPL9ca5HIIQAAAsTkM36RPb5Akb5U9/v6af9LfL6qSNcnGNcdqLCEIAALBQA9qznam8mzu9sEt5fY9sOHikVL0/V14UJd7RmhnlKtwoZwEAADCFdo60Jp6/f1B5cZd8rlqY2ZVN2y5/HyV2dzdOChKCEAAALBwjeqaXEOLJErLk9AJxwwgxqI3RUpDQNQoAAFbB04mC3KhLK8o5Z+ShpGgRAgCApTM8F1wczfwcpKnbRUlR/nG30RpyaBECAIBFO1KqTsq58lzQQaAfY8R1heonh4023x5BCAAAlktS6JFt8o/Rf42OcRBoSbSYdUbdW2KcPlJ0jQIAgOXiAm0YwcVrB8c4irRiuCgaacSMuVuEBQUFlZWVZr6oZSosLCwvL9e6CotQVFR0+fJlrauwCMXFxRcvXtS6Cotw8eLF4uJirauwCJcvXy4qKtK6Ci01BF55eXlhYeF1B2+fuYPwySef3Lp1q5kvapleeOGFzMxMrauwCG+88cbSpUu1rsIizJs376uvvtK6Covw1VdfzZs3T+sqLMLSpUvfeOMNrauwCJmZmS+88ILRT4tnhAAAYNcQhAAAYNdud7CMJEn79+9v+usvX7588OBBFxeX27yuDSguLj58+PCmTZu0LkR7RUVFx44dw60gooKCgtLSUtwKIjp16lRFRQVuBREdO3asqKgIt4KIDh8+XFxc3Kxb0b9/f1dX15u/hqnqbQ0/Xbx48ccff8x5UwO1srLSycmp6a+3YVVVVQ4ODjqdTutCtFdTUyOKooODg9aFaK+2tpYx5ujoqHUh2qurq1NV1cnJSetCtFdfXy/LsrOzs9aFaE+v19fX198y2K725Zdfdu3a9eavud0gBAAAsGp4RggAAHYNQQgAAHYNQQgAAHYNQQgAAHbNrKM3L1++nJeX5+/v36FDB3Ne19IoirJr166TJ0/6+PgMHTpUEOz3z5Hy8vJdu3ZduHChXbt2ERERGDhKRHl5eVVVVb169dK6EM1cvHjx9OnTDd92797dzidc/f777/v373d3dx8yZIibm5vW5Wjj1KlTJSUlDd+Koti3b19jndx8QThq1Ki1a9dyzv/9738//fTTZruuBRo8eHB1dXXPnj0PHDjg5OSUk5Njt7/cjz322Llz53x9ff/444+SkpItW7b4+PhoXZSWCgoKQkNDW7VqdXUS2Jvly5e//PLLvXv3Nnz76aef3nL4uw179tlnFy1aFB4eXlZWtmvXrpdeeknrirSxbNmytWvXGr4+ceKEs7Pzb7/9ZrSzq+Zy9OjRurq6xMTE9957z2wXtUyHDx82fKHX63v06PHRRx9pW48lUBRl2LBhb731ltaFaCwxMfGZZ54JCAjQuhAtffbZZ2PHjtW6CouwYsUKf3//CxcuaF2IZQkPD587d64RT2i+TrmgoCB0fBl0797d8AXn3NfXt7a2Vtt6LIQkSW3bttW6Ci0tWLDAy8srLi5O60K0V15evnHjxkOHDimK0TZftUaLFi2aOXNmfX39li1bSktLtS7HIuTl5W3fvn3y5MlGPKf9Pp2yBNu3b9+5c+e4ceO0LkRL6enpY8aMufvuu/v27Ttt2jSty9FMUVHRm2+++c4772hdiEUoKCiYO3duQkLCkCFDrn4yZG9OnDixefPmpKSkt99+OygoKDs7W+uKtPfFF1+MGDHCuM9QEISaOX78+Pjx4z/55JNOnTppXYuWevfuPWPGjKlTpy5fvvyXX37RuhzNPPbYY3PmzGnfvr3WhWhv2rRphw8fXrt27YkTJ1q1avWvf/1L64o0U11dXV1dvWvXrrS0tDfeeOOxxx7TuiKNSZL0zTffPPjgg8Y9Ldb81Mbp06djY2NfeeWVSZMmaV2LxgICAgICAhITEysrK99///2hQ4dqXZEGjhw5kpWV5enpuWnTpjNnzly6dOnhhx/+z3/+4+HhoXVpGmhYgNfR0XHChAlffvmltvVoyMfHJywsTBRFIoqJiXn44YdramrsedHRrKwsWZYTExONe1oEoQbOnDkTExPz7LPPzpgxQ+taLEhpaWmz1tK1JV5eXu+//77ha2dn5x07doSEhOCZOhHt37/fz89P6yo0ExkZeeTIEcPXeXl5Hh4e9pyCRPTVV1898MADRt+rwHyLbi9dunT37t1Lly4NDAwMDQ2dOHGiEWeBWJfg4ODy8vKxY8cavh00aNDo0aO1LUkrcXFxgwcP9vT03Ldv348//piTkxMSEqJ1URpbu3btzJkz7Xn6xPTp0zt06ODj47N79+6ffvopNzfXbn8rzp8/HxwcfP/993fu3Pntt99+/PHHn3nmGa2L0syFCxc6duy4b9++hvGGxmK+IMzNzT169GjDt9HR0XY7N+i7776rrq5u+Pbuu++2z/5AIsrOzt6+fXtpaWlAQMD48eN9fX21rkh7Z86c2bx587333qt1IZrZunXrxo0bS0tLO3bsOH78eHtuERJRYWHhN998U15ePnz48JiYGK3L0VJeXt6OHTtM8b8GtmECAAC7hlGjAABg1xCEAABg1xCEAABg1/4/ObScI1GmplMAAAAASUVORK5CYII=", "text/html": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ], "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ] }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "p" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "The `info` object passed to the callback contains not just the densities\n", "but also the complete Bloch wave (in `ψ`), the `occupation`, band `eigenvalues`\n", "and so on.\n", "See [`src/scf/self_consistent_field.jl`](https://dftk.org/blob/master/src/scf/self_consistent_field.jl#L101)\n", "for all currently available keys.\n", "\n", "!!! tip \"Debugging with callbacks\"\n", " Very handy for debugging SCF algorithms is to employ callbacks\n", " with an `@infiltrate` from [Infiltrator.jl](https://github.com/JuliaDebug/Infiltrator.jl)\n", " to interactively monitor what is happening each SCF step." ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.5.3" }, "kernelspec": { "name": "julia-1.5", "display_name": "Julia 1.5.3", "language": "julia" } }, "nbformat": 4 }