{ "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 = [ElementPsp(el.symbol, psp=load_psp(\"hgh/lda/si-q4.hgh\"))\n", " for el in load_atoms(silicon)]\n", "positions = load_positions(silicon)\n", "lattice = load_lattice(silicon);\n", "\n", "model = model_LDA(lattice, atoms, positions)\n", "basis = PlaneWaveBasis(model; Ecut=5, kgrid=[3, 3, 3]);" ], "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 - info.ρin))\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 log10(ΔE) log10(Δρ) α Diag\n", "--- --------------- --------- --------- ---- ----\n", " 1 -7.846663040835 -0.70 0.80 4.2\n", " 2 -7.852276630944 -2.25 -1.54 0.80 1.0\n", " 3 -7.852621336747 -3.46 -2.53 0.80 3.2\n", " 4 -7.852621949560 -6.21 -3.35 0.80 2.0\n", " 5 -7.852621959209 -8.02 -4.64 0.80 1.4\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+gvaeTAAAgAElEQVR4nO3deUATZ94H8GeOJMihgIgIKKIVRTksl6CgUjnK6QW0bltta6s97Npr1X3b+lpXd7Xbbl/t4bFb27W7tlVBxVtQARFUKodaRUTF+0IROSJkjvePWEoVNZCQSTLfz1/JNJn5EWm+PDO/5xlKFEUCAAAgV7TUBQAAAEgJQQgAALKGIAQAAFlDEAIAgKwhCAEAQNYQhAAAIGsIQgAAkDUEIQAAyBqCEAAAZM1IQSiKYnl5+dmzZ41zOAAAAB1Rxlli7cUXX2xoaKivrw8MDFywYEGbr2lubj558qSvr6+O++R5nmEYw9UI7SMIAk3jjII08OFLC18+EhIEgaIoiqIMuE9j/L9UWlpaWVm5bt26LVu2rF+//sqVK22+7PLly8nJybrvtrGx0UAFQkfg85dQY2MjVgmWEH75JdTU1CQIgmH3aYwg/Pnnn0eMGEEIYRgmKCiorKzMCAcFAADQhTGCsLa21sbGRvvYzs6upqbGCAcFAADQRScGYXNzs/ZBjx49bt68qX1cXV3ds2fPzjsoAABAu+gVhHl5eSkpKYMGDXrmmWdabz98+LC3t3ePHj369euXn58fGRmZnZ3d1NR0+/bt4uLikJAQ/WoGAAAwGL2CkOf52NjYsWPHXr58uWWjKIrPPffcG2+8UVtb+5e//OXZZ5/t1avXG2+8MXLkyOjo6EWLFtna2nb4iHUacpdvY/uNux3eJQAAyJoBpk/885//XL169b59+7RPDxw4EB8ff/36dZZlCSEeHh7Lly+Pi4t77H6qqqoCAwNHjRrVsmXp0qX29vatX7P1Ev3taeY/IzRWDKmvr9dm6nenmdzr9LdhGj1/EGiXhoaGlku/YGQNDQ3W1taG7SAH3bV8+YDxqdVqpVKp+/QVa2vrx841YvWu6n6VlZVeXl7aFCSEDB48+NSpU7oEISFEpVI9++yzLU8dHR2trKxavyB1AGkg4vMFdPoYpksX3traelWFuO2KsH4MY8UoDPhTwGMJgmBtbS11FTKl/fARhFLheR6//FKhKKpdQajLjFvDB+Ht27dbDxTa1SaqUqnS0tIe/ZqXBxKKolL2CP8eRq8/RTacEzOiFVaY22p0NE1jTrdUtB8+glAq+OWXEP0rA+7T8EHo5OR0586dlqc1NTXOzs6GPcRLXrRISNQelaOVsPNpVoUUBACAjjL8HzXe3t4nTpxQq9WEEEEQjhw54u3tbfCj8AK5y5HCa+KsIv52s8F3DwAAcqHXiLC6urq0tPTEiRO3b9/Ozs52dnb28/Pz9/f38fGZO3furFmzVqxY4ejo2Lr/xSC+OSlsPCfsj2368bLNZ0f4HyqF93yZd3xpJc5VAICxZGdnr1y5UuoqZCQ2Nnbq1KmdsWe9grCysnLx4sWEEBcXl8WLF0dERPj5+RFCfvrppxkzZgQFBQ0aNGjjxo2GvZLxzUkho0pIj2I1jeStwbQtS749KeReEVZVCAuC6FRPhCEAGMOhQ4c4jps0aZLUhchCQUFBVlaWKQZhaGhoVlbWg9s9PDw2b96sz54fZuM5YeM5ISOKVTFEO1tCe71w10Xxq+H0uwf5FSeEz0IZf0c0EQBAp/P29k5NTZW6ClkQBGHDhg2dtHPDN8t0qihXOs6dvq875mUveqwH6a4ixePYVRXC09u5aDf6k2GMSxeJqgQAAPNhZicSbRWkzR7R7ipCCGFpMm0QfSpN0a8r8VmvmVfMt7kMDQAAQAszC0Jd2CrIvACmeDx75g7xWsutPiXgvm0AAPAwFhiEWn1sqdWjmbVjmOUnhNBN3P5rSEMAMAmP/dtcJARfWMZksUGoFepM7U9m3/ejn8/h03bzVXX47QIAKX1+TEjZzWsefot1kZC3Cvh3DuC6jvFYeBASQihCUj3pXyaygU7UsExuThF/B6tzA4BEnnKl9l0VxmdxTW0lnUjIzEJ+VYWQ2Nvyv5xNh1w+a2uWzPanS8azNU1kyHpuZbnAY3AIAEbn70hlx7MHb4gTs+/PQm0K/uukkBnNRrl1ZA6Yh4dHy13Q9SQIQnZ2Ns/rOzBtampycnJqamp6xGv8/f1Pnz6t54H0IZcg1HK1plaEMxujmf9WCn7p3I6LCEMAMLY2s1D/FCSENDU16X9nPa3m5ubo6OhHB5guWJadO3duy/2IHnYsQXj4yeLOZ2bzCA0i0InKTWQ3nxdmFPCeduT/QpkhDpiADwDGo83CqG3cxGwuPYpVMgZIwfuIovjTTz+Vlpb2799/8uTJKpVKuz07Ozs/P9/Z2XnSpEkODg6EkL179/bo0cPHx4cQcvHixaKiovHjx2/cuJEQsmrVKqVSmZyc7OLiostBGxsbV6xYUV1dHRYWlpiYqN3Ycje9jRs3BgQE5OTknDp1asyYMaNHjzbIT6o/eY0IW0vqQ59IYRN706O3ctPzedzjHgCMqfW4cEaBgVOQEPLWW2+tXLnyiSee2LRpU3JysnbjggULZsyY4ezsXFZWNnTo0Fu3bhFCvv322z179mhfcPz48b/97W8dO6JarY6IiCgvL+/Tp8+sWbMWLVpECNFoNNOnT9doNISQxYsXP/vssxUVFd26dZs4cWJOTo7+P6ZByHFE2EJBk5k+9AsD6I+L+SHrNVi5GwD0lLCTa27PSb7+XakdF0WRiP6O1OIj/OIj7XjvujGsvbLt/3Tu3Lk1a9acO3fOzs5uypQpnp6e+fn5/v7+CxcuLC4u1t4RKCkp6euvv/7www/b3MO4ceMIIS+//LLutyBeu3atg4PDihUrCCHh4eGhoaHvvPPOfa+Ji4v76KOPCCG3bt3auHGjiQwKZR2EWo4qsiSMeWMw/f5BHit3A4A+3vNldF/CQyTky+MCTYsKitiw1Hu+NNueAaH1w7+/jx496u3tbWdnRwhRKBSBgYHHjh2ztbW1sbFpuS9eREREcXFxO45HCCHk+vXr7733HiGEZdlvv/229X86cuRIaGio9vGQIUMoiqqqqvLw8Gj9mqFDh2ofuLm5FRQUtPfonQRBeM/AbtTmGDb7kvjeQX75CeEfWLkbANrvKVddvze03TFZl4RtMWyPLiRqG/flL2J6lGHuNN7Q0HD37m/Xe9RqtbW1tbW19d27d0VR1N4RqLGxUTvaYximpTtUeyvZR7C1tZ04cSIh5MF7xDc0NLQ8FgShubn5wdEkw5jijdQx9PmdKDeqZDw7ZQAdt4ObnMNffcyvBABAR9zXI/qIORUdVlZWdvz4cULIhQsXDhw4MGLECE9PT0dHx4yMDEKIWq1et26d9sxk3759tUNDURTT09O1b1epVF26dNFeRGzN2tp63Lhx48aNa7nu2FpmZmZjYyMhJD09vXfv3m5ubob5YToZgvB+NEUmD6ArUn9buVvNSV0TAFiQNmdKGDwL/fz8XnrppeTk5ODg4A8++KB///4KheLbb7/94x//GBcX5+Pj8+STTz733HOEkKlTp+bk5IwaNSooKKhl9gVFUTNmzAgMDAwKCiotLdXxoAMGDBg+fHhiYuKbb765YsWKB0eNpgmnRtumXbl76kD6gyJh4DpuQRD9wgAap0oBQE+PmC9435wKPc+ROjk5bdmypbS01N3d3dXVVbtxzJgxlZWVJ06c6NGjR+/evbUb3d3dy8vLjx8/7u7u7uTk1HJO9ZNPPlmwYEFDQ4P2WqMuwsPDX3/99YqKCm9v765duxJCrKysrl27pp1BsWvXrpapFK+88sqLL76o109oOAjCR+ltQ60ezRy8Lr57kP/yuPCPYUy4C9IQADpufrHwiJkS2iwcs5V7fT+/aqS+l9OUSmVISMh9G7t06RIQEHDfRhsbm+DgYO1jhULReg9K5UM6Ux/CwcFh2LBhrbc4OztrH7QOVJVK1TK1UXIIwscb5kzlJ7HrzwqTc/kgJ2pxCO1phzgEgI6YMoAa2YuN7PXQ7xB/R2p3Ast3aKGV3r17a7tRvL29n3nmmQ4X2TGxsbH29vYdeKObm1vr9DU+BKFOtCt3J/Smv/hFCN7IPfcE/ZcgpquU/3AAYJb62lF9H3eiscMt60VFRdoHfn5+fn5+HdtJh40fP75jb8zOzjZsJe1lHlcyTYR25e4jE9m7PPFaq1lyDCt3AwCYPQRhu2lX7t4Wy244J/ilc9svIAwBAMwYgrCDApyonAR2UQj9x0I+ejv3Sw3iEADALCEI9ZLUhz6ewib2piOxcjcAgHlCs4y+Wlbu/uQI75uueceHeduHNsgiSQBgyjw9Pd9///2dO3dKXYgs1NTUtCxkanAIQsNwVJFFwcxLXvRHPwt+GRxW7gaweJMmTRo4cKChboQLj+Xu7t5Je0YQGtLAbtTaMczuy+K7B/jlJ4TPhjFDu2PGIYDFenBmOpgjjFoMb4zrvZW743di5W4AAFOHIOwUWLkbAMBcIAg7kXbl7pIJ7Jk7ZOA6bvUp3W/YCQAARoIg7HTalbvXRzEryoVhm7j8q0hDAAATgiA0kpAeVH4S+yc/+oVcPmkXd7YOcQgAYBIQhMajXbn7+EQ2vCcdvJGbWcjf0UhdEwCA7CEIja1Lq5W7B2DlbgAAqSEIpaFduXt7LLvxnOCbzm3Dyt0AABJBEEopwInam8AuDWP+dJCP3s4dw8rdAABGhyCUXpQbVTqBTfWkY7Zj5W4AAGNDEJoEBU2mDaKPTVQ4qIhvumZxmdDES10TAIA8IAhNiHbl7rxE9nC16JvBrTsrSF0RAIDlQxCaHK9u1NoxzLIRzMIS4amtXMlNXDgEAOhECEITNcaVKh7PvuhFJ+zkJufwVxqlLggAwEIhCE2XduXuU2mKfl2Jb7pmThFfjwn4AACGhiA0dTbsvZW7LzcQ7/XcynIBS3cDABgQgtA8aFfuTo9iVp8Shm3i9mHlbgAAA0EQmpOQHtS+JHaWPz05l0/axZ3Byt0AAHpDEJqZ1it3D9vEzSzka5ulrgkAwJwhCM2SduXuEykKQojXOs2SYwKHOYcAAB2CIDRjTlZkSRizN4HdeUnwzeC2YuVuAID2QxCavcH21LZY9oswZvYhrNwNANBuCEILEeVGlYz/beXu62qpCwIAMBMIQsuhXbn7l4kKBxXxy8DK3QAAOkEQWhoHrNwNANAeCELLpF25e/kIZmGJEJbJHbiOC4cAAG1DEFqyp1yp4vHs6970xGw+bTd/vh5xCABwPwShhdOu3F2Rxg52IAEbOKzcDQBwHwShLGhX7i7Fyt0AAA9AEMqIuw21ejSTEcV8XymEbOLysHI3AACCUIaCe1B5iexsf/pFrNwNAIAglCftyt3lqWyUK1buBgC5QxDKl5ImM32wcjcAyB2CUO60K3fnYOVuAJArBCEQQoj371fuPnoLcQgAcoEghN+0rNwduwMrdwOAXCAI4XceXLn7LlbuBgCLhiCENmhX7t6XyB6uFgeu41afujf/vq6tVWnUHEGXDQCYLwQhPNSAbtTaMcx3o5jPjwnDM7mC6+L4LC7r0u8uHzZyZGwWt/8arikCgLlCEMJjRPaiDo9jX/emU7N5OwX1P0V8S2dpI0fGZ3EvDKBH9aKkLRIAoMMQhPB4LSt3+3cnZ+rE1/K59WdFbQo+P4B+4Qn8FgGAGcNXGOhKu3J3yXg2xJl+bg8XukPxAlIQAMwfvsWgffrYUt+PYny6U5fV1NcnxKuYYgEAZg5BCO2jPSP6tg99MllzuUH0Xq9ZWY6eUQAwYwhCaIfW1wW7q8TSCWxfW+qTI0LiTu5yIxpHAcAsGSkI58+fHxQUFB8fb5zDQWcQCUnZzbW+LmivJNnxrJMV8bCjntzAYWgIAObISEGYkJDwzTffXLx40TiHg85AEbIinHn+990x3VVkWyz79xAmK45ddkJI2sVdaZSqQACAjjBSEAYGBjo4OBjnWNB5etu0MV/QUUWsWeLnSB0ay4b3pIM2cuvPYmgIAGaDNezucnNzm5qaWm8ZOXKklZWVYY8CpklBk9n+dKQr9WIuv/aM+PUIxgn/8gBg8gwchPn5+fX19a23DBs2DEEoKyE9qOLx7Lxi3i9D8+VwZkJfNGQBgEnTNQi3b99+4MCBkydPvv3226GhoS3bd+7cuWTJErVaPWnSpGnTpn3wwQedUyeYEyuGLApmxnnQL+XxP57G0BAATJquf61/8803dXV1ubm5Fy5caNl47Nix1NTUKVOmzJs3b+HChT/88MPD3v6Pf/zjjTfeOH/+fFpaWk5Ojp5Fg1kIdaZKxrP9uhL/DG5DFa4aAoCJokSxHdO/fH19586dm5qaqn06Y8YMnueXLVtGCPnXv/61atWqgoKCNt947dq1hoYG7eMePXrY2dk9+JqqqqrQ0NA//vGPLVtef/11W1vbhxVTV1fX5n7AOOrr6x/xr9PagRvklXzi60C+DCPdVZ1dlyzU19fb2NhQFNY6lwa+fCSkVquVSiXDMDq+nmXZx/6fotc1wtLS0ldffVX7+L4Mu0/Pnj112aEgCLdu3Wp5ynEczz/0trA8zz/iv0Jn0/3zD3YkB+LJwqNUwCZ6SQif1LuzS7N82g8fQSgVfPlIqL0fPsMwnRuE169ft7e31z52dHRUq9W1tbXdunXr8A5tbGw+/fRTHV+s0WjQhiMhjuN0//ytCPk0jEzoJ76UR607Ty0PZxwxNNSD9sNHEEoFXz4SEkWxXSNCXejV0WdnZ9fYeG/6dENDA8MwOp4rA3ka3pMqHc/260r8MrjMc7hqCAAmQa8g9PDwOH36tPZxZWWlu7u7YVMaLE8XliwKZn6MZN47KKTt5muaHv8WAIBOpVcQTpo06fvvv6+rqxNFcfny5ZMmTTJUWWDZwl2osglsv67EN4Pbch6rdQOAlHQNwujoaIqijh07lpaWRlFUXl4eIWTChAlhYWH9+/f39PSsqamZPXt2Z5YKFsWaJYuCmR8imXcO8NPz+TqN1AUBgFzp2iyTlZX14EaGYb777rvq6uqmpiY3NzeDFgayEOFCHR7P/ukg75fB/SuCGeOK7g8AMDYDLLHm5OSk/05AtroqyIpwZtclcWoeH+tOfTaMsVVIXRMAyAnWgQSTEONGHZnIEkL8Mri9V3DVEACMB0EIpkI7NFw2gnkxl5+ez9fjqiEAGAWCEExLrDt1ZAJLCPHP4HIwNASAzocgBJPTTUlWhDNfDmcm5/DT8/kGTuqCAMCiIQjBRMX1vnfV0D+Dy8XQEAA6DYIQTJe9kqwIZ5aGMc/n8NPz+UYMDQGgEyAIwdTF96aO/tpQuu8qhoYAYGAIQjAD2qHh/4Uyk/byMwv5JtwABwAMB0EIZiOxD3V0AlvTRAI2cEU3MDQEAMNAEII5cVCR1aOZeYF08i5uThGGhgBgAAhCMD+pnnTZBEXlHRK4kfu5GkNDANALghDMknMXsn4M878BdNJObk4R34y7/AJARyEIwYyletKlExQVtSRwA3cYQ0MA6BAEIZi3nl1IRhQzN4BOwNAQADoEQQiWQDs0PHmbBG/kSm5iaAgA7YAgBAvh0oVsiGY+fJKO28HNKeI1GBoCgG4QhGBRtEPDE7dJ8EauFENDANABghAsjUsXsima+eBJOhZDQwDQAYIQLFOqJ106XvFLjRi+mTtxG0NDAHgoBCFYrF7WZHMM++ZgetQWbnGZwCMNAaAtCEKwcJMH0EXj2OzLQvhmrhxDQwB4AIIQLJ+HLbUrjn3Jix6JoSEAPABBCLJAETJtEH1oLLvrkhCxmTtZizAEgHsQhCAjfe2o7Hj2RS86YjO3uEwQkIYAgCAEudEODQ+OZXdeFCK2cBUYGgLIHoIQ5MjTjtqdwE4ZQIdjaAggewhCkCnt0PDAWHb7RWHkFq7yDsIQQKYQhCBr/eyovQns5AH08EwMDQFkCkEIcqcdGhYms9suCKO3cqcxNASQGQQhACGE9O9K7U1gn3+CDt/MrSzHyBBARhCEAPfQFJk2iN6dwH5zUnh6O3ehAWkIIAsIQoDfGWxP7U9in3KlAzdgaAggCwhCgPuxNJntT+9NYP9ZLsTt4C5iaAhg0RCEAG0b4kAVJrORvejAjdzKctzVEMBiIQgBHko7NNwdz64sF+J3cpcwNASwRAhCgMfwcaAOJLOjXOgADA0BLBGCEODxtEPD7Dh2RbmQsJO73IihIYDlQBAC6MrXkTqQzI50oZ/cgKEhgOVAEAK0g4Ims/3prDh2+QkhaReGhgCWAEEI0G5+jtTBsWx4T/rJDdz3lRgaApg3BCFAR2iHhrvi2M+OCGm7+Rt3pS4IADoKQQjQcf6O1IGxbL+uxD9Dk34WQ0MAs4QgBNCLFUMWBTMZUeyHh4W03Xw1hoYA5gZBCGAAoc5UyXjt0JDbUIWhIYA5QRACGIZ2aLg+ivlzkZC2m7/ZJHVBAKAbBCGAIYX9OjT0S+c2nsPQEMAMIAgBDKwLSxYFM+vGMLMPCWm7+VsYGgKYNgQhQKcY3pMq1Q4NM7hMDA0BTBiCEKCzaIeGPz3FvH9ISNvN12BoCGCSEIQAnWvEr0ND3wxuy3ksyQZgchCEAJ3OmiWLgpkfIpl3DvCTc/g6jdQFAUArCEIAI4lwocomsA4q4pfB7b6MoSGAqUAQAhiPNUuWhDErwpmpefz0fAwNAUwCghDA2GLcqCMTWUKIfwa3B0NDAKkhCAEk0FVBVoQzy0YwL+Xx0/P5egwNAaSDIASQTKw7dWTCvaHh3isYGgJIA0EIIKVuSrIinPlqBDMlh5+ezzdwUhcEID8IQgDpPe1+76qhXzqXi6EhgHEhCAFMgr2SrAhnvhjOPI+hIYBxIQgBTEh8b+rorw2leVcxNAQwBgQhgGnRDg2XhDHP7eWn5/ONGBoCdDIEIYApSuh9r6E0eCN36AaGhgCdCEEIYKIcVGRFODMvkB67i5tTxJfcFB+cYsGLZNkJATkJoA8EIYBJS/WkyyYoKu+QSXv5Pxby2Zd+Sz1eJC/m8reaCCVhfQDmD0EIYOqcu5D1Y5i/BNLXG8Upedz2CyL5NQUH2VMfDMX/xQB6YaUuAAB0kupJj3Shp+7jJ+7m/uLPlNYKSEEAg0AQApiNnl3Ilhjmq+PUzEIyxFFcGsZIXRGAJcCfkwDmhBfJgeviNC/haqPY7yfN58eEZkHqmgDMHIIQwGy0XBf85EnN0Qmshx21sUocsJZbfQqNowAdhyAEMA/3dcf0sCJZcWyTIL41mF7yixC6iduHlWgAOsQYQSgIwt///vexY8empKRs377dCEcEsDwF10Rfx991x/SwIptj2BO3xaJx7Pt+9JRcPmkXV3kHcQjQPsZoluE4TqlULlmypK6ubuzYsTt27PDy8jLCcQEsSYQLFeFy/4zBHlbkm5EMISTVkx7rQS87LoRlchP60guCmB5WUlQJYIaMMSJUKpUzZ87s27evr6+vn59fVVWVEQ4KIDdKmsz0octTFQ4q4puuWVwm3OWlrgnAHBg4CDUaTfPvieJvJ2pKS0srKysjIiIMe1AAaNFdRRYFM3mJ7OFqceA69NEAPJ6BT41OmTLl7Nmzrbf8+9//1p4IraysnDx58o8//tilSxfDHhQA7uPVjVo7hjlwXXzvIP/lceHTYczIB06rAoCWrkFYV1dXUlJSUVExYsQIb2/vlu3V1dWrV6+ura1NSkoKCgpas2ZNm2+vqqpKSUn57rvvfHx8DFA1AOgg1JnKT2LXnxVezOWHOJDPQ5knuiIOAe6n66nRiIiIGTNmzJkzJzc3t2VjXV1dSEhIaWmpSqWKiYnZuXNnm+9Vq9VPPfVUcHDwzz//vHLlyoqKCgMUDgA6oAhJ9aTLU9koVzosk5uez19XS10TgInRdURYVFSkUCgiIyNbb/z+++/d3NxWr15NCLG3t1+4cGFsbOyD76Uoas6cObocheO4w4cPtzz18/NTKBQ6VggAD6Pto3lhAP3JEX7wes2MIfQcf8YKC7QBEEJ0D8I2A2nv3r0tyff000/PmDGjublZqVTe9zIrK6tp06bpcpS6urpXX3215emGDRu6d+/+sBc3NDRQFM7zSKahoUHqEuSroaFBFMX2/v4rCfnQmzzrTs0/wg74SZg1hJvSj6fx/1D74ctHQmq1WqlUMoyuf8dZW1vT9GPOferVLHP16tXo6GjtYxcXF1EUr1692qdPnw7v0MHBobi4WMcXi6Joa2vb4WOB/vD5S8jGxqZj38VDbUlGL3Lguvj+QWr1WcWnw5hRvfCd3j748pEQwzDtCkJd6DV9gmEYQbi34i/HcYQQlsXtLADMQKgztS+JneVPv5zHR2/nfqnBJAuQL72CsFevXleuXNE+vnz5MsMwzs7OhqgKADqdto/mRCqb2JuO3MpNz+evoY8GZEmvIExISMjMzNSOBTds2PD0009jRAhgXlqvRzNkvWZeMY/1aEBudA3ChQsXRkdHl5WVLVmyJDo6ev/+/YSQ1NRUpVIZFRX1yiuvfPrpp3Pnzu3MUgGgsziqyKJgpjCZPV5DvNZyK8sFLEgD8kG1XgLtEY4ePXrt2rWWp35+ftqzoE1NTdu3b6+trY2KinJzc9OnlKqqqsjIyPsWpnmEuro6Ozs7fY4I+qivr0e/gFTq6+s73CzzWAevi+8f4tUc+XQYMxp9NG3Bl4+E2ts1qgtdz2T6+vr6+vo+uF2lUo0bN86ABQGAtIY5U3mJ7Pqzwiv7eE878nko4+OAOARLhhvzAsD97vXRpLCpnnTMdvTRgIVDEAJA2xQ0mTaIPjbxtz4aNSd1TQCdAEEIAI+i7aMpHs+euUO81nErywUefTRgWRCEAPB4fWyp1aOZ9CjmP5XCsE3c3isIQ7AcCEIA0FVIDyovkf3fAPrVfXz0du7oLcQhWAIEIQC0T1Kfe300sTu46fn8VfTRgJlDEAJAu2n7aH6ZqHBQEZ/1mjlFfL1G6poAOgpBCAAd5PBrH83lBuK9Hn00YK4QhACgF20fTUYU859KwT+D23YBYQhmBkEIAAYQ3IPKSzwpx9MAABebSURBVGT/FkzPLOSjt3NH0EcD5gNBCAAGk9SHPq7to9nOTc5BHw2YBwQhABiSto/meIrC1eZeH00d+mjAtCEIAcDwtH00JRPYyw1kMPpowLQhCAGgs/S2uddH899KwS+d24o+GjBJCEIA6FzBPajcRHZRCP3OAfTRgClCEAKAMST1oX+Z+FsfzZVGqQsC+BWCEACMRNtHcypN0a8r8U1HHw2YCgQhABiVrYLMC2BKJrA1TeijAZOAIAQACfS2oVaEMxuimTWnBd90bst5hCFIBkEIAJIJcqJyEtilYcycIj56O1eGPhqQAoIQACQW5UaVjGdTPemn0UcDUkAQAoD0HuyjuYM+GjAWBCEAmAptH03pBLamiXit1Sw5hj4aMAYEIQCYFncbakU4szWW3XhO8E3n1p0VpK4ILByCEABMUaATtTeBXRrGzC8WordzpTcxNoTOgiAEANPV0kcTv5NL282fq0ccguEhCAHApLE0mTaIrkhVDHYggRs49NGAwSEIAcAMPNhHw+HSIRgIghAAzIa2j2ZbLLvpnOCbgT4aMAwEIQCYmQAnak8C+0UY85cSIWobV4I+GtAPghAAzFKUG1U8jk3rRyfs5NJ281V1iEPoIAQhAJgr9tf1aAKdqGGZ6KOBDkIQAoB5s2HJbH+6eBxb00QGoI8G2g9BCACWwM2GWhHO7Ilnd10SfLAeDbQHghAALMcQB2prLPvlcGZBiTBmG1dcjQuH8HgIQgCwNNr1aKYMoJN28eijgcdCEAKABaIpMnkAXZHGtvTR1DZLXROYKgQhAFis1n00XuvQRwNtQxACgIXT9tHsTWCz0EcDbUEQAoAsDLantvzaRzM8kyu8jguHcA+CEABkRNtH85o3nZLNp+3mz6KPBhCEACA3rftogjdyMwvRRyN3CEIAkCNtH82Riexdnnit0ywuE5px6VCuEIQAIF+u1tSKcCYngc2/Jviij0auEIQAIHfe9tTmGPar4czCEiEskyu4hguH8oIgBAAgRHtfp/Hs69502h700cgLghAA4J57fTSpv/XR3EYfjQwgCAEAfseaJbP96fJUBSFk4DrN4jKhiSdvH2hjjPifSgGXFS0AghAAoA1OVmRJGJOTwB6uFv0yOHcbMjGbP33ntyz86YzwbYUQ545vUbOHf0IAgIfytqfWjmG+Gs78t1LkRRK/814W/nRGWFkubIpmbRVSlwh6QxACADxGlBt1eBw704eubSbBm7i/HWORgpYEQQgA8Hg0RV72ok8/wwZ0pxYfZ57pRyMFLQaCEABAV1vOCyIhS4P4GQX8/BJe6nLAMFipCwAAMA8t1wXFu+q+jqqJ2fyZO+KqkSxNSV0Z6AcjQgCAx7uvOyahN70rjk0/Kybs5NSc1MWBfhCEAACPV6chm2N+1x0z0oXancBqRGrMNq76rnSVgd4QhAAAj/fKQNr6gUtJIT2orDgmxp0Ky+RO1WJJNnOFIAQA6DiKkHkBzGx/euQWbj9W6zZPCEIAAH29MpD+bhQ7PotbewYrrpkfBCEAgAHEulPZ8eyfDglLf0EWmhkEIQCAYfg5UgVJzKqTwsxCXsBZUvOBIAQAMBg3G2pfEltRK07M5hsxrcJMIAgBAAzJTkE2xbB2ChK5lbuulroa0AGCEADAwJQ0+fdoJq43FZbJld/GSVJThyAEADA87bSKj56kR2/lcq8gC00aghAAoLO86EX/N5JN28OtOY1WUtOFIAQA6ERjXKnd8ez/FAnzinG3ChOFIAQA6Fw+DlRhMrv5vDg1j9dgZGh6EIQAAJ2ulzXJTWCvqcWEndwdjdTVwO8hCAEAjMFWQTbFsE90pcI3cxca0D5jQoxxY15BEDZt2nTq1ClbW9uUlBRnZ2cjHBQAwNQwFPl6BLPkmDA8k98cwwztjlv6mgRjjAg5jisrK/Pw8Lh7925oaGhDQ4MRDgoAYJpm+tCfh9Ix27ntFzAuNAnGGBEqlcp58+ZpH//4448XLlwYNGiQEY4LAGCaUjxpNxtqQhb3cSAzbRAuUUnMGEFICOE47u9///vJkydDQ0ORggAAYc7UviQ2fid/pk78WzCDk6QSMmQQNjc3jx49+r6NBQUFhBCapgMDA21tbVetWnXr1i1HR0cDHhcAwBw90ZUqTGbH7uJezOX/GcEoMTKUSDuCsKGh4dq1a25ubiqVqmWjRqM5ePAgTdMhISFKpXLv3r1tvpem6ZiYmJiYmOLi4j179qSkpOhbOACA+euuItnx7JRcPm4Hlx7F2iulLkiWdApCjuOefPLJiooKjUZTWFg4bNgw7fabN2+OGjXKzs6O53me5/fu3du1a9cH337p0qUrV64MGjTozJkz+fn5f/7znw35EwAAmDMrhqyJZN49wIdv5rbGMh62OEtqbDoNxWmaXrZs2c2bN+3s7Fpv/+KLLzw9PQsKCg4ePNi9e/fly5e3+XZBEJYuXZqYmLhgwYJly5Z5eXkZoHAAAEvBUGRJGPPqQHrEZv5wNVpJjU2nESFN0+Hh4Q9uz8jI+PDDDymKIoQ899xzX3/99axZsx58We/evVevXq3LgW7evOng4NDy9PDhw4+YdNjQ0KA9NEgC02Ak1NDQIIoifv+l0klfPlP7kp4K+unt7LJh/NOuWJi0bWq1WqlUMgyj4+utra1p+jFDPr2aZS5cuODh4aF97OHhcfHiRX32RghxdHQsKSlpedqtW7dH/ACiKNra2up5RNAHPn8J2djYIAil0nlfPs8OJP0cxXFZ1AdDVW8ORvNMGxiGaVcQ6kKvIGxqalIoFNrHKpVKrdb3ZswURbUeEQIAyE1ID2p/Ehu/k6+oFT8PZWj8tdP59PqLw8XF5ebNm9rHN27c6NWrlyFKAgCQNU87an8SW3pTTNvDqzmpq5EBvYIwNDQ0JydH+zg3NzcsLMwAFQEAyJ6jiuyKY1U0eWobd+Ou1NVYOl1Pja5cubKmpqapqen777/Pycl57bXXunXr9vbbb0dHR7u7u/M8/69//Wvfvn2dWisAgHyoGPKfSObjYj4sk9sWy3h1w0nSzqLriPDOnTs1NTVvv/22ra1tTU2NIAiEkJCQkG3btpWUlBw7diwrK8vPz68zSwUAkBeKkHkBzP8MpUdu4fKvYlpFZ6FE0VQ+3KqqqsjIyLNnz+r4+rq6uvvmNYIx1dfXo2tUKvX19egalZDxv3yyLonP53BLw5hn+sm9lbS90yd0IffPFADA9EW7Udlx7KxDwrxizC80PAQhAIAZ8HWkCpOZTefEafk8J0hdjWVBEAIAmAdXayovkb3YIKbs5hsxrcJwEIQAAGbDTkEyo9meXcjordw1fZcwgXsQhAAA5oSlyYpwZmJfOiyTO3HbVLodzRqCEADA/Mz2p+cF0JFbuZwryEJ9IQgBAMzS5AH0mkj2mT3cfyrRPKMXBCEAgLl6ypXaE89+dBjTKvSCIAQAMGNDHKiCJHbLefHlPF6DkWGHIAgBAMxbL2uSm8jeuCvG7+Rqm6WuxgwhCAEAzJ4NSzZGswO7URFbuAsNaJ9pHwQhAIAlYCjy5XBmqhc9PJMvuYksbAcEIQCA5ZjpQ/9fKB27ndt2AVmoKwQhAIBFmehJZ8awr+zjlp9A84xOEIQAAJYm1Jnal8h+fkyYWchjYPhYCEIAAAvUvytVkMwW3xSf2c3fxSTDR0IQAgBYpu4qkh3H0hSJ2sZV35W6GhOGIAQAsFgqhvzwFBPlRo3awlXV4Sxp2xCEAACWjCJkXgDztg8dsYX/uRpZ2AYEIQCA5Xt1EP3lcDp+B7f5PFpJ78dKXQAAABjDWA/a1Zoal8WfG0pmDMYo6Df4LAAA5CK4B5WfxHx9XJhZyAs4S/orBCEAgIx42lH7k9gjt8SU3byak7oa04AgBACQFwcV2RnH2rAkcht3XS11NSYAQQgAIDtKmqwezTztTg3fzJ2slftJUgQhAIAcaadVfDCUHrWF23dV1lmIIAQAkK+XvOj/jGZTdnM/nJbvtAoEIQCArEW5Udlx7JwiYV6xTNckRRACAMidryNVmMxknhNf3cdz8hsZIggBAIC4WlN5iezlRjFxF1enkboa40IQAgAAIYTYKsimaNbTjgrfzF1skFH7DIIQAADuYWmybATzshc9PJMvuyWXLEQQAgDA78z0oT8LpZ/ezu29IossRBACAMD9Uj3p9VHs83v57ystv3kGQQgAAG0Y0ZPKimfmHrb8aRUIQgAAaNtge6owmd16Xnwpj9dY7sgQQQgAAA/l0oXkJLI375K4HVxts9TVdA4EIQAAPIoNSzZEM972VPhm7ny9BbbPIAgBAOAxGIp8MZx5ZSA9fDNfXG1pWYggBAAAncz0oZeG0U/v4LZesKgsRBACAICuJvSlN8ewr+7jlp2wnOYZBCEAALTDMGcqP4ldckyYWcgLFjEyRBACAED79LOjCpLZkpviM3t4NSd1NXpDEAIAQLs5qkhWHKugyZhtXPVdqavRD4IQAAA6QsWQ/0YyMe5UWCZ3qtaMT5IiCAEAoIMoQuYFMLP96ZFbuP3XzDULEYQAAKCXVwbS341ix2dxa8+YZSspghAAAPQV605lx7N/OiQs/cX8shBBCAAABuDnSBUkMatOmt+0CgQhAAAYhpsNlZPIHr0lTszmG81nWgWCEAAADMZeSXbEsXYKErmVu66WuhrdIAgBAMCQlDT592gmrjcVlsmV3zaDk6QIQgAAMDDttIqPnqRHb+Xyrpp6FiIIAQCgU7zoRf83kk3dza05bdKtpAhCAADoLGNcqd3x7P8UCfOKealreSgEIQAAdCIfB6owmd18XnxlH68xyZEhghAAADpXL2uSm8BebRQTd3J3NFJX8wAEIQAAdDpbBdkUw/bvSkVs5i42mFb7DIIQAACMgaHI1yOYl73osEy+9KYJZSGCEAAAjGemD/15KB2zndtx0VSyEEEIAABGleJJb4xmX8rlVpabRPMMghAAAIxteE9qXxL76VFhThEv+cAQQQgAABJ4oitVmMzmXxVfyuWbJR0ZIggBAEAa3VUkO55V8yRuB3e7WbIyEIQAACAZK4asiWR8HKjwzdy5emnOkiIIAQBASgxFloQxrw6kR2zmD1dLkIUIQgAAkN5MH/rL4XTcDm7LeWNnIYIQAABMwjgPekssOy2f++q4UZtnEIQAAGAqQnpQ+5PYL48LMwt5wVgjQwQhAACYEE87an8SW3pTTNvDqzljHNGoQXjhwoW6ujpjHhEAAMyOo4rsimNVNHlqG3fjbqcfznhBmJ+fP3DgwIyMDKMdEQAAzJSKIf+JZGLdqbBMrqK2c0+SGikIm5qaFi5c+OyzzxrncAAAYO4oQuYFMH/2p4dv5l7Zd/8N7u/y5M0CvskQ9703UhB++OGHb775pq2trXEOBwAAlmHqQPq/keya00Lyrt8uGDYLJG037+tAqRgDHII1wD5+debMmW+++ab1Fg8Pj2nTph06dOjGjRuJiYm7du0y4OEAAEAOYt2owmQ2fDM3YjOXHUWaBfLMXj6+N/Wat2HGcgYIQo1Go1AoCCH29vYjR45s/Z8cHBwIIfPnz3dxcZkzZ05BQcGZM2d8fHwCAwP1Py4AAMiEvyN1IoV9cgMXso3xsBETPWhDpSDRMQirqqo+/vjj4uJitVpdUVHRsv369et/+MMfDh06pFQqFy9ePHXq1NjY2AffPn/+/Fu3bhFCTp486eXl1atXL0NVDwAAMuFuQx2dqBiSrvGwEQ2YgkTHINRoNP7+/uHh4a+99lrr7bNmzXJ1da2pqTl27FhERMTo0aP79+//4NsDAgK0D/bv3+/r6+vq6qp/3QAAICvNApm2j/9fP6G8jpl1iP8kxBCXBwkhhFCiqGtb6tGjRwMCAjQajfapWq3u3r17UVHRkCFDCCHPPPOMt7f3vHnzOlxKVVVVYGDgqFGjWrYsXbrU3t7+Ya+vr69H942EGhoabGxspK5CphoaGqytrSmKkroQmcKXj/E1C+T5fEWsq/AHt3qFUjmrVGXNiAuGPr5n1NramqYfM3zs+DXCS5cuNTU1DRo0SPvU29u7srKyw3vTUqlUradYODo6WllZPezFPM9bW1vreUToMEEQ8PlLRfvhIwilgi8fI2sWyKTdQlJfavogSq0WlErlsgjmrULh418Ui4MfE3KPTUGiTxDevn3bysqKYe4NTrt27aq9EKgPlUqVlpam44tpmtblJ4ROgs9fQtoPH0EoFfzyG9mVBvHZ/vRzT9Dk1w+foemvRtD/OCpoRFr/GRQd/7d0cnJqbGxsOVNaU1PTs2dPfctpj0uXLjU3S3dLY3lTq9VXr16Vugr5unbtWmNjo9RVyFRzc/OlS5ekrkJePO0obQoSQqqrq7VLdVKEvOdrgBQk+gShq6trt27dSktLtU9LS0tbTpMaR2pqanl5uTGPCC3y8/OnT58udRXy9frrr+fl5UldhUydOnVqwoQJUlchXx999NH69esNu0+dgrC5uTk7O/vAgQOiKGZnZ+/bt48QolQqp0yZMnfu3KtXr2ZmZubm5k6ePNmwxQEAAHQ2na4RNjY2Ll68mBASGRm5ePFiJyeniIgIQshf//rXd955JyQkpGfPnuvWrcMEQQAAMDs6BaG9vX1WVtaD221sbFauXGnokgAAAIynHfMIO9vNmzcTExOvX7+u4+tramrs7OxY1pDLpYKONBpNY2Njt27dpC5Epmpra7t06aJUKqUuRI54nr9z5452/Ugwvrq6OoVC8YiZdffZsmWLt7f3o19jQkFICLlx4wbu3AsAAIbi7u7+2D8ZTSsIAQAAjAxzQgEAQNYQhAAAIGsIQgAAkDUEIQAAyJpZzj24du3azz//fOnSpTFjxrR5B0ToPI2NjTt27Dh8+DBN09HR0SNHjpS6InnJy8vLy8u7deuWq6vrCy+8YOQFfkGroqIiJycnKSkJq4gY04EDB44cOdLy9OWXXzbU9DmzHBGOHDnyr3/96+zZs4uKiqSuRXa+/vrrL774wsbGRqlUTpgwYcmSJVJXJC/r1q3TaDT9+vUrKSnx9fW9fPmy1BXJjkajeeGFF2bOnFlRUSF1LfKSnp7+3XffnfmVIAiG2rNZTp8QBIGm6aFDh86ZM6f1/QvBCO7evdsylfX7779fuHAhlj6XytChQ999912s8Wtk8+fPFwRh2bJla9eubX0jcehsf/rTn5RK5cKFCw2+Z7McEeJOYBJqvaDD3bt3cZ9uqVRVVV2+fHnIkCFSFyIvJ06cSE9PnzNnjtSFyFRZWdnixYt/+OEHtVptwN0iUaCDqqurP/7449mzZ0tdiOzMnz/f3d3dy8vrww8/DAwMlLocGeF5/qWXXvrqq690X98LDKhXr17u7u537tz54osv/Pz89L8VfAuzbJYByd25cychISElJSU1NVXqWmTn/fffnz59+sGDB6dOnern5zd69GipK5KLzz77LCgoKDw8XOpCZOrdd9/VPhBFceTIkUuXLp03b55B9owRIbRbfX19fHx8cHDw559/LnUtcmRtbd2zZ8/k5OSJEyemp6dLXY6MrFmzZs+ePUFBQUFBQbdu3Zo2bdqqVaukLkqOKIoaPnz4mTNnDLVDjAihfRobG5OTkwcOHLh06VKKoqQuR154nuc4TqVSEUI4jispKUlJSZG6KBlZs2ZNy6WpmJiY999/Py4uTtqSZEWtVnfp0oUQ0tTUtGvXrkmTJhlqz2bZNfrWW28VFhYeP37cxcXF0dFx+fLlQUFBUhclFwsWLJg7d+7QoUO1LUtKpbKgoEDqouSitrb2iSeeGD58eLdu3fbv3+/i4rJz5070K0miZ8+e6Bo1sv79+w8aNMjBwSE/P79Pnz47duywtrY2yJ7NMghPnTp1586dlqdeXl52dnYS1iMrly9fvnLlSstTiqICAgIkrEduLl++/PPPP6vV6v79++PvPwmVlZX169cP3zzGdOnSpcOHDzc2Nmp/+Q14RsosgxAAAMBQ0CwDAACyhiAEAABZQxACAICsIQgBAEDWEIQAACBrCEIAAJA1BCEAAMgaghAAAGQNQQgAALKGIAQAAFlDEAIAgKz9P91b3mJA4OdMAAAAAElFTkSuQmCC", "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" ], "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" ] }, "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.7.3" }, "kernelspec": { "name": "julia-1.7", "display_name": "Julia 1.7.3", "language": "julia" } }, "nbformat": 4 }