{ "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 Input and output formats 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.846721755867 -0.70 0.80 4.6\n", " 2 -7.852284510527 -2.25 -1.54 0.80 1.0\n", " 3 -7.852621481395 -3.47 -2.53 0.80 3.2\n", " 4 -7.852621950730 -6.33 -3.38 0.80 2.2\n", " 5 -7.852621959566 -8.05 -4.61 0.80 1.6\n" ] } ], "cell_type": "code", "source": [ "scfres = self_consistent_field(basis; tol=1e-8, 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+gvaeTAAAgAElEQVR4nO3daUAUV9o+/FNVp9lBUXZQIwiiAiqCLBIXFg0C7hjNG9ckmnXIY+LoTDIZoklGsxNjEsw7JiGTZB4URHBFIKi4h0VUQDGuATdUkJ2uqv5/6DzEKGoDTVd31/X71JTdVTcN9sVd55wqRqVSEQAAALlipS4AAABASghCAACQNQQhAADIGoIQAABkDUEIAACyhiAEAABZQxACAICsIQgBAEDWEIQAACBrOgpCURSPHz9eXl6um8MBAABoiNHNJdbmzp2rUCju3Lnj6en5wQcfdPichoaG999/f9WqVRruU6lUKhQK7dUIncPzPKVU6ipkCm++tPDhIyFBEFiWZRhGi/vURUdYWFhYXV2dkpKSnp6+ffv2qqqqDp9WU1Pz/fffa77blpYWLRUIXYH3X0ItLS24SrCE8Msvoba2NlEUtbtPXQRhcXFxcHAwIYRl2VGjRpWWlurgoAAAAJrQRRDW1dVZWlqqH1tZWdXW1urgoAAAAJrowSBsDzwHB4eamhr14xs3bjg5OfXcQQEAADqlW0GYm5sbFxfn5uY2efLku7cfOXLE3d198ODBbm5ueXl54eHhe/bsaW5uvnnzZnFxcVBQUPdqBgAA0JpuBaGJicns2bMXLlxYX1/fvlGlUs2bN2/lypXXrl379NNPn376aQcHh2XLlk2YMCE6Ovrjjz+2sLDo8hHrlaRF6GD7DQxdAwBAl2hh+cTXX3+dkpKyf/9+9ZcHDhyYMmXKtWvX1NO7Bw4c+Pnnn8fExDxyP2VlZQEBAeppNWrz58+fNWvW3c/ZXsV+8yv3nzFKM440NDRYWVkRQr79ldt7nf0mRNnNbwQ6pbGxsX3oF3SssbHRwsJCuzPIQXPtHz6ge83NzSYmJhzHafh8MzOzR6410v5SpHPnznl5ebUf2Nvb++zZs5q80NzcvFevXn/729/atwwePPie9jHekzQS1dMH2bQIztxcsLCw2HhGteOKuDmCM+OwrEenRFHsTnMP3aF+8xGEUhEEAb/8UmEYplNByLKPPvGp/SCsq6u7+1fE2tpaw2miDMOYmZlFRUU9/GmLBxOGYWblid8FsZsryZaLqvQohZmm7wloDcuymvyGQU9Qv/kIQqngl19C7P/R4j61H4T29vZ1dXXtX9bW1jo4OGj3EIu8WBUhQdmmA6zE3U9QU6QgAAB0lfb/qBk2bFh5eXlTUxMhRBCE4uJiHx8frR9FEIkokkPXVCuPCXVtWt89AADIRbc6wpqampKSkvLy8tra2pycHAcHBz8/Px8fn5EjR65cuXL58uXJyclOTk5hYWHaKlft36fFjIviL9Gt/622/KhU+OGs+MYI7qWhLMW5CgCQTl5e3ldffSV1FUZr0qRJzzzzTE/suVtBePbs2bVr1xJCnJyc1q5d+/jjj/v5+RFC/vd///d//ud/wsPDhw0blpGRod2RjH+fFtMviGmRVNlEXhnKWlGy8bS44zfxqwrxwyAuph9GTQBAGkePHlUqlU899ZTUhRihgwcP7tmzRx+DMDg4eM+ePfdvd3V1TU1N7c6eHyTjophxUUyPpKYcUa+WUI8XZv+mWh7CvnpY+PQk+SSY87FFHAKABLy9vePj46WuwgiJorhly5Ye2rmB3ckl0oWNdmPvmR2z2IudOoD0NSXF0+k3Z8SJO/m4/uw7AZy9mURVAgCA4TCwUTUrBelwjmhfU0IIUbBkiTd7cqbC1pT4pinXHhfbtHyzDgAAMDYGFoSa6GNK1gRye2NpwTXRN43fdB5hCAAAD2SEQag2uBeTNZGuD+VWF4uRO/jSW7iLKQAAdMBog1At0pUpmkZnu7OTdvLz84XrzVIXBABwn0aePGQcR0VILVZL9yQjD0JCCGXJEm+2bJbCxZL4pSvXHhdbO7p/BQCAVObnC7NzhQ6zUEVIwiEhfDsv4qxWjzH+IFSzNSVrArn9sbSwRuWbjoFDANAj//RnD1wTp+/h7/kzXUXIq4eEDRXi6gCO7dKiMHd39+vXr2ulSFEUc3JyBKG7nYRSqbSzs1NffexBRo4cWVlZ2c0DaU4uQajm2YtJjeC+GsO9WyyGb+dLbuJPLACQnl8fJncyPXpDNSPnjyxUp2ByhZgWSbt8qZDW1tbu32tPjef5qKio5ubujjBxHPfWW2+ZmJg85DltbW2iqLt2xcDWEWpFuAtTNJ3+56w4eTcf6cK+H8Q5mUtdEwDImzoLI3bwM3L49EhqwmkhBe+XmppaVFTk7u6+YMECU1NT9cbc3NyCggI7O7u5c+f26dOHEJKfn9+3b19fX19CSFVV1ZEjR2bMmJGRkUEI+eabb0xNTePi4pydnTU5YlNT04YNG27cuBEcHBwXF6feaGb2+yrvzMxMX1/f/fv3V1ZWhoeHT5gwQVvfaafIqyNsxzJkvid7Jl7hbkN8NisTi4QOb3wPAKAzd/eFrxzUfgomJCR88cUXgwYNysrKio2NVW987733XnzxRXt7+5MnT44cOfLmzZuEkO+++y43N1f9hIqKivfee69rR2xpaRk3blxZWVn//v1Xrlz57rvvEkJ4nl+6dGlraysh5P33358+ffrp06d79+4dHx+fl5enhe+z8+TYEbazUpBEf26xF/vmL6JXKv9OADvPs2vn4QEAOpBTpVpb2rm/sj1smN1VKpVK5deX+fSk8OnJTrx2vic7b1DH7c3ly5dTUlIuXrxoY2OzYMECDw+Pffv2+fv7v/POO8eOHRs2bBghZNq0aevXr3/rrbc63MO0adMIIYsWLbKystKwns2bN1tbW2/YsIEQMnbs2MDAwNdee+2eWwlOmzYtMTGREHL79u2MjIzw8HBNv1vtkXUQqvW3YlLGc/lXVMsOC1+Ui58EcyEOSEMA0IJRdswKv07cMVVFyBdlIktUCo5YcsxrviztzKfRkN4P/KcTJ054e3vb2NgQQhQKRUBAwIkTJ3r37m1qaqpOQUJIWFjY0aNHO3E8QgghN27cWLZsGSGEZdnvvvvu7n8qLS0NDg7+vbYhQyil58+f9/DwuPs5I0aMUD9wc3Pbu3dvZ4+uFQjC3413Zn6ZRv9zVpyVI4xxZD4MYvtbIQ4BoFtsTUmkq6afJOrZMTt/E7dE0X6WJGIHv+6USn2Pge5rbGxsaWlp/7K5udnCwsLCwkI9m0Z9jyD1RkIIx3Hts0MfPr2TEGJpaTlz5kxCyP03GmpsbGyfFCOKYmtrq3r/d+M46W+tLtMxwg79PnA4mw61Jf5b+JXHhAal1DUBgDzcM0e0w3mk3XTixImTJ08SQn777bcDBw6EhYU99thjDg4OmzdvJoQ0NzenpqaOHz+eEDJgwICioiJCiEqlSktLU7/cxMTEwsLi1q1b9+zWwsJi2rRp06ZNmzp16v0HzczMbGxsJIRs2bLFxcWlX79+2vlmtApBeC9LShL9uZIZtLqRDNnMp1RiGSsA9KwOV0poPQt9fX2fffbZKVOmBAYG/u1vf/P09KSUbty48dVXX42Ojvb19fX19Z03bx4hZPHixfv37x83blxAQMDdyxheeeWVwMDAgICA4uJiDQ86ePDgMWPGxMbGvvDCC8nJyfcMEOoJnBrtmJslkzKeO3pD9T+Hhc/LxE+CuTGOOFMKAD1i2WHh69Ni1kQa9efzqO1rKp7ME9Iju7imvl3fvn137NhRUlLi5ubm4uKi3hgeHn727Nny8nJ7e/v2ds3V1bWiouLUqVNubm729vbtZ0fXrFmzevXqhoYGa2trDQ8aGhr60ksvnTlzZsiQIeoRShMTk2vXrllaWhJCdu7c2b6UYtGiReoY1j0E4cOMtmcK4ujm8+LT+UKgHfP+aPYxa8QhAGhZiAMzud+9KaimzsK9V1RamdFuYmIyevToezaam5v7+/vfs9HCwiIwMFD9WB1gagqFwtbWtlMHtbW1DQoKunuLg4OD+sHdgWpqatq+tFHHEISPwBASP5CN6ceuOyUGZfKLvNg3RnDWCqnLAgAjMtv9YScM/fowfn26GINubm7q2Sje3t5z5szp2k66LCoq6u4Q1ZyLi4tCobvPWQShRiwoWTGcfXoQs6pYHLqZ/8dI9tnBLJYcAoCeO3LkiPqBeghQx0dXLz3sgj179mi3kofTx3FLveVqySSHcVuiuO/PiqO38vuvYhoNAIDBQxB2WoAdsz+W/tOfXbBXiMvmz9UjDgEADBiCsIvi+rOnZtIwR3Z0Bp9wSLiDFYcAAIYJQdh15pSsGM6emKloEYhnqjLppCigOQQAMDSYLNNdzhYkOYxb6s0uOyJsqBA/CuaecMMsGgA5opR+++23Op7oIRO3b99uv2yp1iEItcPfjsmPoVmXxJcPCkN6k0+DOQ8bxCGAvCxdulSqO+rJgZubWw/tGUGoTXH92Sfc2C/KxOBM/ikPdtUortfDbsIMAEbF2tp61KhRUlcBnYYxQi1TsCTBh62IVxBCvDYpk06KvPjIFwEAgGQQhD2irylJCuHyY+juKtE3nd9xGbNoAAD0FIKwBw3pzeyYRNeFcMuPCFE7+VO3EYcAAHoHQdjjIl2Zkhk0fiAbtZNfWiDcaHn0SwAAQGcQhLqgYMkSb/bkTIUZR4ZtVq49LrZh4BAAQD8gCHWnjylJCuH2x9GCa6JvGr/pPMIQAEB6CEJdG9yLyZpI14dyq4rEqJ38iVsYOAQAkBKCUBqRrkzxdBo/kJ20i19aIFxvlrogAAC5QhBKhrJkiTd7aqbC1pT4pSvXHhdbBalrAgCQHwShxGxNyZpAbl8sLaxR+aZj4BAAQNcQhHrBqxeTGsF9NYZ7p1gM384fx8AhAICuIAj1SLgLUzydLvRio3fx8/OFqxg4BADoeQhC/cIyZL4neyZe4W5DfDYrE4uEFgwcAgD0JAShPrJSkER/rmg6PXeHDN7Ep1SKOFUKANBDEIT6q78VkzKe+24c98lJMSSTP3QdaQgAoH0IQn033pkpnEZfHMLOyhFm5wqXGhCHAADahCA0AL8PHM6mQ22J/xY+sUho5qWuCQDAWCAIDYYlJYn+XMkMeu4O8cLAIQCAliAIDYybJZMynkuL5JIrxKCt/IFrSEMAgG5BEBqk0fZMQRxd7sc+nS/MzhUu1CMOAQC6CEFoqBhC4geyp2bSUXZMUCa/8phQr5S6JgAAA4QgNGwWlKwYzhZNo7dbydDN/IYKESOHAACdgiA0Bq6WTHIYtyWK+/6sOHorv/8qwhAAQFMIQuMRYMfsj6X/9GcX7BXisvnzGDgEANAAgtDYxPVnT82kYY5sYAafcEi4g4FDAICHQhAaIXNKVgxnS2fSFoF4pSqTTooCmkMAgAdAEBotFwsmOYzbMYmmXxADM/i9VxCGAAAdQBAaOX87Zm8sfXsUu3ifEJfN/3oHcQgA8CcIQlmI68+Wx9NIFzY4k084JNS1SV0QAIDeQBDKhQlLEnzYingFIcRrkzLppMiLUtcEAKAHEITy0teUJIVw+TF0d5Xom87vuIwzpQAgdwhCORrSm9kxia4L4ZYfEaJ28qduIw4BQL4QhPIV6cqUzKCx/djx2/mlBcKNFqkLAgCQAoJQ1hQsSfBhT8crzDgybLNy7XGxDQOHACAzCEIgfUxJUgi3P47uvyr6pfGbziMMAUBGEITwu8G9mG2T6Oeh3KoiMWonf+IWBg4BQBYQhPAnka5M8XQaP5CdtItfWiBcb5a6IACAHoYghHtRlizxZk/NVNiaEr905drjYqsgdU0AAD0GQQgdszUlawK5fbG0sEblm85vOi+qCFlVLDbz9z4z/YJ49AbOowKAoUIQwsN49WJSI7gvx3DvFIuRO3hTVhWbzTfelYXpF8TPTolDejPS1QgA0C0IQni0CBemeDpd4MkmnRKVAnli1+9ZqE7BrInUWiF1iQAAXYUgBI2wDJnvyZ6JV4S7MiU1Kt805Y8XOKQgABgBBCF0gpWCJPpzxTMoS8jLR+nrvhxSEAAMHYIQOq30lsrNilkySJiVyz+zX2i6b/oMAIABoVIXAAamfVyQaW3xdzT5R6Hgc0X13VjucSfMlwEAg6SjjrCpqengwYPZ2dm6ORz0kHtmxywezH4UxFlyqvhcfuUxAcsNAcAQ6agjfOaZZ5qams6fP19aWqqbI4LWqQg5fF11z+yYOR5sm0j6mjIbz4ijMvhvx3EBdmgNAcCQ6CgIf/rpp0uXLsXGxurmcNATGELeH83dv32+J0sIienPbTovxu7mF3qxq0ZxJhh9BgADgY8r0Jr4gezxGYrTtSQggy++iWvNAIBh0HJHOG/evNra2ru3fPPNN3Z2dto9CugtR3OyJYrbdF58Yhe/yItdPYpT4G8tANBvnQvC27dvW1pampiY3L2xsrKytbV16NChLMt++OGHgvCnKRO2trZaKBMMSvxANsyRXVogBGbw343nhvfBqCEA6C9NgzA+Pj4nJ6e2tjY1NTU+Pl69UalUxsfHl5SU2NjYMAyTk5Pj6OjY4csbGxvr6uoEQegwSsH4OFuQzIlcSqUYtYN/zZd73Y/lkIYAoJc0PW+1cOHCoqIiHx+fuzempqaePXu2vLy8tLTU29t77dq1D3r5xx9/vGzZMhcXl9mzZ+fn53enYjAg8z3ZY9NodpUYlsWfrsOoIQDoI007wpiYmPs3pqamPv300+bm5oSQZ599dvHixR9++GGHL//HP/7xyEM0NTVVVVUNHDiwfcuLL774wgsvPOj5jY2NDIMuQzKNjY2aPK0vIRljybe/cmGZ4l+8xQRvnsUPrdsaGxtVKhV+/6WCDx8JNTc3m5iYcFwHM9g7ZGZmRukjkq5bk2UuXbr01FNPqR+7u7tXV1crlUqFootXn7SwsHB0dMzNzW3f4ujoaGlp+aDnq1QqKyurrh0LtELz9/+V4STWXbV4n7DrKv12LOfZCx8i3WVpaYnPYqngw0dCHMd1Kgg10a0gbGpqMjMzUz82MzMTRbG5ubnLQUgIoZS6u7t3pyTQWwOtmbwY+nWFGJrFv+7LLfdj0RoCgD7o1tx2R0fHW7duqR/fvHnTwsLCxsZGG1WBcWIIWeLNHp5Ct18Wx2/nf72DUUMAkF63gnDUqFEHDhxQPz548GBAQIA2SgIj52HD5MfQmY+xIZn8hgoRYQgA0tL01Oj27durqqpu3bqVk5Nz+/btuLg4Z2fn559/PjAwcPTo0U5OTm+//fYXX3zRo7WC0WAZkuDDRroyC/cK6RfE//9xzs0S50kBQBqadoSVlZWFhYXqi4UWFhY2NDQQQgYPHrxz5878/PxvvvkmKSlp+vTpPVgpGJ1htsyhKXSCM+u/hd9QIUpdDgDIFKNS6cupqQsXLkyYMOH8+fMaPr++vt7a2rpHS4KHaGho0NbEuRO3VAv3CY7m5OswzhWtoQYaGhowa1RC+PCRUGeXT2gCF4IE6fn2YQ5PoeOcWP8MtIYAoGsIQtALCpasGM7uiaZflotx2fyVJqkLAgDZQBCCHvHrwxydSkfZMaMylGnn0RoCgC4gCEG/KFiS6M+lR9I3C8XZuUJNi9QFAYCxQxCCPgp2YIqnU3cb4peuzLiI1hAAehCCEPSUGUfWBHKbI+iKo+LsXOFWq9QFAYCRQhCCXgt1ZEqmU3cb4pvGZ11CawgA2ocgBH1nTsmaQO5/w7llh8X5+UK9UuqCAMC4IAjBMIQ5MYXTqTklful8brW+XAUCAIwAghAMho2CJIdxyWHc4n3C0gKhAa0hAGgDghAMzERXpnQGJYT4pfM/X0FrCADdhSAEw9PLhCSHcetDuQX5wtICoZGXuiAAMGQIQjBU0f2Y0pmUEDI8nd93Fa0hAHQRghAMWG8TkhzGJYVwT/0sJBwSWgWpCwIAA4QgBIMX048pnk6rm4j/Fv7YDbSGANA5CEIwBvZmZFMElziKnZLNrzyG1hAAOgFBCMYjfiBbMkNxpo4EZPBFNWgNAUAjCEIwKo7mJD2Se8ufnbybX3lMaMNF2QDgURCEYITUrWF5LQnM4EtuojUEgIdBEIJxcjInW6O4N0eyk3bxiUWCgDQEgAdAEIIxix/IHptKD1xTjcniK2oRhgDQAQQhGLn+Vkx2NF3sxY7dxq89LqI1BIB7IAjB+DGELPFmj06lu38TH8/iz9QhDAHgDwhCkIvHrJncGLrQiw3L4tceF0WkIQAQQhCEICvq1vDwVLrjsjh2G3/2DsIQABCEID/u1szPMXS+JxuSySedRGcIIHcIQpAjliFLvNm9sfSHX8UndvKXG5GGAPKFIAT5GtqbORhHw13YUVv4DRW4CA2ATCEIQdYoS1YMZ/Ni6IYKMXoX/xtaQwD5QRACEB9b5vAUOt6ZHZWB1hBAdhCEAIT8X2u4J5p+VS7G7uarm9AaAsgFghDgD359mCNT6eNO7Mgt/Pdn0RoCyAKCEOBPFCxZMZzNmkj/VSLOzhVqWqQuCAB6GIIQoAOj7Zmi6dTdhvilK7dcQGsIYMwQhAAdM+PImkAuLZL+7Zg4O1e42Sp1QQDQMxCEAA8T4sAUq1vDNH7rRbSGAEYIQQjwCOaUrAnkUiO45UfF2bnCbbSGAMYFQQigkTGOTMl06m5DfNP5bZewuALAeCAIATRlQcmaQO6bsdzLB4WlBUK9UuqCAEAbEIQAnRPlypTOpIQQv3T+5ytoDQEMHoIQoNNsFCQ5jPtqDLcgX1haIDSgNQQwZAhCgC6a5MacmEkJIcPT+b1oDQEMFoIQoOt6mZDkMG5dKPd0vrC0QGjipS4IADoPQQjQXZP7/d4a+qXzBVfRGgIYGAQhgBb0NiHJYdy/Atn4XH7lMaFVkLogANAYghBAa+IHssdnKCrryKgM/pcatIYAhgFBCKBNDuYkLZL7pz8bu5tfeUxow0XZAPQeghBA+9StYUUtCcjgi2+iNQTQawhCgB7haE4yorh/jGSf2MWvPCYo0RoC6CsEIUAPih/IlkxXnLqtCsviy2vRGgLoIwQhQM9ytiBZE+lLQ9lx2/i1x0UBaQigZxCEALow35M9No1mV4lhWfzpOoQhgB5BEALoyAArJmcyXeTFPp7Frz0uikhDAP2AIATQHYaQJd7skal012/i49v4SrSGAHoAQQigawOtmbwYusCTDUVrCKAHEIQAElC3hoen0O2XxUm7+EsNCEMAySAIASTjYcP8HEMjXdiADH5DBTpDAGkgCAGkxDFkxXD25xj6dYUYvYv/rRFpCKBrCEIA6Q2zZQ5NoROc2VEZ/IYKXIQGQKcQhAB6gbJkxXA2J5omV4gxu/kqtIYAuoIgBNAjvn2Yw1PoWCfWH60hgK4gCAH0i4IlK4az2dH0y3IxPle40SJ1QQDGDkEIoI+G92EOTaEeNmR4ujLtPFpDgB6EIATQU2YcWRPIpUfSNwvF2blCbrVq032J2MyTxCJcxxugWxCEAHot2IEpmkbdbci8fOHdEuGHs39kYatAnswTBlozjIT1ARg+BCGAvjOnZE0gtzmCa1Iyrx0RkstFQkirQOJzhZkDmQWe+F8M0C1U6gIAQCOhjszxGfTNQiHhsFBaQy83i0hBAK1AEAIYDHNKPgriwhzZJ/NU7tZipAv+/wJoAf6cBDAkrQL55oz4r5G8UiTem/jEIqFFkLomAAOHIAQwGO3jgks9xdIZdKQdk1ul8krlUyqxvgKg6xCEAIbhntkxFpTseoKaUTLPk/nohBi+nT9xC8soALpCF2MMKpUqOzv7wIEDlNK5c+d6enrq4KAARqbkpmrWQGb+XbNjLCjJiKKrioSi6fSHs+KkXXxcf/bdAM7OTMIyAQyPLjrC1tbW1NTUgIAADw+PiIiI6upqHRwUwMgEOfwpBdUsKVk7muMYMt+TPTVTYcaRoZuVSSdFrLEH0JwuOkIzM7N///vf6sdbtmw5fvy4i4uLDo4LICu2piQphHthCPvqYSG5Qvw0hJvoiqX2AI+m0zHCqqqq0tLSoKAgXR4UQFa8ezO7nqBrR7PPFwhx2fyFevSGAI+g5Y4wPj7+3Llzd2/58ccfBw8eTAipra2dNWtWUlJSnz59tHtQALhHXH820oX97JQYlMkv8mLfHMFZKaSuCUBfdSIIa2trf/311/79+9vb27dv5Hn+559/rq+vHzt2rJ2d3aZNmzp8bX19/ZQpU5YvXx4dHd3dkgFAA+aUrBjOPj2I+dsxcchm/t0Adp4ni1OlAPfT9NTopEmTHB0dw8LC0tLS2jcqlcqIiIg33ngjNTV1yJAhJSUlHb62ra1t8uTJwcHBAwYMKCwsvHnzphYKBwANuFoyKeO5tEjui3Jxwnb+OJZYANxH0yBMSkqqr68PDg6+e2N6evrNmzcLCgr++9//vvTSS2+//XaHr21tbR06dGhdXd2GDRs2bNhw5syZ7lYNAJ0x2p45GEcXe7FP7OTn5wvXm6UuCECfaHpq1Nvb+/6NmZmZ06ZNMzExIYTMnj37vffeEwSB47h7nmZtbZ2cnPzIQzQ1NV28eJG565YyiYmJr7322oOe39jYyOD+M9JpbGyUugT5amxsVKlUnf39n+FMIiYzH5Wxvmncy4PFlwbzJriiRpfgw0dCzc3NJiYm9wfNg5iZmVH6iKTr1mSZqqqqxx9/XP3Yzc1NqVReu3aty0sjLCwsBgwYcP78eQ2fr1KprKysunYs0Aq8/xKytLTswmexFSEfjyFLfVTLDgs/7KafBHPR/fCB3mn48JEQx3GdCkJNdOsPQp7n26tRR25bW5sWigKAnjS4F7N9Ev0giH3lkBCXzZ/DEguQt24FoZOT040bN9SPr1+/zjCMs7OzNqoCgB4X158tn0UjXdjRGXzCIaFeKXVBABLpVhA+/vjjubm56sc5OTlBQUGmpqbaqAoAdEHBkgQftnQmbRHI0M18SqWI3hBkSNMxwh9//LG0tPTcuXMZGRkXLlyYN2/esGHDFixY8MEHH7z44ove3t6rVq1qv44aABgQFwsmOYz7pUb1l4PC5y7hZb0AABm8SURBVGXiZyFcsAMGDkFGNA1CR0dHd3f3N954Q/2lpaUlIaR3795Hjx7duHHjxYsXt2zZ0j5xBgAMToAdc2AK/b5SnJHDR7qwHwRxjuZS1wSgE4xKpS/nQi5cuDBhwgTNZ43W19dbW1v3aEnwEA0NDZg4J5WGhoauzRrVRCNPPigVPj8lvjyM/dtwzlSbs/OMBD58JNTZ5ROawDIiAPgTS0oS/blDU2hhjco3nd92SV/+VgboIQhCAOiAZy8mayL9IpRbeUyI2smX1SIOwWghCAHggSJdmeLpNLYfO24bn3BIuIMlFmCMEIQA8DDqJRYV8QpCiGeqMumkiDUWYGQQhADwaH1NSVIIt3MS3XReHL2VP3ANYQjGA0EIAJryt2P2x9EVw9n/72dhdq5wqQFxCMYAQQgAncAQEj+QLZtFh9oS/y18YpHQIkhdE0D3IAgBoNMsKEn054pn0HN3iG8av+m8KHVFAF2HIASALupnyaSM55LDuNXFYuQO/uRtnCkFg4QgBIBuCXdhiqbR2e5s1A5+aYFQ0yJ1QQCdhCAEgO6iLFnizZ6apTDjyLA0ZdJJUUBzCIYDQQgA2tHHlCSFcLufoOkXxMAMfv9VhCEYBgQhAGjTiL7M3lj69ih2/l4hLpu/iCUWoPcQhACgfXH92VMzaZgjO3orn1gkNPNSFwTwYAhCAOgRFpSsGM4WTqPn7pDBm/iUSiyxAD2FIASAHuRmyaSM574fz318QpywnS+9hTOloHcQhADQ48Y5M0XT6SIvduJOfn6+cANLLECfIAgBQBdYhsz3ZMtnKWxNybDNyqSTIo9zpaAfEIQAoDu2piQphNsXS3f9Jvqm87t/w5lSkB6CEAB0zbs3s/MJ+v5o9sUDQlw2f74ecQhSQhACgDTi+rPl8TTShQ3O5FceE+qVUhcEcoUgBADJmLAkwYctnk6rG8nQzXxKpYjeEHQPQQgAEnOxYFLGc2mR3JflYvBW/sh1pCHoFIIQAPTCaHvmQBx9aSg7PYefny9cb5a6IJANBCEA6Av1EouyWQoXS+KXrlx7XGzDEgvoeQhCANAvvU3ImkBuXyzdf1X0TeN3XMaZUuhZCEIA0EdevZhtk+j6UO71I0LUTr68FnEIPQVBCAD6K9KVOT6DxvZjH8/iEw4Jd7DEAnoAghAA9JqCJQk+7ImZihaBDN3Mb6gQscYCtAtBCAAGwNmCJIdxGVHcd5Vi0Fb+EJZYgPYgCAHAYATYMQVx9K/D2bl5wvx84SqWWIA2IAgBwJAwhMQPZE/Nou42xGezMrFIaBWkrgkMHIIQAAyPJSWJ/tzhqbTsNvFN5zedx3pD6DoEIQAYqkE2TGoE9+UY7u0iMWonf+o2Bg6hKxCEAGDYIlyY4uk0th87fjufcEioa5O6IDA0CEIAMHjqJRYV8QpCiNcmZdJJUUBzCBpDEAKAkehrSpJCuJ2T6OYL4uitfMFVhCFoBEEIAEbF347ZH0sT/dl5e4XZucKlBsQhPAKCEACMUFx/9tRMOtSW+G/hE4uEFiyxgAdDEAKAcbKgJNGfK55Bz90hXql8SiWWWEDHEIQAYMz6WTIp47nvxnMfnRAjdvAnbuFMKdwLQQgAxm+CM1M4jT7pzk7axS8tEGpapC4I9AmCEABkgbJkiTd7cqbCjCNDN2OJBfwBQQgAMtLHlCSFcHtj6c7fxIAMfh+WWACCEABkaEhvZtcTdNUoduFeIS6bv1CPOJQ1BCEAyJR6iUWYIxuUya88JjQopS4IJIIgBAD5MqdkxXC2aBqtbiRDNvMplSJ6QxlCEAKA3LlaMinjuf+M5z45KU7Yzh/HEguZQRACABBCyDhnpnAaXezFPrGTn58vXG+WuiDQFQQhAMDvWIbM92TL4xUulsQvXbn2uNiGy9HIAIIQAOBPepuQNYHc3li696rol8bvvKwihGReFPn7QvFcvarkJs6jGjwEIQBABwb3YnZMoh8EsS8fFOKy+WM1qrk/C8q7svBig2pWjoBV+UYAQQgA8EBx/dnyeBrpwn5ZJl5vVs3O/T0LLzaopu8Rvn6cG2XHSF0jdBeCEADgYUxYkuDDls6k3r2Z/CuqMVn8rw0MUtCYIAgBAB7NxYJJDuO2T+IuN5KQXSbvj0YKGg8EIQCAplwtibM58bMl0/YIxTUYHjQSCEIAAI20jwvmRLTOGsiEZPHpF7C6whhQqQsAADAAd8+Oqa8n347j7M3IvHyxxFe1ahQndXXQLQhCAIBHe7tI/GYcN7zPH+OCHwRxVgrhu0qxqpF8GcaZ4PyawcKPDgDg0TaO/VMKqv3TnzsxU1HbRsK387gkm+FCEAIAdJ0lJZsjuUhXJjSLP3Ub02cMEoIQAKBbGEIS/bl3A9jwHXzWJUyfMTwYIwQA0IIn3dkBVsysXKHsNlkxHD2GIcFPCwBAO4IdmCNTuNTz4nP7/3RVUtBzCEIAAK1xtWT2xtCaFhK+g7/RInU1oBkEIQCANlkpSHoUF+HChGby5bWYPmMAEIQAAFqmnj6zahQ7fju//TKyUN8hCAEAesRcD3ZLJF2yX/jsFAYM9RqCEACgp4Q6MvvjuA0V4tICTJ/RXwhCAIAe5G7NHJpCrzSRmN18bZvU1UBHEIQAAD3LWkG2RHH+dszorfzpOgwZ6h1dLKjnef7NN988ffp0r169nnvuuTFjxujgoAAA+oNjyJpAbpCNOG4b/+MEGu6Cm/rqER11hOPGjVu3bt1zzz335JNP1tTU6OagAAB65dnB7OYI+nQ+v74MA4Z6RBcdIaU0OjqaEOLq6mpra9vQ0GBnZ6eD4wIA6JswJ2Z/LI3LFkpvqdaHchTDU3pARz+E1tbWqKgoHx+fOXPmPPbYY7o5KACAHvKwYY5MpVWNqthsvg7TZ/SANjvCtra2kJCQezYWFhYSQkxMTFJTU8+cObNo0aI5c+Z4eHho8bgAAIbFWkG2TqTLDgthWXzmRG6gNYYMpaRpEN65c6e4uPj69esRERF9+vRp337x4sXNmzezLBsfH+/m5qaOvfsxDGNraxsUFBQUFFRSUoIgBACZ4xiSFMJtqBBDM/n/htNxzshCyWgUhG1tbfb29t7e3mVlZQUFBUFBQertlZWVQUFB8+fP53l+5MiRv/zyy4ABA+5/eUVFRU5Ojre396+//pqfn/+vf/1Lm98BAIDBWuLNDrRmnszj1wRyC70wYCgNjYJQoVDU1taam5v36tXr7u2ffPLJnDlzPv30U0JIU1PTunXrPvzww/tf7uTkRCnNycmxt7cvKChwcnJ60IF4nj937lz7l46OjpaWlpp+KwAABijKlSmIo3HZQvFN1SfBHIvOUOc0CkKGYczNze/fnpOT88knn6gfT548efXq1R2+vHfv3s8///wjj9LU1HTt2rWIiIj2LQkJCc8+++yDnt/Y2Mgw+JWRTGNjo9QlyFdjY6NKpcLvv1S0/uHjxJLscGbeARq9k/82lLemWHT/QM3NzSYmJhzHafh8MzMzSh+RdN2aLFNdXd3e3jk5OVVXV3dnbxYWFq6urufPn9fw+SqVysrKqjtHhG7C+y8hS0tLBKFUeuLDx4qQnBjy6mFhUh6XNZEbYIUfbsc4jutUEGqiW6ekWZYVxd+XhYqiqN3KAABkhbLk81Du5aFs8FZ+/1U0hbrTrSB0dna+cuWK+vGVK1ecnZ21URIAgHwt8Wa/HUfjc/mUSlx9Rke6FYTR0dFbt25VP87MzFRfPgYAALpjkhuzL5a+WyImHBJEdIY9T9MxwoSEhCtXrjQ3N7/xxht9+vRZv369vb39q6++GhQUtHDhQp7n8/Pzjx071qO1AgDIhFcv5uAUOiuHn50npIzjLHRxNUz50vTdffLJJ5uampYsWaL+Uj1Q/Nhjj504cSIrK4tl2aSkpL59+/ZUmQAAMtPXlOyOps8XCKGZfOZErj+mz/QYTYMwNDS0w+0ODg7PPPOM9uoBAIDfmbBk41gu6aQ4JktIj+QC7ZGFPQIXMgAA0GsJPuyGMC42m//hLKbP9AgEIQCAvovux+RE038UiiuPYfqM9iEIAQAMgG8f5uhUeuiaak6e0MRLXY1xQRACABgGOzOSHU3NKQnL4i83ojHUGgQhAIDBMOXId+O4BZ7smEzhlxpkoXYgCAEADEyCD/tpCBuzm99yAdNntACrNAEADM+Mx9hBNszUPcLxW6p/+nNYV9Ed6AgBAAySXx/mYBzdeVk1N09oxvSZbkAQAgAYKmcLsi+WKlgSvoO/1ix1NQYLQQgAYMBMOZIynps2gA3M4IswfaZLEIQAAIaNIWTFcPbjYDZ6N7/1IqbPdBomywAAGINZA3+fPlN8U5Xoj9ukdwI6QgAAIzGiL3NoCrf9kuqpn4UWQepqDAeCEADAeLhYMPvjKENIxA7+OqbPaAZBCABgVMw48p8JXJQrE5LJn7qN6TOPhiAEADA2DCGJ/tx7gWz4Dj7rEqbPPAImywAAGKcn3dkBVszMHKFsGFkxHG3PA+GtAQAwWsEOzNGpXOp58bn9ghKd4QMgCAEAjJmrJbM3hta0kPAd/I0WqavRSwhCAAAjZ6Ug6VFchAsTksmX12L6zL0QhAAAxk89fWb1KHb8dn77ZWThnyAIAQDkYq4HuyWSLtkvfHYKA4Z/QBACAMhIqCNTEMdtqBCXFmD6zO8QhAAA8jLQmjk8hV5pIjG7+do2qavRAwhCAADZsVKQLVGcvx0zeitfIfvpMwhCAAA54hiyJpD7qx87fjufVy3rLEQQAgDI17OD2c0R9Ol8fn2ZfAcMEYQAALIW5sTsj6Xry8SlBQIvyzREEAIAyJ2HDXNkKq1qVMVm83Xymz6DIAQAAGKtIFsn0hF9mdFb+TN18hoyRBACAAAh/zd95jVfdtw2Pv+KjLIQQQgAAH9Y4s1+P57OyeO/PSOXAUMEIQAA/EmkK1MQR9eWigmHBFEGnSGCEAAA7jXIhjkQR0/eVsVm83eUUlfTwxCEAADQgT6mZPcT1MOaCcviL9Qbc2OIIAQAgI5RlqwL5V4eyoZk8vuuGm0WIggBAOBhlniz342ns3P5lErjnD6DIAQAgEeY6Mrsi6Xvlhjn9BkEIQAAPJpXL+bgFFp6SzU7T2jkpa5GqxCEAACgkb6mZHc0tVGQMZn8pQbjaQwRhAAAoCkTlmwcyy3yYsdkCUdvGEkWIggBAKBzEnzYrx/n4rL5/5w1hukzCEIAAOi0J9yYnGj6VqG48pjBT59BEAIAQFf49mGOTqWHrqnm5AlNhjx9BkEIAABdZGdGsqOpBSVjsvjLjYbaGCIIAQCg60w58u04bqEnG5op/FJjkFmIIAQAgO5K8GGTw7iY3fxPvxre9BkEIQAAaMHkfsyeaPr3X8TEIsGwGkMEIQAAaIdfH+ZgHN15WTU3T2g2nOkzCEIAANAaZwuyL5YqWBK+g7/aLHU1mkEQAgCANplyJGU8N20AOzqDLzKE6TMIQgAA0DKGkBXD2Y+D2ejdfMZFfZ8+Q6UuAAAAjNOsgewgG2bqHqHkpirRn5O6nAdCRwgAAD1lRF/m0BRu+yXVUz8LLYLU1TwAghAAAHqQiwWzP44yhETs4K/r5fQZBCEAAPQsM478ZwIX5cqEZPKnbuvd9BkEIQAA9DiGkER/7r1ANnwHn3VJv6bPYLIMAADoyJPu7AArZlauUHabrBiuL52YvtQBAAByEOzAHJnCbTovPrdfaNOPzhBBCAAAOuVqyeTH0JoWErGDv9EidTUIQgAA0D0rBUmP4iJcmJBMvrxW4ukzCEIAAJCAevrM6lHs+O389stSZiGCEAAAJDPXg90SSZfsFz47JdmAIYIQAACkFOrIHJnKfXtGXFogKKVIQwQhAABIzM2S2RdLrzaTmN18bZuuj44gBAAA6VkpSHok52/HjN7KV+h2+gyCEAAA9ALHkDWB3Ao/dvx2Pq9ad1mIIAQAAD3yzGA2LZI+nc+vL9PRgCGCEAAA9MsYR2Z/LP2iTFxaIPA9n4YIQgAA0DseNszhqbS6SRWbzdf18PQZBCEAAOgjawXJiKIj+jIjtvCvHrr3rr4tAnnpoNCqjZv96i4Ib926FRMTs2fPHp0dEQAADJp6+szK4ezXp8W5P/8Rem0imZ0r+NoyppwWjqK7IPzrX//a0NBQXV2tsyMCAIARWOrNZkTRzIti9G6eENImklk5wuR+zPNDtBNhOrof4Y4dOxwcHCwsLHRzOAAAMCZRrkzJDBq0lQ/dyblaqmL7s9pKQaLdIGxpabmn4TMzM3Nxcblz585HH320bdu2FStWaPFwAAAgH542zMmZisGpSguq0mIKEu0G4dmzZ99+++27t3h7e69evfrNN9986qmnrly5UldXd+PGjfr6emtray0eFwAAjF6bSJ4vEN71Fyvqub8eFd4frY3hQUKIhkHY2tp67NixkpKS5ubm5cuX3/1P//3vf3ft2uXk5PTyyy/7+Phs2rTp/pc/9thjR48ePXr0aGFh4aVLl8aMGRMSEqKd8gEAQAbaxwUXPCYqTGjCEaLFLNSouzxy5MgLL7ywc+fOv//973dvX79+/d///veoqKjm5uYxY8Y0Nzd3+PJly5YlJycnJyeHh4cvXLgQKQgAAJq7Z3YMQ8jnoVwjT/56VBuLJwhhVCpNr+d24sQJf39/pVKp/lIUxUGDBq1bty4mJoYQEhQU9MILLyxcuPAhezh9+rSNjY2zs3OH/3r+/PnQ0NCUlJT2LV5eXv369XvQ3nCKVVoNDQ1WVlZSVyFTDQ0NlpaWDMNIXYhM4cNHx87Xqw5dJ095MISQ5uZmExMTjuNUhHxyUvXSkEesoGDZR/d7XQ/C6upqV1dX9X9IQsiKFStu3769YcMGDfd2v7KysoCAgODg4PYt8+bNi4+Pf9DzP/vsswULFvTq1avLR4Quq66u3rVr1+LFi6UuRKY2btwYFRX1kD8ToefU19dv3LgxISFB6kJkKj093cvLy8fHR8Pnm5mZUfqIQcCuT5a5evWqubm5OgUJIfb29mVlZV3eGyHEwsLC0dExLy9Pw+f/9NNPsbGxrq6u3TkodM3Fixe3bt36l7/8RepCZGrbtm2enp5DhgyRuhA5unjx4g8//PDGG29IXYhM5eXltbW13d0ydV/XZ6CamZkplcr2hrKlpcXc3FxLVQEAAOhI14PQxcVFEIT2hYOXL192c3PTUlUAAAA60vUg7N27d3h4+A8//EAIqaur27Zt2/Tp07VXGAAAgC5oNFnm5s2bo0ePbmtrq6qqGjhwoKOj48GDBwkhx44di42NDQwMPH36dHBw8Pfff9+dUiorK318fDRvK6uqqhwcHBQKRXcOCl3T0tJSW1vr5OQkdSEydfXq1V69emEwQhJKpfLatWs4ASaVmpoaMzMzzaesP/XUU6tXr374czQKQkEQLl682P4lpbR///7qx3V1dUePHrW3tx8xYoSGZT3EuXPnNH9ya2urqalp9w8KXYP3X0J486WF919CSqWS4zhNFkWoOTs7P/JPxk4snwAAADA+uDEvAADIGoIQAABkDUEIAACyhiAEAABZ09Ed6rXr2rVrv/zyS1VVVUREhIeHh9TlyEtTU9OuXbsKCwtZlo2Kiho7dqzUFcnLvn379u3bd+vWLRcXl3nz5jk6OkpdkRydOXMmPz8/Li7uQbcQgJ5w+PDh0tLS9i8XL178yIuIasggO8KxY8e+9957K1asOHbsmNS1yM4XX3yxbt06S0tLExOTGTNmJCUlSV2RvGzatEmpVLq7uxcXF/v6+rZf2gl0RqlUzps3LyEh4cyZM1LXIi9paWnffvvtuf8jiqK29myQyydEUWRZdsSIEStXrpwzZ47U5chLS0uLmZmZ+vH333//7rvvVlRUSFuSbI0YMWLZsmXz58+XuhB5WbVqlSiKX375ZWpq6rhx46QuR0aWL19uYmLy7rvvan3PBtkRar6UErSuPQUJIS0tLbgloVQuXLhQXV09bNgwqQuRl/Ly8rS0tJUrV0pdiEwdP3587dq1P/3004PuA981SBToopqamrfffnvFihVSFyI7q1atcnNz8/LyevPNN0eNGiV1OTIiCMKiRYvWr19/95+DoDPOzs5ubm537txZt26dn5/frVu3tLVng5wsA5K7c+dOTEzMrFmzHnLnZOghr7/++tKlS48cOfLMM8/4+fmNHz9e6ork4qOPPgoICAgLC5O6EJlatmyZ+oFKpRo7duxnn32WmJiolT2jI4ROa2homDx5cmBg4CeffCJ1LXKkvoX1lClTZs6cmZaWJnU5MvLjjz/m5eUFBAQEBATcunVryZIlGzdulLooOWIYJjQ0tFPXpn44dITQOU1NTVOmTBk8ePBnn33GMIzU5ciLIAg8z6sv98zzfHFx8axZs6QuSkZ+/PHH9qGpiRMnvv7669HR0dKWJCvNzc3qy2e3trZmZ2fPnTtXW3s2yFmjr7zyyqFDh8rKypycnPr06fPVV18FBARIXZRcvPPOO2+99daIESPUU5ZMTEzU9+QCHairqxs0aFBoaGivXr0OHDjg5OS0e/duzFeShKOjI2aN6piHh4e3t7etrW1BQUH//v137dplYWGhlT0bZBBWVlbeuXOn/UsvLy9ra2sJ65GV6urqK1eutH/JMIy/v7+E9chNdXX1L7/80tzc7OHhgb//JHT8+HF3d3d88uhSVVVVYWFhU1OT+pdfi2ekDDIIAQAAtAWTZQAAQNYQhAAAIGsIQgAAkDUEIQAAyBqCEAAAZA1BCAAAsoYgBAAAWUMQAgCArCEIAQBA1hCEAAAgawhCAACQtf8Hv+mgHE0nUZAAAAAASUVORK5CYII=", "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.8.3" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.3", "language": "julia" } }, "nbformat": 4 }