{ "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.844807180103 NaN 1.97e-01 4.3 \n", " 2 -7.850313458268 -5.51e-03 2.93e-02 1.0 \n", " 3 -7.850616387864 -3.03e-04 2.79e-03 1.5 \n", " 4 -7.850646669310 -3.03e-05 1.29e-03 2.5 \n", " 5 -7.850647322268 -6.53e-07 6.18e-04 1.0 \n", " 6 -7.850647503335 -1.81e-07 3.67e-05 1.0 \n", " 7 -7.850647511548 -8.21e-09 6.38e-06 2.5 \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/8PfnnDMMA0KyCAwgiknlgjuaiiKCCApqhls3TVv0pt20bjf7WnmtR920+lm2XrW6ZamV5gYKogLumribkgsKiqKCLLINc5bfH2OE5gI4zJnl9fwLzpyZeTOO85rPOZ/PeTNFUQgAAMBRcWoXAAAAoCYEIQAAODQEIQAAODQEIQAAODQEIQAAODQEIQAAODQEIQAAODQEIQAAODQEIQAAODQEIQAAOLQ7BmF1dfXrr78+evToX3755S73P336dEVFRf2fT5blBlRn1/BS1JJlGZf6M1EUBS+FCV6KWoqi4OOiVlO8FHcMwpdeesnV1fXjjz+eP3/+zp0777TbM888k5WVVf/nq6ysxJvbpLq6Gm9uE4PBIEmS2lVYhZqaGlEU1a7CKoiiWFNTo3YVVkGSJIPBoHYVVkGW5erqarM/rHCnG1atWpWbm+vs7Dx16tRly5b17dvX7M8NAACgutuPCEtKSnQ6nbOzMxG1atUqLy/PslUBAABYyO2DUBCE2kNVoihqNBoLlgQAAGA5fx4azc7OXrVqlU6nGzdunJ+fHxGVlpY+8MADJ0+ebNu2rXoVAgAANKEbI8KsrKxevXpVVFScOHGie/fuV65cmTJlyowZMzZu3Pj5559PmjRJ3SoBAACayI0R4bx586ZPn/72228T0YgRIxYvXvz666//9NNPv/7665IlS9q1a3en+8uynJube+LEibobvby8fHx8an/NvKQsPyP/N5xndfY5Vaq8tEdKGizU3QgAAGBhzLSYwdPTc8OGDY8++igRLVq06Mcff0xPT6/P/bt27VpUVKTT6epuHD9+/PTp0+tuefeYcKmKfdLDWFVZ4eLiklPOPbXbaVGvmvYPOO5SiqqqKq1Wy3G4pgFVVVVpNBpBuOMcZsdhMBg4jsNZeSIyGo2yLGu1WrULUZ8oikaj8ZaPWccky7LBYGjQS+Hi4sLz/N33EYioqqqquLjY19fXtMnX1/fixYv1fA53d/ePP/44IiLi7rvN7U2z90v/PKyZ34UVSC6T9spLB/IdPZzq+Sx2ied5Z2dnBCERCYKAIDRxcnJCEJogCGshCGvJsqzRaFxcXMz7sBwR8TzPGKtd3C3L8j3zsxHe7s776ejpXcLodPmHAXxHDxwTBQAA9XFE5OTk5OXldenSJdOmS5cu6fX6pniyp0K4XYX8uXIF11MBAAArceNgVFxc3Jo1a8LDw4lozZo1Q4YMMfsznSpVRqdLawbUvH9C22ut+J/u/PRQDqNCALADmzZtOnPmTBM9uCzLkiQ5+AFzb2/vxMTEJnpwYdCgQd99991rr73Wv3//kpKSwsLC3Nzcp59+2rxPY0rB7wfwrTXKiih++h557hFpwwV5yQDBD8e9AcDGLVy4sKysrE2bNmoXYp+KiopOnz7dhEE4c+bM5s2b+/v7Hz16NDU11dnZOT4+3s3NzYzPcapUGZMu/TCA7+DBysuJiD7pzTcTKC1f6bba+HU/Ia4lRoYAYNsmT57cdJ/UDu7gwYNmH57VJURHR5t+0uv1TbRw3lVDP0Ty7ZvflHb/CeMj/RWBo6cypcdasw968U6YPgkAABZnifDxd2G3pKDJoAAWqWdHHhcuV1GfdeLJUsddUwgAAGpRfxTW3Il+HMi/2IHrmyQuOIb5pAAAYFHqB6HJhBBuW7zw7Sk5cYt0DR0oAQDAUqwlCImoXXO2e5jQ1p26rRa3FeAwKQAAWIIVBSEROfM0N4xf3I//W4b02j7JiAOlAADQxKwrCE0GBbD9I4Sj15R+yeKZMgwNAcBuZZcoJTX32Ke0hk6U4JOwCVljEBKRj46SBwvj2nC914lLT2NgCAD2Kec6JaSJd8nC0hqKTxNzrluwJsdjpUFIRIxoekcufagw97A8IVMqN6pdEACAuQ1pyV7vwg9JFW87SbC0hhLSxFdCuaGNverIJ598cvXq1fsq8b5lZ2f/9NNPt73pzJkz3377rWXLuQ3rDUKTjh7s1+GCh5Y6rRJ3X8HBAQCwN7GBbHY3Pn7jrVloSsF/hnLDWzX+g3rx4sWFhYX3U97YsWMPHDhwP49w5syZ9evX3/ams2fP/vjjj/fz4GZh7UFIRDqBFvTmP+zFjdgkzjkgyUhDALAvf81Cs6TgLSorK8vKym7ZWFZWZmrPfifHjx+/fr1hR2YVRSkqKhJF0fTr0KFDlyxZUnvr9evXDQbrWiRnA0FoMrI1lzVCyLioDEoR8ysQhgBgV+pmodlTUFGUV199tU2bNp07d+7Xr19BQQER7d27t127dj169NDr9R999BERrVixYsyYMbX38vb2fv3110+fPj158uQePXosX768Ps917NixDh06hIeH6/X6L7/8kohWr149YsSII0eOdO7c+YUXXujatauPj897771nlj/NLGypJ3hLV5YxVPj0Nzlsrfjfvvww831RAgAwu4xLSpXYsLvEBHJ9k0QiGt2GaTi24XwDvvTrBIrU3/5U4vLlyzMyMnJycnQ63ZQpU1599dXFixePGTPmrbfeeuqpp86dO9e9e/ewsLCamppyU2MEIiIqLi5++eWXk5KSPv3004iIiPrUoCjKqFGjZs+ePW7cuNOnT3fv3j0yMtJoNFZWVoqieOTIkZkzZ3722WenT5/u1KnT5MmT6//XNSlbCkIi4hhN78j19mVPZEgrzyr/DeddbOwvAABHseSUfLmqYYevRJkuVigK0a7LtO+q1KD7+jizSD1/25vWrFnzzDPPuLi4ENGLL744YMCAY8eOVVVVTZgwgYhat26dmJi4fv36jh071vO5KisrT548SUQajaZDhw6127Ozs69duzZ27Fgiatu2bWxsbGpqqr+/v+lWT0/PJ554wnRTQEBATk5Og/7ApmOTMdKzBTvwmDBtp9Rjjbh8IN/ZE12cAMDq/K//7WPpTkxHRJcM4LU8e/uAlDxY8NSap5LLly97eHiYfvbw8CgpKSkqKvLw8GDsxoenl5dXUVFR3bsoinKXc4e5ubn/+Mc/iMjb23v16tW1269cudK8efPah/X09CwuLq4NQnd399o9tVptTc29VlBaiq0eXXTX0PcD+Le6c7Ep4oJjmEADALat7nnBO80jbTRFUY4ePWr6+ciRIw8++GBISEheXl5JSYlp46FDh0JCQjw8PGrj8NSpU6Yg1Gg0knTr2LRdu3bbt2/fvn173RQkIlmWz58/X1paavr18OHDbdu2Nc/f0JRsckRYa1QwF+bN/pYpbTgvf4dm9wBgm/46OyY2kBHx8RtFc40Lly5dGhIS4uPj89JLL82YMSM4OHjo0KHjx4+fMWPGzp079+/f//333xPRiRMnPv300+Dg4EWLFpkGdqGhoV9++WVBQUFYWFhISMg9n8jJyenZZ5+dNm3ali1b8vPzExMTk5KSzPAHNCVbHRHWau3Gtg4VevuybquNKQ05sQwAYA3uNEfUvOPCt9566/jx40uXLn3zzTeff/55Ilq2bFlUVNQ333xTUVGRlZXl5eXl5eWVmpp68ODBtWvXzps3b+bMmTqdbsGCBZGRkUePHr18+XJ9nqh169bPPvvs4sWLr1+/vnPnTp1O165du7Fjx/r6+v7973+v3e2ZZ54JDAw0wx9mDrY9IjQROJrTjY/QcxMypZFodg8AtuPuKyXMOC584IEH5s2bV3eLRqOZMWPGLbv17t27d+/epp9rVzhMnTq1Qc81ePDgwYMH1/4aGhoaGhpKRDNnzqzd+NJLLxHRqVOnGvTITcR+EiNSz46i2T0A2JSVZ+V/dbrbesHYQDarC7/qHC653ITsYURYy9TsfskpuW+S+EYXfnpH+4l5ALBLzzx874+p+CBG1Mi58fPnzw8MDNy2bVvj7t5QkZGRhw4dqv/+nTp1euONN5qunnqyqyA0mRDChbVgT2RI2y8ri8J5c80/BgCwOYMGDVK7hLvx8fHx8fFRuwo7OjRaF5rdAwBAPdlnENIfze4Xodk9AADcld0GoUkMmt0DAMBd2XkQEprdAwDAXdl/ENIfze5TY4W3D6LZPQAA3MQOZ43eSTdvdugx4bV9UqdV4tJIvrcPLtUNAGZz4cKF3377Te0q7FNT96lwoCCkP5rdR+jlEZvE59txs7vyHNIQAO7b+PHjZ82atXjx4iZ6fEVRals6OKZ27do13YM7VhCajGzNhbVgT2ZI2wvEJRF8gKtDv70A4P4NHz58+PDhTfTgoigajUadDl0FmopDnCP8q5auLH2oEBPAha0V1+ViBg0AgONy0CAkIp7RzM7cmkHCy3vlCZlSpah2QQAAoAbHDUITU7N7hajHGvHwNSw0BABwOI4ehPRHs/vXOnPRG9DsHgDA4SAIb5gQwu0bLvx8Vo5NEQuq1K4GAAAsBUH4JzS7BwBwQAjCm5ia3S+NFCbvkKbvlmownxQAwN4hCG8jUs8OPiacK0ezewAA+4cgvD1vZ1o7iH+xA9c3SVxwDANDAAC7hSC8mwkh3LZ44dtTcuIW6ZpB7WoAAKAJIAjvwdTsPsCFuq0Wt6PZPQCA3UEQ3pszTwt684v68U+g2T0AgN1BENZX3Wb3OdcxNAQAsBMIwgaobXb/6Fo0uwcAsBMIwoYxNbtPQbN7AAB7gSBsjO7e7NBjgoeWOq0Sd1/BYVIAABuGIGwkU7P7D3txIzaJcw5IuFY3AICNQhDel5GtuawRQsZFZVCKmF+BMAQAsD0IwvuFZvcAADYNQWgGaHYPAGC7EIRmY2p2LxOFodk9AIDtQBCak7uGfhjAz0SzewAA24EgNL/aZvdxqWh2DwBg7RCETcLU7P5RHzS7BwCwdgjCpoJm9wAANgFB2LTQ7B4AwMohCJtc3Wb3r+6VZu2TbtmhpIYmZEoihowAAGoQ1C7AUUwI4cJasCcypOtGKjVKn/fhTdtLa2hYmvhKKCfgOwkAgBrw6Ws5pmb3Q1uy70/JY9MlIiqtoYQ08ZVQblgr/EMAAKgDI0KLMjW7HxLEJW4W+5bwrhr51U5IQQAANeEjWAWDA9i+EZqTpVxOmRIVgH8CAAA14VNYBaU1NHm7+HUf0cuZ2v5svIJF9wAA6rlbEJaWlh44cKCystJi1TiCP88LtlQOPca3cWMPrTBiZQUAgFruGITffvttZGRkbGxsdna2JQuyb7fMjmFEOxKEsBas62px12VkIQCACu4YhBMnTjxw4EBoaKglq7F71wzKrC583dkxjCgtTpjSjnt8k5ichywEALA0nCO0qGA3FhvIbtnIiOb34tcNFibvEP97AuvqAQAsSiCiNWvWbNiwoe7W+Pj4YcOGqVSSgwprwXYkCHGp0rly5b0w/ta0BACApiEQUc+ePVu2bFl3q16vV6keh9bGje0aJgxLEydtlRb34zUYrgMAND1hy5YtUVFR/v7+alcCREReWkqLE8ami0M2ir9EC+4atQsCALB33MSJE0VR/OsN27dv79GjR2lp6eTJk//+979bvjKH5SrQmkFCW3fWL0nMr8D0GQCApsX8/f2XLFkSFRXVuPt369atpqbG3d297sbExMTnnnvutvtXVFS4uLgwhlNgVFVVpdVqOe6OB0A/OiF8c4ZfFVET4mbncVhVVaXRaAQBF/wjg8HAcZxGg0MBZDQaZVnWarVqF6I+URSNRqNOp1O7EPXJsmwwGBr0Uri4uPA8f/d9hJCQkJycnEYHoU6nmzRpUpcuXepuDAwMdHNzu+3+jDFXV1cEIRHxPO/s7HyXIJzdk1p5yEMy2MooIdzPnl8xQRAQhCZOTk4IQhMEYS0EYS1ZljUajYuLi3kfVmjWrFlpaWnj7y8InTp16tevnxlrglpPhXABLmzkZvGLvnxiMCbPAACYH1dcXOzr66t2GXBH0QFsQ6wwfbf8xXEsMQQAMD/uyJEjnTp1UrsMuJse3mxHAv/Jb/L03ZKdny0EALA4rnPnzp07d1a7DLiHYDe2a5iwv1CZuFUyYmQIAGA+3MqVK9WuAerFU0ubhwiVIsWlimVGtasBALAXnJ+fn9o1QH058/TjQP6hB1g4lhgCAJgJJiLaGJ7RF335Zx7i+iVL2SXIQgCA+4UgtEnTO3JzunED1ovbC5CFAAD3BUFoqyaEcEsjhcc3iz/nYPIMAEDjIQhtWJQ/2zJE+Nev8udYYggA0FgIQtsW6sm2DuU/Oy5P3y3JOEoKANBwCEKb19qN7UwQDhYpT22VajAyBABoIAShPfDU0qY4oUamuFSxtEbtagAAbAqC0E5oeVoWybdrzsKTxAtYYggAUG8IQvvBM/qsD//sw1zvddKRa8hCAIB6QRDam+kduffCuJgUcRuWGAIA1AOC0A492ZZbFimM2iL+hCWGAAD3giC0TwNNSwz3yh8cQRYCANwNgtBudfRgu4fxP5zGEkMAgLtBENqzAFeWOVQ4VKSMTpeqJbWrAQCwSghCO+ehpbQ4QWA0BEsMAQBuB0Fo/7Q8LR/I92jBwpPE81hiCABwMwShQ2BE7/fkn32Y67NOOowlhgAAdSAIHcj0jtxHj3KxKeLWS8hCAIAbEISOJTGYWzZQGJ0u/ngGyyoAAIgQhA4oUs+2DBFm7pPfxxJDAAAEoWMyLTFchiWGAAAIQofl78Iy44XD15RRW7DEEAAcGoLQcTV3oo2xghNPURvEIoPa1QAAqARB6NBMXQzDfVn/JDGvHAdJAcARIQgdHSOa15Of/AjXJ0k6VIQsBACHgyAEIqLpHbkFj3IxKWJaPrIQABwLghBueDyY+zlKGJ8pLsMSQwBwJAhC+NMAPUsfIszaJ885gImkAOAoEIRwkw4ebNcwfm2ugiWGAOAgEIRwK38Xti1eyC5RErdIVaLa1QAANDEEIdyGm4aSBgs6ngZuEAur1a4GAKApIQjh9pw4+iGSHxzIIpKxxBAA7BmCEO6IEc3pxk9px/VJkg5iiSEA2CkEIdzDix24T3tzg1PE1AvIQgCwQwhCuLfHWnNrY4Snt4lLT2OJIQDYGwQh1EtvH7Z5iPB6FpYYAoC9QRBCfbVvznYPE9blKv/YJUk4SgoA9gJBCA2gd6Ft8cLpMiVxM5YYAoCdQBBCwzTT0NoYwVWDJYYAYCcQhNBgThx9P4AfHMh6rxNPl+EgKQDYNgQhNIZpieGLHbgB66UDhchCALBhCEJovH904D7rw8WmiinnkYUAYKsQhHBfRrTi1sUIT28Tv/odSwwBwCYhCOF+PerD0ocK7x7CEkMAsEkIQjCDds3Z7mFCUp7y7HZJxMgQAGwKghDMw09HW4cKFyuVxC1SJZYYAoDtQBCC2TTT0LpBQgtnGrhBvIolhgBgIxCEYE4CR4v68bGBrPc68VQpppICgA1AEIKZmZYY/l9nLnKDtB9LDAHA6iEIoUk88zD3eR8uLlXcgCWGAGDdEITQVIa34pJihGe2iYuzMZEUAKwXghCaUC8ftiNB+OAolhgCgPVCEELTetCdbYsXkvOUp7dhiSEAWCMEITQ5Px1tjRcuVylDN4rXjWpXAwBwMwQhWIKrQGsHCa3cWNQG8UqV2tUAANSBIAQLEThaGM4Pacn6JGGJIQBYEQQhWI5pieGsLly/ZHHXZWQhAFgFBCFY2tMPcYvC+YQ0cVzGrVNJDRJN3SlV4VKlAGBBCEJQwbBW3IZYISlXjt7wZ+jVyDQmXerowXSCiqUBgMNBEII6erVghx8X9hcqPdeKRFQj0+gtUkwAm9oe70kAsCh89wbVPOjGTozSdP7F2COZf9CdYlsiBQFABQhCUJOfjn5L1LT5yXj6utJcq3g5yzEBnIdW7bIAwJHc8Qv4//73v549e3bt2nXUqFFXr161ZE3gOGpkena79HYX+ckH6Uq18stZpc1Pxh5rxDkHpP2FCuaVAoAF3DEIH3744YyMjIMHD4aGhr711luWrAkcRO15wSkPyf/tyx5+gAW60qW/aeaG8dUSPZkp6ZcaJ2RKK87KuB4NADSdOwZhnz59XF1diSgkJKSiosKCJYFDuGV2DCP6uDcvKTRrnxQdwOaG8ScShV3DhO7ebFG2HLjMOChFnHdYzi7BKBEAzEwgoitXruTl5dXd6uPjExQURETFxcXvvPPOjz/+qE51YL8uVymJwezJtn9+FTNl4SfH5CqRTCso2rix6R3Z9I5chUjpF+XkPCU6RdZyFB/EEoK4/nrmhLk1AHDfBCI6evTozz//XHfrgAEDgoKCysrKEhIS3nnnndDQUJXKA7vV0pU92ZbdspERTe94m3BzFSghiEsIIiL6rVhJzlPmHJCOXFMi9CwhiIsPYv4utz4UAEA9CYWFhVFRUVFRUbfcUFlZ+dhjj02dOvWxxx5TpTKA2+rgwTp4sJmducJqyrgkJ+UqM3+V/F1YQisW35Lr48s4ZCIANISQkJCwe/fuv94wderUvLy8rVu3bt269eGHH3755ZctXxzAXXg706hgblQwSQq/+7KSfF6esUfKK1cG6Ln4IDasFdfcSe0SAcAWME9Pz6SkpD59+txyw+nTp0tLS00/u7m5PfTQQ7e9f+/evSMjIx9++OG6Gzt06NC1a9fb7l9eXu7q6soYvrRTZWWls7Mzx+E0F1VVVWk0GkEww6rWc9eVTRdpy0XadFHp0JziW7LoANbN6/4f2EIMBgPHcRqNRu1C1Gc0GmVZ1mqxqpREUTQajTqdTu1C1CfLcnV1tYuLS/3vwnHcPRNH6NWrV1ZW1l+DsG3btvV5DqPRuHfv3rNnz9bdqChK+/btb7u/wWAQBAFBSEQGg4ExhiAkIoPBIMuyJN16De5G0DvRhNY0oTVVibSnkNuQzx7fwvGMBvrJkX5yjF5pZt0RYwpCWZbVLkR9piBUuwqrYApCfFYQkSzLBoOB5/n638XZ2fneQejt7V1QUNDoslxdXWfPnh0REVHP/WVZdnFxQRCaYERowhgz14iwlgvRUHca2oY+J8q5riTlKt+dlaftVXr6sPiW3GOtWVAza3wT8jyPEaEJRoS1MCKsJcsyx3ENGhHWh1BZWWlaLwhgr2qXYZQbKeOSnJynzFsn6XhmWoYRoWcafBsBcGDCuXPnhg8frnYZAJbQTHNjGcaXCn+wSNmcr8w5IGWXKAP9uegAlhDE6c38RRMAbICQnZ0dHx+vdhkAFsUx6u7NunuzmZ25q9WUcl5OzlNe22ds48aiA1h8S66vHw7fAzgKYf78+R4eHmqXAaCaFs40IYSbEEKizO+5cmMZxvkKZXAAl9CKxQRwD2AZBoBdEyZPnqx2DQBWQeAo3I+F+/EURjnXlc35ypJT8qStUhcvlhDEJbRi7ZtjlAhgh9CPEOA22rixyY+wyY9wlSLtuqwk5cmxKbKGo+gAFh/EYgI4bQPmbwOAVUMQAtyNi0DRASw6gF/Q+8ZlTj/5TX4yQzItwxgZzFq6YpgIYNsQhAD1VXuZ0yIDpV+UN+crc9dInlqW0IpF+3MD9EzAMgwAG4QgBGgwL+2Ny5x+0Zc/VKQk5cmv7ZNyy5VIPRcdwIa34nyx9BnAdiAIARqPv7EMg5/TjS5X0cYLcnKeMvNX44PuN1brd/PGMgwAa4cgBDAPX92NZRjVEr+jQNl8UR6fKV0zKDEBXEIrFhvIueG6aQBWCec0AMzMmafoADY3jD+eKOwaJoT7sSWn5MBlxkEp4rzDcnaJcsv+lyrpjSzplq2SQv/aK5XUWKxqAMeFESFAE6pdhlEhUvpFOTlPiU6RtX8swxgcyDlxpHchhWjKDmlh+I01GZJCT2+T2roztFQEsACMCAEswVWghCBuYTh/YZywLoZv48bmHZa9vzcmpImLsuWp7bgWzjRlh6TUScE3u+K/J4AlYEQIYGm1yzAKqynjkpyUq8z8VfJ3YQ9o6fFMzsOJHm6OFASwHAQhgGq8nW8swxBlfsdlZf15+bPfFFmhxGBlySk5JpDzwzIMgKaHIARQn8BRPz/2v5P0akeluIYdL1M259NLe4zNnZjpbOKgAM4ZF3UDaBoIQgD11Z4XfLW9wnFszmF2tZquPKk5VKRszlc++U0enymFtWDR/lx0AMPaRADzQhACqKzu7BiDgYjo3R7861nS8zulheG8qWlihUi7LytJeXLiFrlaVAYFcKbrunlo1a4ewPYhCAFUVlhNPVuwae1vmh3zbg9+/lG5rIZM3RBd61z729QiakWOMmXHjU7C0f5chJ5pML0GoFEQhAAq89XRLSlo8nLo7ZPtj7WJJMr84Ws3rnR6slSJ0LOEIG5wIGvVDIdOARoAQQhgqwTuzyudXq2mzEvy5nzl7YM3FuxHB7DBgZw7rusGcC8IQgB70OKPlRhElHNdScpVFmXLk7ZKXbxYQhCm2ADcDYIQwN60cWPTO7LpHblKkXZdVjZflKfskM9XKBF+nGkxhr8LMhHgTwhCALvl8scUGwqjS5W0KV9OzlNe2yfpdTeaCffzY1osTwSHhyAEcAh6lxtdoiSFNy1PnHdEStyi1C5P7O6NYSI4KAQhgGP5o5kwm9mZKzJQ+kV5c74ycrMssBtTbAYFcOh6AQ4FQQjguLy0N02xMS1PnLzd+KD7jeWJA/RMwPJEsHcIQgAgqrM8sVridxQomy/Kr+2TcsuVSD0XHcDiWrKWrjh2CvYJQQgAN3Hm/5xic6WKthbIm/OVOQckHY8rgIN9QhACwB356G4cO/1S4Q/iCuBgpxCEAHBvXJ0pNqYrgG++KE/IlIsMSn8/Lj6IxQdxnrgCONgmBCEANEztFcDnht2YYpOcp7y0xxj8xxXA++uZE6bYgO1AEAJA49VOsaldnjjngHTkmtLLh8W35Ia3Yq3dcOgUrB2CEADMoO7yxMJqyrgkb85XPjwqO/1xBfCYAO4BLE8Eq4QgBAAz83a+dXniklPypK3SI81ZfBBLCOK6ejHuDgPFciM1+0vHjNtuBDAXHMgHgCbUxo1NfoRLihGuPqmZG8ZXSzRlh+S31Dh6i3VTqe8AABg7SURBVLQoW86vUOruXFJDEevFM2U3bcyvUAasFy9VWrZucCQYEQKAJejqXAG8oIrSLty4AriH041hYj8/1tyJfhjAj02Xlg/kW+mIiC5UKCM2SZ/14fUuav8BYL8QhABgaX66u10BfGZnbly69H1/cmI0eqv0WR/+UR/MuIEmhCAEANXUnWJTbKAtF+WNF5QvTiiiTH2SyUPLLYvkeyEFoYkhCAHAKnhoKTGYSwwmIsq8JI9Jl/Iq6Gq1QoQghKaFyTIAYF0uVCiv7JV/iaQPesgjN0tpF2S1KwI7hxEhAFiR2tkx3T3EMC8yEJeQJqXFUYQe39qhqeC9BQDWoqSGRmySPq8zO2ZmJ/6fodzgVOlIsbqlgT3DiBAArEVzJ/olmm/V7KaTgv8J46tlGp8hZg4VPHBdb2gCGBECgBW5JQVN5vfiYwLYkI1iudHyFYH9QxACgA14vxcf6ske2ywaJLVLAbuDIAQAG8CIvuzLe2nZ2HRJxDRSMCsEIQDYBp7R9wP4Gll5epskK/feH6CeEIQAYDM0HK2MEnLLlem7cYQUzAZBCAC2RCdQ0mBh9xVlzgFkIZgHghAAbIy7hlJihZ9zlA+P4mwhmAGCEABsTwtnSovjvzguf/U7shDuFxbUA4BNCnRlaXH8gPWSu4ZGt8F3emg8BCEA2Kq27ixlMB+TIrppWFxLNKmARsLXKACwYaGebM0gYeI2cVsBVlRAIyEIAcC29fJhyyOFUVvE/YXIQmgMBCEA2LyB/mxROJ+QJp4oQRZCg+EcIQDYg+GtuOtGGpwibYvnW7vhfCE0AIIQAOzEk2250hoalCJtTxD8dGpXA7YDh0YBwH5Ma8+ND+FiUsRrBrVLAduBIAQAuzK7KxcbiOaF0AAIQgCwN/N68p092YhNYjUuRwr1gCAEAHtjal7YQofmhVAvCEIAsEMcoyURvCgrk9C8EO4FQQgA9knD0YooIa9ceRHNC+GuEIQAYLdMzQv3XFH+vR9ZCHd0x3WEWVlZ6enplZWVERERkZGRlqwJAMBc3DWUGitEJIvNNPK/OuGrP9zGHd8Wp06datWqVZcuXV599dWVK1dasiYAADPydqa0OP7LE/LibMycgdu444hw3Lhxph8uXLhw4MCBxMRES5UEAGBmAa5sUxwfsV5yd6IxaF4IN7vbJdaWLVu2Z8+e7Ozsr7/+2mIFAQA0hQfrNC8cguaFUAdHRLNnz+5+szVr1hBRaGjooEGDRFHcsWOH2nUCANyvUE+2NkaYhOaFcDOBiN5+++233377r7eFhoaGhoa6u7t/+OGHtUdKAQBsV88WN5oXbhgsdPfGuBCIiDjT4O+vkpOT8/Pzc3Nzv/rqq549e1q4LACAJjLQny0O5xPSxONoXghERCS88sorw4YN47hbzx7//vvvCxcuFAQhMjLy+eefV6U4AICmMKwVV2ak2BRpazwfjOaFDo81a9Zs+/btXbp0adz9O3funJ2dfUuOTp06dc6cObfdv6KiwsXFhTG886iqqkqr1f71K4gDqqqq0mg0goDumGQwGDiO02g0aheiPqPRKMuyVqttuqdYfIr//KSQOrDGT2fVQ0NRFI1Go06HFosky7LBYGjQS+Hi4sLz/N33EfR6/YULFxodhM2bN09KSurXr1/djXf5UGOMubq6IgiJiOd5Z2dnBCERCYKAIDRxcnJCEJpYIAhf7kaVTH58B5c5VPBswue5XwjCWrIsazQaFxcX8z6s4OTkZDDcVwtLrVaLfyEAsEVvdOWuG5UhG8XNcUIzfP1wVNyVK1f0er3aZQAAqGNuT76LJxuO5oUOjDMYDJ06dVK7DAAAdTCiL/ryPjo2Bs0LHRU3derUZs2aqV0GAIBqTM0LJTQvdFTcu+++q3YNAAAqMzUvPF+h/APNCx0Ph1mLAABEpBNoXYzw6xXlTTQvdDBIQQCAG9w1lBIrrDqrvH8EZwsdCIIQAOBPpuaFC0/Ii9C80GFgCTMAwE0CXFlaHD9gveSuobEPYrRg/xCEAAC3etCdpcTygzaIbk5sKJoX2jt82QEAuI2OHmxtjPD0NnHrJayosHMIQgCA2+vZgv04UBidLmYVIgvtGYIQAOCOIvXsq358wkY0L7RnCEIAgLtJCOI+7MXHpkhnryML7RMmywAA3MPf2nJlRhqUIm2PF/RmbgEE6sOIEADg3p5vxz39EBeTIhbdV9s6sEYIQgCAepnVhRsaxIakiteNapcCZoUgBACor/fC+K5ebHgamhfaFQQhAEB9mZoX+rmwMemSEZdgsxcIQgCABuAYfRfBywqaF9oPBCEAQMNoOPp5oJBfobywC0dI7QGCEACgwUzNC7MKldezkIU2D0EIANAYbhpKiRXW5irzDuNsoW1DEAIANJKXltLi+EXZ8kI0L7RluLIMAEDj+buwTUP4iGTJXUPj0LzQNiEIAQDuSxs3lhrLR28Q3TQsPgjNC20Pvr8AANyvDh5sXYzwzHYxE80LbRCCEADADMJasFXRwrh0cd9VZKGNQRACAJhHX1/23QBhxCbpt2JkoS1BEAIAmE1MAPukNxeXKuWgeaHtwGQZAABzejyYKzNSTIq0LZ73d8HcGRuAESEAgJlNeoh7oT0XkyKheaFNQBACAJjfjI5cQhCLQ/NCW4AgBABoEu+F8X190bzQBiAIAQCayvxH+dZubPQWNC+0aghCAICmwogW9+OdeBqXIUmYRmqtEIQAAE2IZ7Qski83Ks9tRxRaKQQhAEDTcuJoVbRwslR5eQ/OFlojBCEAQJNzESh5sJB5SZmL5oXWB0EIAGAJzZ1oY6zw7Un542PIQuuCIAQAsBAfHaXF8R8fk/93ElloRXCJNQAAywlqxtLi+Mj1kruGHg/GUMQq4J8BAMCiHnqArR/MT9slpeVjGqlVQBACAFhaFy+2Klp4MkPceRlZqD4EIQCACvr4sh8ihZGbxUNFyEKVIQgBANQRE8C+7MsP2Sj+XoosVBMmywAAqGZka66shganSNvi+aBmaF6oDowIAQDUNPEh7qWO3KAU6XKV2qU4KgQhAIDKpnfkEoNZbKpYUqN2KQ4JQQgAoL53e/CRejZ0o1ghql2K40EQAgBYhf/3KP/IA2zkJtGAS3NbFoIQAMAqMKJF/Xh3J/YEmhdaFoIQAMBa8IyWRvIVovIsmhdaEIIQAMCKmJoXni5Tpu68zRHSYgMhIM0OQQgAYF1cBFo5UFhySp6266YsvFipxKSIZ8oQhWaGIAQAsDq+LrQjQfjmd/mlP5raX66iEZukD3vxbd2x7t7McGUZAABr1NWLZcbz/ZMlF55NbstGbRM/6MlH6JGC5ocgBACwUr1acCmDKTZF+vYkvyKa7+OLFGwSODQKAGC9HmnOWjejqwZ6+6BcWK12NXYKQQgAYKVM5wX/25cdiBcPXZPbrTQuOSWrXZQdQhACAFijy1WUkCZ+0JPv70cPutG+4YKXM/v3ATkhTbxYiYmj5oQgBACwOlUixW8UP370z9kxLV1ZWizvoaX2zVnX1eKibAwNzQaTZQAArI5OoOTBgq/upo1BzdjGWMFLS0+05SZtldbmygvD+UBXzKC5XxgRAgBYo1tS0KSFM3GMOnuyvcOF/n5ct9XiomwZx0nvE4IQAMD2aDia2ZnbMkRYlC0PSRXPVyANGw9BCABgq0I92Z5hwgA91x1Dw/uAIAQAsGECRzM7cxlDha9+lwcki6dxJdKGQxACANi8Dh5sV4IwpCXXa60477CMsWGD3C0IZVnu16/fwIEDLVYNAAA0jmlouHe4sOG8HLFePFmKMKyvuwXhggULevToYbFSAADgPrV1ZxlDhfFtufAkcd5hGe196+OOQXj27NnMzMy//e1vlqwGAADuE8do8iPcnuFC6gW5f7KYXYIwvAeBiJYtW1ZWVlZ36+jRo2fMmPHRRx8VFxerVBgAADReGzeWPlRYnC33Txb/Gcq/0onjsfL+Djgi4v8iLS3t8uXLK1as+O677/Ly8pYvX652nQAA0DCMaPIj3N7hwqZ8OTxJPIGh4R0IRqNxzJgxt2zNzc3led70s7Ozs7+/v8ULAwAAMwh2Y5uGCIuz5YhkcWp77vUuvAbLBW7GzZ8//69bW7VqNWrUqFGjRkVHR/v4+ERERFi+MgAAMAvT0HDfCGHXZSVsjXiwCEPDm7A2bdqcOXOm0ffv2bNnQEBAQEBA3Y0REREJCQm33b+8vNzV1ZUxHKumyspKZ2dnjsN3M6qqqtJoNIKAS8CTwWDgOE6j0ahdiPqMRqMsy1qtVu1C1CeKotFo1Olud+3RhluVx2b8yiY8qMzurDjZ2sePLMvV1dUuLi71v4uTk9M9P2aFs2fPVldXOzs7N64sjuP8/f3btGlTd2OLFi1qj6zewnQOEkFIf7wUCEKqc5Za7ULUZ3pL4KUgIlmWGWN4KYhIURRZls31UowKpnBfemEP65tCX/Wlbl5meVQLMb0lGvRS1CduBEVRioqKbhnS1Z9Wqx09enT9j51qNBqNRoMgpD9eCgQhEYmiiBGhiSzLGBHWkmUZLwX98VFuxpci6AFaN5hWnJWHb5EmPsS91Y3X2sj3DVmWJUky+7uCIyJPT0/zPigAAFi5UcHcoZGak6XUY42476pDnzXkAgICzHXoGQAAbIifjlZF87O7ccPSxNf2SQZJ7YJUwk2cOFHtGgAAQDWjgrnDIzVnyqjjL+K2AkccGnKzZs1SuwYAAFCTj45WRPH/CeNGbxGn7JAqRLULsiyuQfNQAQDAXo0K5k4kaoio8yox85IDDQ0xZREAAG7w0NLCcH5Bb358pjRlh1RuVLsgi0AQAgDATYa2ZEcfF4io0yox/aL9Dw0RhAAAcKvmTrQwnP+8Dz9pmzRlh3TdroeGCEIAALi9uJbsyMgbQ8PN+XY7NEQQAgDAHT3gRAvD+f/25Z/dLk3IlIoNahfUBBCEAABwD4MD2ZHHBQ8tha4S1+bKapdjZghCAAC4N3cNLejNL4/k//WrPHqLVGRHQ0MEIQAA1Fc/P3boMaGNO3X6RVx9zk6GhghCAABoABeB5obxP0fx/7dPHr1FKqxWu6D7hiAEAIAG6+vLDpqGhquMK8/a9tAQQQgAAI2hE2huGP9LtPDmfnn0FumqzQ4NEYQAANB4vX1uDA3brzQuyrbJoSGCEAAA7oszT3PD+E1xwpcn5IQ0Mb/CxpbeIwgBAMAMunixX4cL4b5ctzWibQ0NEYQAAGAeGo5mduY2xQkLs+UhG8ULNjI0RBACAIA5dfJke4YJEX5c9zXiomzZ+sMQQQgAAGZmGhpujhMWZ8txqWJeuVWnIYIQAACaRKgn2z1MiNRzPdaIC47JVjs2RBACAEBTETia2ZnLGCosOyMPWC+eLrPGMEQQAgBA0+rgwXYmCI+35nqvE+cdtrqhIYIQAACanMDR9I7cnmFCygW5f7J4stSKwhBBCAAAFvKgO8sYKkwI4fomifMOy5J1pCGCEAAALIcRTX6E2ztc2HhB7pckZpeoH4YIQgAAsLQ2bmzLUGHiQ1z/ZPWHhghCAABQgWlo+OtwYfNFuW+SeFy9oSGCEAAAVNPajaXFCU8/xPVLEl/bJxnVuEYpghAAANRkGhoefVw4Xkxha8QDhZYeGiIIAQBAff4ubF0M/3pXbshG8bV9Uo0Fh4YIQgAAsBajgrlDIzW/l1CPNWKWpYaGCEIAALAifjpaPYh/sysXv1F8bZ9kkGjf1dsk4u+lSpnRPM+IIAQAAKszKpg7PFJzqpS6rRb/uVe+pdPvb8XKuHTpapV5howIQgAAsEa+Ovolmp/TnTtZKs87LH987EYWZpcoT2ZKSyP5B92ZWZ4IQQgAANZrVDB3ZKSmqzd7I0v+5175ZBkblyEti+TbNTdPChKRYK4HAgAAaAo+OloZxS87Iz+9TfrlHL9lCB/ygNlSkDAiBAAAm9DZk4W4sTbNKOOSmWeTYkQIAADWznRe8MeBLMBJnLibF2V5anuzDeQwIgQAAKuWXaLUnhd04ujnKD4tX/niuNmW3CMIAQDAeoky/X2n9PPAP2fHOHH000A+9YJysMg8x0hxaBQAAKyXwNGWIQJ/8+QYLU+rB/G8mWbMWHpEeP78+fLycgs/qXXKz88vKytTuwqrUFBQUFxcrHYVVuHq1auFhYVqV2EVCgsLr169qnYVVqG4uLigoEDtKtRUG3hlZWX5+fm3bLx/lg7CF198cceOHRZ+Uus0a9aslJQUtauwCu++++6KFSvUrsIqLFiw4JtvvlG7CqvwzTffLFiwQO0qrMKKFSveffddtauwCikpKbNmzTL7w+IcIQAAODQEIQAAOLT7nSwjiuLhw4frv39xcfHRo0ddXFzu83ntwNWrV48fP75161a1C1FfQUHBqVOn8FIQ0fnz50tKSvBSENG5c+euX7+Ol4KITp06VVBQgJeCiI4fP3716tUGvRQ9evRwdXW9+z5MUe5r+uny5cs/++wzQahvoJaXlzs7O9d/fztWUVHh5OSk0WjULkR9VVVVPM87OTmpXYj6qqurGWNarVbtQtRnMBgURXF2dla7EPXV1NRIkqTT6dQuRH1Go7GmpuaewVbX119/3bZt27vvc79BCAAAYNNwjhAAABwaghAAABwaghAAABwaghAAAByaRWdvFhcX5+TkBAYG+vr6WvJ5rY0sy1lZWWfPntXr9eHh4RznuF9HysrKsrKyrly54uXlFRERgYmjRJSTk1NRUREaGqp2IaopLCzMzc2t/bVdu3YOvuDqxIkThw8fbt68ed++fd3c3NQuRx3nzp0rKiqq/ZXn+S5dupjrwS0XhMOGDdu4caMgCO+8885LL71ksee1Qn369KmsrOzYseORI0ecnZ0zMjIc9s09bdq0S5cu+fv7//7770VFRdu3b9fr9WoXpabz58+HhYU1a9asbhI4mlWrVr355pudOnUy/frll1/ec/q7HXvllVeWLl3av3//0tLSrKysN954Q+2K1LFy5cqNGzeafj5z5oxOp/vtt9/M9uiKpZw8edJgMMTFxc2fP99iT2qdjh8/bvrBaDR26NDh008/VbceayDLcr9+/ebOnat2ISqLi4t7+eWXg4KC1C5ETQsXLkxMTFS7CquwevXqwMDAK1euqF2Idenfv/+8efPM+ICWOygXEhKCA18m7dq1M/0gCIK/v391dbW69VgJURQ9PDzUrkJN3377rY+PT0xMjNqFqK+srCw9Pf3YsWOybLbmq7Zo6dKlkydPrqmp2b59e0lJidrlWIWcnJzdu3ePHz/ejI/puGenrMHu3bv37ds3atQotQtRU3Jy8siRI9u3b9+lS5dJkyapXY5qCgoK/vOf/3zwwQdqF2IVzp8/P2/evNjY2L59+9Y9M+Rozpw5s23btvj4+Pfffz8kJGTTpk1qV6S+r776asiQIeY9h4IgVM3p06dHjx79xRdftGrVSu1a1NSpU6fnnntu4sSJq1at2rt3r9rlqGbatGlz5sxp0aKF2oWob9KkScePH9+4ceOZM2eaNWv273//W+2KVFNZWVlZWZmVlZWUlPTuu+9OmzZN7YpUJorid9999/TTT5v3YXHNT3Xk5uZGR0fPnj173LhxateisqCgoKCgoLi4uPLy8o8++ig8PFztilSQnZ2dmprq7e29devWCxcuXLt2bcqUKe+9956np6fapamg9gK8Wq12zJgxX3/9tbr1qEiv1/fu3ZvneSKKioqaMmVKVVWVI190NDU1VZKkuLg48z4sglAFFy5ciIqKeuWVV5577jm1a7EiJSUlDbqWrj3x8fH56KOPTD/rdLo9e/Z0794d59SJ6PDhwwEBAWpXoZoBAwZkZ2ebfs7JyfH09HTkFCSib7755qmnnjJ7rwLLXXR7xYoV+/fvX7FiRevWrcPCwsaOHWvGVSC2pVu3bmVlZYmJiaZfH3300REjRqhbklpiYmL69Onj7e196NChn3/+OSMjo3v37moXpbKNGzdOnjzZkZdPPPPMM76+vnq9fv/+/b/88ktmZqbDvisuX77crVu3J598Mjg4+P3333/hhRdefvlltYtSzZUrV1q2bHno0KHa+YbmYrkgzMzMPHnyZO2vAwcOdNi1QT/88ENlZWXtr+3bt3fM44FEtGnTpt27d5eUlAQFBY0ePdrf31/titR34cKFbdu2PfHEE2oXopodO3akp6eXlJS0bNly9OjRjjwiJKL8/PzvvvuurKxs0KBBUVFRapejppycnD179jTFfw20YQIAAIeGWaMAAODQEIQAAODQEIQAAODQ/j9sh4o3lsJ7LQAAAABJRU5ErkJggg==", "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 }