{ "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.844911905950 NaN 1.97e-01 4.0 \n", " 2 -7.850326616771 -5.41e-03 2.93e-02 1.0 \n", " 3 -7.850617259203 -2.91e-04 2.81e-03 1.5 \n", " 4 -7.850646790884 -2.95e-05 1.26e-03 2.5 \n", " 5 -7.850647351232 -5.60e-07 5.91e-04 1.0 \n", " 6 -7.850647505110 -1.54e-07 4.22e-05 1.0 \n", " 7 -7.850647511556 -6.45e-09 5.46e-06 2.2 \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+gvaeTAAAgAElEQVR4nO3deVwU5/0H8OeZGY5dDhERueRUREFRARURRQVELiVeMTHGajRnJZdV+2tatTGJsa1JbGxijiaxtkaNiiAqhxA1gCheiRo1KuKBGBGQG+b4/bGGEEVdYNnZ3fm8/8hrdpid+T4Y9rPPzPPMUEmSCAAAgFIxchcAAAAgJwQhAAAoGoIQAAAUDUEIAACKhiAEAABFQxACAICiIQgBAEDREIQAAKBoCEIAAFA0BCEAACianoIwLS1t7Nix48ePz8/Pf9A2TU1Np0+f1n6fgiDoojRjhebLXYKcRFGUuwQ5oflylyAnURR1fmdQqod7jZaXlw8fPrywsLCmpmbcuHGnT582Nze/f7Pi4uKIiIji4mItd1tdXW1jY6PLQo0Kmq/k5tfU1FhZWVFK5S5EHrW1tSqVimEUekKrrq7O0tJSsc2vr683NzdnWVaH+9THrzIjIyM6Otre3t7d3d3f37+goEAPBwUAANCGPoKwtLS0V69emmVnZ+fS0lI9HBQAAEAb+ghCKyurxsZGzXJdXZ21tbUeDgoAAKCNTgVheXl5SkrK8uXLP/roo9br6+vr33777enTpy9btqy6utrf3//48eOaH508eXLAgAGdOSgAAIAOcZ158//+978tW7YQQkRRfO6551rWz507t7y8/MUXX/ziiy+mT5+enp7e1NS0dOnSysrKQYMGeXl5dfiI1c3EjCGW910l/bmB9LTs8F4BAEC5dDBq9JNPPvnqq68OHDigeVlSUtK3b99r1645ODjU1dU5OjoWFBT06dMnOzvb3Nx83LhxDxrtU1xcPHLkyLVr17ZeaW1tHRUV1fJyZ4m0/kdx63jWkv113OCnZ6V9pdJ/I5Q1hkrhwyYV3nyMGsWoUcU2v72jRrX5RXWqR9imI0eO+Pn5OTg4EELUanVwcHBBQUFAQEBcXNwj39vQ0PDVV1+1XtOrV6/w8PCWl9GOpKyGmZwh/i9cEBobzczMvrzApFyh/w0XGhp03hSD1tjYaGZmJncVslF48xsaGliWVWwQNjQ0UEoVmwQNDQ1Eu893k9TQ0CCKovZBqM2XBt0H4Y0bN3r06NHy0sHBQfthora2tikpKQ/f5tkAYm4uzspjvxyu+t8Vy7Tr4o4J3P0nS02eIAhqtVruKmSj8OaLoqhWqxUbhJIkKblHSLT7cDdVlFKdzyPUfRCqVKqWMaKEkIaGBp1/YP3Ol5EIicyx6Gkppk/gLJSXggAAoCu6/07h5uZWUlLS8rKkpKR37946P4ogktom8l2ZtLNE0XcbAgCATtJ9j3D06NGNjY3Z2dnjx48/evRocXFxTEyMbg/x2Vlxx2WxIKZx7QX1nG+FHZekT0azat03BQDggbZu3XrPmAb9EEWRUqrAE+Ph4eGLFi3qij13Kj327ds3f/786urq6upqHx+fmJiYDz/80MLC4r333psxY0ZISMiRI0dWr15ta2urq3IJIZ+dFbcVi99Ecs11ZEUQ66SmK48JwdulTePZQfaK+z8DAORy7tw5CwuLp556Su5CFOHgwYNHjhzpop13KghHjBiRmZnZ8tLKykqz8MQTT0RGRp4+fbpfv37Ozs6dKvC3dlwWd1wWt0VyFixpJoQQ8kJ/xpIl68+I49P5Pw1mkwMUegEZAPSvb9++iYmJclehCPX19ZcvX+6inXcqCNVqtbe3d5s/cnR0dHR07MzO2xTpwkx0Y+4ZHTPXl5nkwdxqkJ7IEXJLpU9Hsz0sdH5kAAAwTUbWf7I2I22OEe1hQfp1owWJXGAPErSd33+jy58tBQAApsHIgvDhzBiybCj7aTg7c5+QnC80YzwpAAA8ikkFoUakKz2WxP10RwpP4y9Wo2sIAAAPY4JBSAhxVJG0CdxMb2ZECv+/C+gYAgDAA5lmEBJCKCHJAczuGO4vR8XZuUItL3dBAABgkEw2CDWCHGjRZE7FkeAd/PFynCYFAJmlX5EeOZrvdKX01XmcytIfEw9CQoiNGfl4FLsiiJm4h191QhSRhgAgn2AHuuiQkFP6wE+iHyulJ3OE4J64PYj+mH4QakzzYgomcTtLxIl7+Bv1clcDAErlqCKp0dySwraz8MdKaWaOsHEsO8CuI0EYGRlZWVnZ6RoJIUSSpIqKis4/sLapqSk8PLypqekh2zz22GOt71Ctf0oJQkKIhzX9No4L7UWHbm/ecxUdQwCQx4OysJMpSAj54YcfeF43AyIaGxvt7e3r6zvbb2AYJiAg4OEPjTpz5kzrZxbpn4KCkBDCMWTZUPY/Edz8A0JyvtCEk/AAIIf7s7DzKXi/48ePf/nllwcPHmy98uLFi19//fW+fftaIvPUqVNXrlzRLJeXlxcWFhJCDh8+TAjJycnJysrSvpcpSVJmZuaXX3559uxZzRqWZZ944gnN4wMPHz5869atwsLCjRs3FhcX66CFOqKsINQY50KPJXEXq6VRqfxPd9A1BAAZtM7CrkjBNWvWzJgx48SJE88///zzzz+vWblhw4awsLCDBw+uWLFizJgxmofdr1q1avv27ZoNioqKXnrpJUJIdnY2IWTbtm1btmy5deuWNkcUBGHq1KkrV648duxYRETEf//7X0JIY2Oj5pFEhJCFCxfOnTt31apVubm5Q4YM+f7773XV2E5S6LOLHCzJzmjukx/FUan834azs/oo8QsBAOjcqwXtu6dVQHc6NZsXRRLvTv91un0nqd4OYa3N2v7Rzz//vHz58jNnzjg7O1dVVXl5eT3//PO+vr7JyckpKSnh4eGiKIaHh69fv37hwoVt7mHJkiXLly9fu3at9k9WT01NLS4uPnz4MMMwjz322PTp06dPn37PNt7e3u+99x4hRK1Wb9iw4d1339W2tV1JoUFICKGELPBjwnrRmTlCxlVpXdgD/5cCANBS326Ub0+c3aiXLBlKWMnegvaxbV93kH3w5kVFRX5+fpqH/3Tr1i00NFSTT5IkhYeHE0IYhpk4cWIHHmxUWVm5bt06QgjLsosXL279o7y8vOjoaM3lwNGjR1dXV5eUlLi4uLTeJiIiQrPQp0+fQ4cOtffoXUS5Qajh350eSuSWHBYGbeM3jmVDHTFkGQA67vn+7Ti99GOlNDNHzIxlHSxoQgY/2ZMZ66ybj6CqqqrWLymlmvGfoijev7JlgRDy8OGdLW9s+W9rd+7csbOza73m/kGn5ubmmgVNKmvRFH3AKUGi4sj7oezfhjOTM/llRwVMNAQAPWh9XfDhcyo65uTJk2VlZYSQ6urq/Pz8kJAQHx8flmXz8/MJIZIk7dmzJzg4mBDi7Oz8008/ad61f/9+zYK5ubmZmVltbe09u7Wzs1u6dOnSpUv/8Ic/3H/Q7OxsTdYePHjQ2traw8NDV83pUkrvEbZ4zJMJdqCzcoWDN/ivIlgXNbqGANBV7h8do8nChAz+nWGsTvqFvXv3njRp0pgxY/bu3ZuUlBQYGEgIWbNmzZQpU2bMmHHy5Mnm5uYFCxYQQmbNmqWZ6ldZWXnz5k3N2xmGiY+Pj4uL8/T0XLFihZ+fn5bHjYmJCQwM3Lhx49/+9jeO43Q1naNLIQh/5W5N98Vybx4Xgrbzn47m4nojCwFA9x40RlS3Wejh4fHFF19kZGTExMS0XJl7+umnw8LCCgoKYmNjIyIizMzMCCEBAQFHjx799ttv+/Tp079//wsXLmg23rp165kzZ0pLS52cnLQ8aFJSUlxc3LFjxzIyMgICAgghFhYWubm5FhYWhJAPPvjAy8tLs6UmpDvZRl1BEP6GZqJhjBvzZI6Q4U7fHca2+RxgAICOefhMCd1moYuLy5w5c+5Z2adPnz59+tyz0svLqyWievTooVlgGMbf39/f379dBw0MDNT0PjUopS2BFxIS0rLezc3Nzc2tXXvuOrhG2IYRjvRoEldWT0am8ueqcM0QAHSmqol8Pe5h8wUdVWRnNNfJ54r36tVryJAhndpF+/n5+bWkqXFBj7Bt3czJpnHsV+fFsFT+T4PZ5AB8YwAAHRiuxdD0XioS7dqR7uClS5dUKhUhJCIiouV0qN68/PLLHXvjsWPHNOdO5YLP94eZ3ZfZH8/9+5w4LVuokPNOeAAAj6ZJQaNjaWl5/2QMfUIQPkJ/O1owiXNRkyHb+e/KcJoUAMDUIAgfzZIl74ey74cyU7L4ZUcFAWkIAGBCEITamuTBHH/MLL9Mikrnr9UiDAEATASCsB2cVGTPRG6SBxOSwqeW4BlOAACmAEHYPpSQ5ABmRxT3SoE4O1eoM4J7JgAAwMNg+kRHDOtJjyZxzx8UQnbwm8axA+1xDxoAJTp06NDq1avlrkIRTpw40XU7RxB2kK0Z2TiW/eq8OC4dEw0BlGjRokV37txpuTmn3jQ3N3McJ+98A/1zdnYePHhwF+0cQdgps/sywx3pzH3CtzekT8NZezmnhAKAXpmZmb3zzjv6P25dXZ2lpaXmsX+gE/hVdla/bjQ/kfPtRoZu5/ffwGhSAAAjgyDUAQuWvBPCfhLOPpkjLDksdPImgQAAoE8IQp2JcqVFk7mTt6XRafylanQNAQCMA4JQlxxVZNcE7nFvZngKv+kCOoYAAEYAQahjmomG6THcG0Xi7FyhFhMNAQAMG4KwSwQ70KNJHCEkeAd/vBynSQEADBeCsKvYmJGvItgVQUzMHv79H0SEIQCAYUIQdq1pXkxBIvf1RXFyplCOJxoCABgeBGGX87Sh++O5IT3IwG+a915FzxAAwLAgCPWBY8iyoeyGCG7eASE5HxMNAQAMCIJQf8a70GNJ3MVqKSyVv3AHXUMAAIOAINSrnpZkZzT3pA8TupPf+BM6hgAA8kMQ6ptmomFOHPfOCXF2rlDTLHdBAADKhiCUh393WjiJ625BgnfwR2/hNCkAgGwQhLJRceT9UHZlMDNxL7/sqICZhgAAskAQymyKF3N4Epd9XZqwhy+tk7saAADlQRDKz92a5sRyYb1o0I7m9CvoGAIA6BWC0CBoJhr+dyz33EEhOV9ownhSAAB9QRAakAhnenIKd6OejNzJn69C1xAAQB8QhIbFzpx8PY5d6M+Ep/Hrf0THEACgyyEIDdHsvsy+OO7D0+L0bKGySe5qAABMGoLQQA2wo4cmcc5qMmQ7n1eG06QAAF0FQWi4LFnyfii7ZgTzWBYmGgIAdBUEoaGb7MEcnszlXJeidvMZ18T7n3fPi2QDblsKANBRCEIj0NuK5sRxie7MU7nC1GzhSKtbsgkSmf2tcL1WxuoAAIwbgtA4MJQkBzA7ozhelCbu4b8rkwghgkSeyhUC7eniQPw7AgB0ED5AjclwR3pyillYLxq5m99cwiAFAQA6D5+hRsbWjOyI4v4UyDxbYF7HS0hBAIBOwseo8REkcqqSzPQQ0q9Kiw8LcpcDAGDc9BSEDQ0Nhw8fLigo0M/hTFjLdcF1w5tzY9kPTokzczCzAgCg4zj9HGb27Nl1dXUlJSUnT57UzxFNUuvRMdXVZGQv5lAiDUvlEzLIjijWDN17AID209Nn5+bNm9etW6efY5mwXSXisJ6/GR0zyJ4emcSdqZAm7uHvNMtYGgCAsUInwpgkejAvB9z7T9bPjp6fwfl2o6NS+Wu1OEkKANA+Oj41OmfOnMrKytZrPv30UwcHB90eBe7BUrIujH3/B3F0mpAew/brRuWuCADAaLQvCCsrK9Vqtbm5eeuVly5damho8PPzo5S+9dZbgvCbcYzdu3fXQZmgheQAprsFGZPGfxPJhfVCFgIAaEXbIHz88cezs7Nv3bq1efPmadOmaVbyPP/4448fOnTI1tbWwsIiIyPDxcWlzbfX1tZWVVUJglBRUWFlZXVPlIKuzO7LuKhpUia/Loyd6oXz3gAAj6btZ+WTTz6Zl5fn7+/feuWWLVvOnDlz9uzZU6dO+fj4rFq16kFv/8c//vHqq6+6uLhMnz49Nze3MxXDw0W60vQYLjlf/OgM7sQNAPBo2vYIExISCCGU/uaE29dffz1r1iy1Wk0ImT9//jPPPLN69eo23/7GG29oc5SKioq+ffu2XtO7d++UlJQ2N66pqdFmn6bqIc3vZ0H2jqNJ35qfuy0sG8Sb5ElShf/r19bWiqJ4z9+jctTV1fE8zzAKPedRX1/f3Nys5Oabm5uzLKvl9mq1+pEbd2qwzOXLl2fOnKlZ9vb2vnbtGs/zHNfxfdra2qalpbVeY25ubmNj86DtH/IjJXhI8wNsSMFkkpjBLzxq/km4aU4xVPK/PqXUyspKsUHIMIxKpVJsErAsa2lpqdjmcxzXriDUap+deXNtba2lpaVmWaVSiaJYV1dna2vb4R2yLNuvX7/OlAQteliQjIncjGw+bi//TSRnYyZ3QQAABqlT3yl69ep1+/ZtzfLt27fVanVnUhB0zoojKdGclw0dn87frJe7GgAAg9SpIAwKCsrLy9Ms5+XlBQUF6aIk0CWWko9HsVM8mdCd/LkqTLcHALiXtqdGd+3ade3atdu3b2dlZVVUVCQkJDg7Oz/33HPDhg0LDQ3t1avX8uXL165d26W1QoctDmQcVWR0Gr89igt1VOiFJQCANmkbhOfPnz9z5kx8fDwhpKioaOzYsYQQPz+/Xbt2rV27tqGh4e9///uUKVO6sFLonN/5Mg6WZFIG/+8xXFxvZCEAwF1UkgzldFlxcXFERERxcbGW21dXVyt53GDHml/4szQ5k18+lJ3vZ9xDzhT+r19TU6PkUaO1tbVKHjVaV1en5FGj7Z0+oQ2F/ioVa1hPeiCeW/29uOwonugLAEAIglCBfGzp/ngutUSat1/gcfMZAFA8BKESOanIt3Fcab00NVuo4+WuBgBAVghChbI2IzujuJ6WZFw6/3OD3NUAAMgHQahcHEPWh7MxbnTkTv6nO4YyZgoAQM8QhIpGCVk2lE0OYCJ2CUdvIQsBQIkQhEBeGsD8cyQTs4ffcxVZCACKgyAEQgiZ7MHsjObm7Rc2/oSBpACgLAhCuGuEI82KZf/vCKYYAoCyIAjhV/3taH4il3JZWpgviDhLCgDKgCCE33BWkwPx3NlKaWq2UI8phgCgAAhCuJe1GUmdwFmyZHw6X94odzUAAF0MQQhtMGfIxrHsaGc6OpUvqcFJUgAwZQhCaBsl5J0QdoEfMzJVOF6OLAQAk4UghIdJDmA+CGWid/MZ15CFAGCaEITwCI95MpvHc0/l8psuYIohAJggBCE8WoQz3RfLLT4svnsSWQgApgZBCFrx707zE9mNP4nJmGIIAKYFQQjaclHT3DjueLk0fZ/QgJvPAICpQBBCO3S3IBkTOY6S2D18VZPc1QAA6AKCENrHgiX/G8cGOdBRqfzVWpwkBQCjhyCEdqOErB7OPtOPCd0pnLyNLAQA44YghA5KDmDeDmGid/MHbiALAcCIIQih42b1YTaO5aZm85svYloFABgrBCF0yngXmjWRW1QofngaWQgARglBCJ010J5+G8f+87SYnC/gJCkAGB0EIeiApw39LoE7Wi49nSs0o2cIAEYFQQi6YW9BMidy9QKZuIe/0yx3NQAAWkMQgs5YsmTTOLZvNzoqlb+GKYYAYCQQhKBLLCX/CmOf9GHC04SzVchCADACCELQvcWBzF+GMmPS+O/KkIUAYOgQhNAlnu7LbIjgpmbxaSXIQgAwaAhC6CpRrnRnNLfgIP/RGQwkBQDDhSCELhTSkx5M4Nb8IC45jCmGAGCgEITQtbxtaF4id+CGNHc/phgCgCFCEEKX62FBMiZyP9dLcXv5akwxBAADgyAEfbDiSEo052VDx6fzN+vlrgYAoBUEIegJS8lHo9jY3nRkKn8eUwwBwGAgCEF/KCHLhrJ/HMyMTReKbiELAcAgIAhB3+b6MutGMhP38OlXkIUAID8EIcgg0YNJm8DN289/8iMGkgKAzBCEII9hPenBBG719+Kyo4LctQCAoiEIQTY+tnR/PJdaIs3bL/DoGQKATBCEICcnFfk2jrteJ03NFup4uasBAEVCEILMrM1IajTnYEnGpfM/N8hdDQAoD4IQ5Mcx5JNwNsaNjtzJ/3QHQ0kBQK8QhGAQNFMMkwOYiF3CsXJkIQDoD4IQDMhLA5h/jmQm7Ob3XEUWAoCeIAjBsEz2YHZGc3P38xt/wkBSANAHBCEYnBGONCuW+78jmGIIAPqAIARDNMCO5idyKZelhfmCiLOkANCVEIRgoJzVZH88d7ZSmpot1GOKIQB0GQQhGC4bM5I6gbNkyfh0vrxR7moAwEQhCMGgmTNk41g23ImOTuVLanCSFAB0D0EIho4SsmoYu8CPCdohrDx+7/CZyiay6BDG1ABAx3FyFwCgleQAxtGS/G6/UNtM3gphNSurmkhCBr9oIL7PAUDHIQjBaMzsw3S3pJMy+MomsmoQqWoi8Rn8ooFMogeCEAA6DkEIxiTGjeYlcqNS+WvVZrd5pCAA6AA+RMDIBDnQ/fHc7uvsmQop/6ZUcFPCREMA6Ax9BCHP80888URgYOCQIUNWr16thyOCCatqIq8eEj4b3tSvG/2xUnopT3Da2Dw7V9hySaxplrs4ADBCeuoRPvvssydOnPjuu+82bNhw+PBh/RwUTE/LdcHJ7mJWLFfHk4X+zOHJ3Cgn+tV50e1/zVG7+fd/EEvr5C4UAIyHPq4Rchw3ZswYQoharXZ2dq6vr9fDQcH0tB4dU11NVBzZEcVNzuQJYRb4MQv8mIpGknVdTL0sLT/W7G1D491pgjsT5EDlLhwADJqOg/D777+/J+cCAgLUarVmeceOHfX19aNGjdLtQUEhTldKSwLZuN6/BpsmC987dfc5Fd0tyDQvZpoXESQ2v0zacklMyhTMGKJJxDHO1AzXxAHgPjoOwh07dpSWlrZe88c//lEThDk5OW+++ebu3bsZBp9G0BGhjm307VQcWRp47/9RLCWjnOgoJ/b9UHKqQkorkZYdFc5WSWOdmXh3muTJ2JjppWIAMAbaBuHhw4fz8/N//PHHxx9/fPTo0S3rjx49umbNmsrKyqSkpLlz577xxhttvj0vL+/1119PTU3t2bOnDqoG0Jp/d+rfnS4OZEpqpD1XpS2XxIX5QkhPGt+bmepFXa1w4hRA6bQNwlWrVqnV6v379w8aNKglCEtLS8ePH//nP//Zz8/vxRdfFEXxmWeeuf+91dXVMTExo0ePXr58OSFkzpw5oaGhumoAgJbcrekCP7rAj6njSfZ1cctFadlRwUVNp3njUiKAomkbhFu3biWEjB07tvXKzz77bMyYMa+88goh5J133lm+fHmbQahSqXJyclpeenp6PugojY2Nn3/+ees1tra2SUlJbW4sCIIgKPcmk2h+h5tvQUmsK4l1JYLEFNwku65KT+YIjSKJciFxvckEVyO4lKhpPqUKDW9N8yVJoRNI0fx2/e0zDPPIv5ROXSMsKioKCwvTLIeFhZ0+fbqurq5laMyvx+C4oKAgbXbY1NS0b9++1mt69uwZGxv7oI0bG5X7bB40XyfND7YjwXbkLwHkTBXdfY1ZdYLO3U8inMhEVzHeTbI1M9DPmsbGRo7jFBuEjY2NDMModrRBY2MjpVTJzZckiWVZLbe3tLTs2iAsKyuzt7fXLPfo0YMQcuPGDW9v7w7v0MbG5j//+Y+WGwuCcH/oKgear9vmB6lJkDP5UzD5uYHsviJuuSS+fFga3INO82KmeFE3A7uUKIqiWq1WbBBKkqRSqRSbBIQQS0tLxTafUmpubq59EGqjU0GoVqtbJktoFmxsbHRQFIB8elqS2X2Z2X3vXkpMK5He2SHYW9AEDxrfmwlzUmr4AJiuTgWhu7t7cXGxZvnSpUsqlUrTLwQwAWqOJLgzCe5kXRibXyalXREXHBTqeDLBjca70wlujLlCv5EDmJpO/SlPnz5969atVVVVhJDPPvts2rRpiu2tgwnTzEp8J4Q9PZXbF8cOsKOrTogOG5oTMvivzouVTXLXBwCdo22PcN68edu3b6+uri4sLPzjH/+4YcOGuLi46OjoiIiIAQMGODk5VVdXZ2ZmdmmtALLztqHJATQ5gLnVQNKviGkl0sL85oDuNMGdSfKkvt1w3hTA+FAtx+DW1tY2Nf361dfa2trM7O7NOS5fvlxVVeXv79/Jq5fFxcUREREt51ofqbq6WsmXJNF8A2l+PU++K5NSS8Qtl8Tu5nq6lFhTU2NlZaXYwTK1tbVKHixTV1en5MEy9fX1sg2WsbKysrKyavNHHh4euqsHwMioOBLpSiNd2TUj2GPlUmqJ+HKBcKVWmuDKTPOm0a6MhS7/YAFA9/CEegDdYCgJcqBBDuyyoeRitZR6WfrglPhUrhDuRBPcmUkeTC+V3CUCQFsQhAC613IpsbyR7Lsupl6WFhc2+3fXJCL1s1PoKU0Aw4QgBOhCPX55MlSDwB68IaWWiJG7RQuGxLvTaV7MyF6UQSYCyA1BCKAPluzdS4maJ0NtuSS+XCCU1EgxbkyCB53oxljjyVAAMkEQAuibf3fq351dNpQUV0sZ16SvzovzD9y9lJjowTjhUiKAfil0AC6AIfC0oQv8mNRo7uIMs9l9mYM3pAFbm4N38MuOCqcr753XdOSWtPWSeM/KBoGsOCYa6K3BAYwEeoQA8rP/5VIiL7IFN6Utl8QJuwVzhsS70wR3JsKZcgwZYEeXHhaaRPKEz93vr40CmZ4txLkrdTohgI4gCAEMCMeQUU50lNPdS4lpJdKyo8LZKmmiG5PgQTdGcE/m8oSQxF6kUSDT9wlx7vRZP5zXAegUBCGAgfLvTv2708WBzOUaae/Vu5cShzrQ5UeFq97MwXIRKQigEwhCAEPnYU0X+NEFfkxlE9l9Rdx8UVp6wszBXAxyYAt/loIdMAcDoFPwdRLAaNiZk8c8GUEibw/mPW3IyQrxhe+EXhubp2cL638Ub9TLXR+AcUIQAhiNRoFMyxbi3OkLvsK+WO5OE3l1IPPDFLN4d5p17e6I0+JrhGUAABioSURBVCWHhYM3JAwkBdAeTo0CGIeWFHzWj6mpIWqOpERxkzJ5QsjsvszsvndHnKZduTtVP8KZ0Qw67W4hd+kAhg1BCGAcjpVLj3nSOb6/nsVRc2R7JLfimDDTh9BWI05JyN2p+mkl0ssFzd42NNJVHw+HAjBSCEIA4zDCkY5wvDfIrM3Iu8PaeM6Tp41mfA26iQCPhiAEMGWtu4mXqqXMa1JaiZSc3+xjSzWJONQB3URQOgQhgFJ4/dJN1DwKI+u6ODtXuNUoTXBlEjxolCtjZy53iQByQBACKE7LozDeCSEXq6Wsa9KWi9L8A8190E0ERUIQAiia9y/dxHqe/a5MyrouPpUr3G6Uol2ZBA8a7cp0QzcRTB2CEAAIIUTF3dtN/Oq8+LtvhcE9aII7E+lK0U0EU4UgBIB7/dJNZOp58l2ZlFoiTskSKSXRrjTSlU5wY2zxGGEwIQhCAHiglm7i+6FtdxODHNBLBKOHIAQArbR0E+t4klcmpZaIj2WJLCVRrjTSlca4MTboJoJxQhACQPuof9tNTL0srf9RnH9ACOlJI12YBA86wA7dRDAmCEIA6DhvG5ocQJMDmFqe7LsuppVIMbtFM0aTlHSiG2ONbiIYPAQhAOiAFUcS3JkEd0JIG93ERA/aH91EMFQIQgDQsfu7idG7RXOGRLrSeHca5cpYtnF7VADZIAgBoKu07iaeqpDSSqQPTolP5QohPWl8b2ayJ/WwRjcR5IcgBAB98O9O/bvTxYFMTTPJKRXTSqTVqaLFL93EaFfGAt1EkAmCEAD0ytrsYd3EJE/q3lY3sVEgVU3EUXXv+iu1Um8rdCuhU5hHbwIA0DU0fcTMidyFGWYL/JjTldKInbzP13xyvpB1TWoSf93yTKUUn8HfqP/N27cVi0/nCoKk56rB1CAIAUB+PSzINC/m41Hs1Zlmm8ezLmq66qTg+J/mhAx+/Y/ilVppcA+6ZgSbsJe/Xnc397YVix+cEndEcyw6hNA5ODUKAAaEoSTIgQY50MWBzK0GklMqZl2T/lIkqDka705n9WESM4RNo8iJm9La09LOaA53PYXOQxACgIFysCTTvJhpXmRdGFv4s5R+Rdzwk3iuShq6y8zbVjwQb4abuoFO4NQoABg6lpJQR/rXIPbIZO69UNaWI6cryLHbuDYIuoEgBACjsa1Y/Oq8eCSe/9NgJnIXv69UfPR7AB4Fp0YBwDhoRsfsjObYpsY/DaaUMjG7hb0TyVhnfKGHTkEQAoAROPyz9MEpMTWaszEjtU2EEPLGEFaSSMwe4cgkOtAeI0eh4/BNCgCMQJAD3TWBu2d0zJ+HsutHsTF7hFMVuF4IHYceIQAYAYYSq7Y+rp7uy1iyJHq3kBnL4jmI0DEIQgAwbjO8GUJIVDqyEDoIQQgARm+GNyNJyELoIAQhAJiCx33QL4QOQhACgIl43IdpEEhkOp8dy/VHFoLWEIQAYDrm+DKEkMh0ISuWRRaClhCEAGBS5vgyErIQ2gNBCACm5ne+DCEkZo+QHcv2sUUWwiMgCAHABP3Ol5EIGZ+OLIRHQxACgGma68sQZCFoAUEIACarJQv3xbI+yEJ4AAQhAJiyub6MJJFxyEJ4MAQhAJi4ef1+PUeKLIT7IQgBwPQhC+EhEIQAoAjz+t0dR7ovjvW2QRbCrxCEAKAUz/RjCCHjdiEL4TcQhACgIM/0YyQJWQi/gSAEAGWZ73e3X5gTx3ohCwFBCAAKdDcL04V9schC0EsQiqK4d+/ec+fO2djYTJ482d7eXg8HBQB4iPl+jIQsBEIIIYwejsHzfG5urr29fWlp6fDhw+vq6vRwUACAh1vgxywNZMalC5eqJblrATnpo0dobm6+atUqzfLOnTsvX77cv39/PRwXAODhFqBfCHq7Rsjz/N///vfi4uIhQ4b4+fnp56AAAI/0rB/TwJOxu4TcONYTWahIugzC5ubmmJiYe1ZmZ2cTQiil3t7elNKtW7dWVVXZ2dnp8LgAAJ2RHMBIhIxNF3JikYVK1I4gbGpqunr1qqurq4WFRctKQRCOHTvGMMzgwYPNzMy++eabNt/Lsuy0adMIIadPn87Kypo6dWon6wYA0KGXAxiCLFQqrYKQ5/nQ0NCTJ082Nzfn5+cPHz5cs76ysnL8+PGSJAmCoFarMzMz2+zqlZaW3r5928/P7/z58/n5+YsXL9ZlCwAAdAFZqFhajRplGOatt94qLS21sbFpvX7t2rWOjo5FRUXHjh1TqVQff/xxm29vamr661//Ghoaunjx4jVr1mCkDAAYppcDmGR/JnqPcK0W40gVRKseIcMwUVFR96/funXr0qVLKaWU0qeeemr9+vWvvfba/Zt5eHhs2rRJmwPdvn3b2dm59Rp3d/esrKw2N66pqdFmn6YKzZe7BDnV1taKokipQnstdXV1PM8zTJfM/prnQRoauTFp4q6IRhd1Vxyhs+rr65ubm7uo+Yavvr7e3NycZVktt1er1Y/cuFODZa5cueLp6alZ9vLyunLlSmf2Rgixs7MrKChovYbjuHu6oa095EdKgObLXYJsKKVWVlaKDUKGYVQqVdclwZIgYmEhJuxncmJZVyuD+yWzLGtpaanYIOQ4rl1BqNU+O/NmTTJrli0sLGpraztZDcMwLi4undwJAEAnvRLw63PtDTALQbc6FYROTk63b9/WLN+6deues5oAAMbr1YF370eaE8e6qJGFpqxTnethw4YdOHBAs3zgwIFhw4bpoiQAAIPw6kBmgR8zdpdwvQ5jZ0yZtj3C9evXV1RUNDY2btiwITc397nnnuvWrVtycnJsbKynpyfP8x9//HFOTk6X1goAoGevDWQIIWN3oV9oyrQNwjt37lRUVLz88suEkIqKClEUCSEjR47ctm3bl19+yTDMrl27hg4d2oWVAgDI4bWBjIQsNGnaBuHrr7/e5vpx48aNGzdOd/UAABic1wfefZbvPmShKVLoAFwAgHZ5fSAzrx8zDtcLTRGCEABAK4sG3c3CUjxT1bQgCAEAtLVoEDO3HzN2F48sNCUIQgCAdvgDstDkIAgBANpHk4Xj0pGFJgJBCADQbn8YxMzpy4xL52/Uy10KdBqCEACgIxYHMjO86dhdyEKjhyAEAOigZUNZZKEJQBACAHTcsqHsdC86DllozBCEAACdsjyIneZFJ+zmbzXIXQp0CIIQAKCzlgexkzzo+HRkoVFCEAIA6MCKIHaSB41EFhohBCEAgG6sCGITkIVGCEEIAKAzf0UWGiEEIQCALv01iI13p5HpfHmj3KWAdhCEAAA69mYwstCYIAgBAHTvzWA2rjey0DggCAEAusSbwWwsstAYIAgBALrKSmShMUAQAgB0oZXB7EQ3GpXO30YWGioEIQBA13orhI1xo5HIQkOFIAQA6HJvhbATkIWGCkEIAKAPb4ewE9xo1G5kocFBEAIA6MnbIWyUK7LQ4CAIAQD0551fsrACWWgwEIQAAHr1dgg7qhey0IAgCAEA9IoS8l4oG4YsNBgIQgAAfdNk4UhkoWFAEAIAyIAS8n4oO7IXjdvL32mWuxplQxACAMhDk4VBDjRmN7JQTghCAADZUEI+GIkslBmCEABATi1ZOHEPX40slAOCEABAZposHNKDxiAL5YAgBACQHyVkLbJQJghCAACDoMnCwT1wjlTfEIQAAIaCEvLPkWwgslC/EIQAAAYEWah/CEIAAMOiycJB9jR2L1+DLOx6CEIAAINDCfkwjB3YnU5EFnY9BCEAgCHSZGEAsrDrIQgBAAwUJWQdsrDrIQgBAAyXJgudLGnfzc13mn7zo7J6kpgh8KJMlZkQBCEAgEGjhGyOZD2sqe+WX7PwZj1JzOBfGchw+BTvNPwKAQAMHSUkL5FzU9/Nwpv1JCGDf2cYO9aZyl2aKeDkLgAAAB6NoaRwMheygx+cZuZmLa4ejhTUGfQIAQCMA0NJ6gSuXiQ/3ZHcreSuxoQgCAEAjMPNepKUyX89WhjnTAO+4bddwjgZ3UAQAgAYgZbrgqMdxa3jmWf6MU/vFxYXCoIkd2XGD9cIAQAMXevRMXV1hBCydiTb3YJ+dlYo/Fn6ehznqJK7RGOGHiEAgKETCXkv9N7RMSuCmH+OZEc70ZAUvuAmOoYdhx4hAIChc1IRJ1UbY0STPJkkTxLcU5ycyS8NZJMD0LfpCPzWAACMW4I7cyCe++yc+FSuUMfLXY0RQhACABi9vt3ooUTOjCFhqfyFOzhN2j4IQgAAU6DiyOej2ef7MyN28tuLMbOiHXCNEADAdCzwY4Y60GnZwqGfpZXBLIubz2gBPUIAAJMS7EAPT+KO3pKi0vmb9XJXYwwQhAAApsbBkuyO4UY705AU/hBmVjwKghAAwASxlCwbyv5zJDMpk3//B1wyfBgEIQCAyWqZWTEbMyseTH9BWFJS0r9//23btuntiAAAoJlZwWFmxYPpLwhfeeUVHx+f6upqvR0RAAAIZlY8ip6mT/z73/8ODw+/ePGifg4HAAD3wMyKB9FlEFZUVBQVFbVeY2dnFxwcXFpaumnTpvT09FdeeUWHhwMAgHbRzKyYmcNHpfOb8MyKX+gyCG/fvp2Zmdl6jaenZ3Bw8JIlS6Kjo3Nycq5cuaJSqa5fv+7i4qLD4wIAgJYcLMmeGO6vx4SQFH7zOHa4IzqG2gVhZWXlzp07jx8/Xl9f/69//atlPc/z77777p49e5ycnJYuXTpkyJBVq1bd//ZRo0ZduHAhKyvrypUrzc3NpaWlCEIAALloZlYEOYiT8MwKQoiWQXj+/PmtW7c6Ozt//vnnrYNw5cqVqampa9euLSwsjIyMvHDhgp2d3f1vnz9/vmZh4cKFQUFBQUFBOikdAAA6LMGdORBPp2QLRbekj0axagXfcJNKkrajab///vuhQ4c2NzdrXvI87+bm9vXXX48ZM4YQMn78+KSkpJdeeukhe7h165aFhYWNjU2bPy0uLh4xYsSSJUtar+zWrdusWbPa3L66uvpBu1ICNF/Jza+pqbGysqJUoSe1amtrVSoVwyi0H1NXV2dpaamr5tfzZOEhcvw22TyWeFnrZJddq76+3tzcnGVZLbfnOO6Rfykd/w5QWlpaVlY2YsQIzcsRI0YcO3bs4W9xcHB4+AY8z587d+6etwiC0ObGgiA86EdKgOaj+YoNQk3ztf8Sb2J023xzSj4aQT47T0em0XUjpMTehv5bbe/fPsuyXRiEZWVlarXawsJC89Le3v748eMd3puGtbX1unXrtNy4ubnZ0tKyk0c0Xmi+kpvP87ylpaWSg1CHXSKjI4qizpv/4kAy3Fmali0craQGPrNCkqR29Qi10fFfpa2tbUNDgyjenZtZW1vbrVs3HVUFAAB6pZlZUaTIZ1Z0PAhdXFwYhikuLta8vHDhgru7u26KAgAAvdPMrFDgMys6HoTW1tYJCQmaQaRXr15NS0ubOXOm7goDAAB9U+YzK7QKwps3b1JKBw0axPM8pdTNzU2z/t133925c6e/v39gYOBLL70UGBjYlaUCAIA+KO2ZFVoNlnF0dGxzhFKfPn3OnDlz6dKl7t2729vb67o2AACQh+aZFS/mCWGp/DeRrLeNAY+f6bTOjjtiGMbHxwcpCABgYlqeWTE8xcSfWaHgewkAAMCjKOGZFQqdiAMAAFoy+ZkVCEIAAHgE055ZgSAEAIBHM+GZFQhCAADQlknOrEAQAgBAO2hmVrAMCUvlL1abwmlSBCEAALSPiiP/NqGZFZg+AQAAHWEyMyvQIwQAgA4yjZkVCEIAAOg4E5hZgSAEAIBOMfaZFQhCAADQAeOdWYEgBAAA3TDSmRUIQgAA0BljnFmB6RMAAKBjxjWzAj1CAADQvWAHWmgkMysQhAAA0CV6GsnMCgQhAAB0FaOYWYEgBACArmXgMysQhAAA0OUMeWYFghAAAPTBYGdWYPoEAADojwHOrECPEAAA9MrQZlYgCAEAQN8MamYFghAAAGSgmVmxNlT+mRUIQgAAkE2ih/wzKxCEAAAgJ9lnViAIAQBAZvLOrMD0CQAAMAhyzaxAjxAAAAyFLDMrEIQAAGBA9D+zAkEIAACGRc8zKxCEAABgiFpmVkzO5P9cdG8cNgjkxTyhUdDBgRCEAABgoDQzK7pb0H+dERYc/DX0mkQyPVsY2J1asDo4CoIQAAAMl2ZmxZvB7IafxKQsnhDSJJKpWUJsb/pcf91EGKZPAACAoXvWjxnqQMen86P3sr1UUrw7o6sUJOgRAgCAUQhxoN8/ZnahmhIi6TAFCYIQAACMQpNIfp8nLAsU3a3pHwp1MUjmFwhCAAAwdC3XBef3FT8YQWt5osMsRBACAIBBu2d0DCXknyNZHWYhghAAAAzatVpphvdvxohqsrCXiupkHiFGjQIAgEHzsqFeNvfegZsS8tpA3fTljLhH+PrrrxcXF8tdhWyWLl167tw5uauQzV/+8pcffvhB7ipk8+abbx49elTuKmTz7rvv5ufny12FbNasWbN//365q5DNhx9+mJWVpdt9GnEQfvfdd3fu3JG7CtkUFBRUVFTIXYVsCgsLy8vL5a5CNkVFRTdv3pS7CtkcO3asrKxM7ipkc+LEievXr8tdhWy+//77K1eu6HafRhyEAAAAnYcgBAAARaOSpI+nPWmjvLx84sSJpaWlWm5fVVVlbW3Nsrq45aoRunPnjlqt5jiFDneqrq5WqVRKbr6lpaWZmZnchcijpqbGwsJCyc03Nzc3NzeXuxB51NbWchxnYWGh5fZ79+4dMGDAw7cxoCAkhNy+fbu2tlbuKgAAwEQ4OTk98juTYQUhAACAnuEaIQAAKBqCEAAAFA1BCAAAioYgBAAARTPKIHz77bcnTJjQp0+ftLQ0uWuRwYsvvujp6WlmZubr67tx40a5y9G3efPmOTo6sizr5eW1bt06ucuRx82bNwcOHLhgwQK5C9G3JUuW+Pyif//+cpcjg08//dTb25thGB8fn+PHj8tdjl6NHDnSp5VFixbpas/GOg1r7ty5b7zxRk1NjdyFyKBHjx6ZmZne3t579+6dMmXKwIEDBw0aJHdR+jNnzpzVq1fb29vv379/4sSJwcHBw4YNk7sofXvppZdsbW0VeJuxn3/+efr06fPnzyeEUHrvXZhN3qZNm1asWLF58+aQkJBLly5ZW1vLXZFebdq0ied5QogoiiNHjhw+fLiu9myUQbh06VJCyNtvvy13IfJYsWKFZiE2NrZv374nT55UVBCGh4drFkaPHu3l5XX58mWlBeHOnTtFUUxMTMzLy5O7FhnY29t7e3vLXYU8Vq5cuXLlyhEjRhBC+vTpI3c5+ubu7q5Z2Lt3ryiKCQkJutqzUZ4aBY1Lly5duHBBaTFACCkqKtq0aVNycrKNjU1sbKzc5ehVeXn5okWL3n//fbkLkc0//vEPZ2fn8PDw9PR0uWvRq+bm5jNnzpSXlw8dOtTX1/eNN95obm6Wuyh5fP75508//bT2N5d5JKPsEQIhpKamZsaMGa+++qqvr6/ctejbkSNHMjMzT5w4MWHCBKXdZe3ll19+/fXXXV1d5S5EHnPnzl28eLGtre2uXbumTJly4MCB4OBguYvSk5s3bwqCsG3btvT09Pr6+vj4eDs7u9dee03uuvStvLw8JSWlsLBQh/tEj9Ao1dfXJyYmDho0qOU0qaI8++yzW7duPXXqVF5e3qeffip3Ofqzf//+wsLCIUOGFBUVXbt2rbKy8sSJE3IXpVdhYWG+vr5OTk7z5s1LTEzcvn273BXpT48ePSilr776qpOTk5eX1wsvvKDM0YL/+c9/Bg8erNvrQcr6Nm0ampqapk2b1rNnz48//liB4wVamJubDxw4sKSkRO5C9KexsdHd3V1zjby4uLiqqmr58uXbtm2Tuy55KO1/fktLSy8vL7mrkN8XX3zx/PPP63afRhmE58+fv3PnTl1d3aVLl4qKinx9fW1sbOQuSk8kSUpKSiorK1u3bp1m8LSbm1uvXr3krktPqqqqUlJSxo4da2FhcfDgwR07dqSkpMhdlP5ERUVFRUVplletWpWXl6e0FFy/fn10dLRard69e3dKSsq3334rd0V69dxzz61ZsyYsLKyhoeGjjz6aN2+e3BXp25EjR86ePTt9+nTd7tYog/CDDz7Iz8+3tbX95ptvvvnmm48++kg51wmam5s1g+ZfeOEFzZrk5OSnnnpK1qL0h2GYlJSUP//5z42NjT4+Pp9//nlERITcRcnDycnJx8dH7ir0bc+ePW+99VZDQ4Ovr+/WrVuVNlLslVdeKSsrGzp0qEqlmjVr1u9//3u5K9K3gwcP/v73v7ezs9PtbvH0CQAAUDQMlgEAAEVDEAIAgKIhCAEAQNEQhAAAoGgIQgAAUDQEIQAAKBqCEAAAFA1BCAAAioYgBAAARUMQAgCAoiEIAQBA0RCEAACgaP8PPcXDCU9tD18AAAAASUVORK5CYII=", "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.6.0" }, "kernelspec": { "name": "julia-1.6", "display_name": "Julia 1.6.0", "language": "julia" } }, "nbformat": 4 }