{ "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", "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 Eₙ-Eₙ₋₁ ρout-ρin α Diag\n", "--- --------------- --------- -------- ---- ----\n", " 1 -7.844807293676 NaN 1.97e-01 0.80 4.0\n", " 2 -7.850320248871 -5.51e-03 2.94e-02 0.80 1.0\n", " 3 -7.850646832409 -3.27e-04 2.96e-03 0.80 3.0\n", " 4 -7.850647500193 -6.68e-07 4.54e-04 0.80 2.0\n", " 5 -7.850647510503 -1.03e-08 1.72e-05 0.80 1.2\n", " 6 -7.850647511580 -1.08e-09 7.09e-06 0.80 4.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+gvaeTAAAgAElEQVR4nO3dd0AU1/428HNmZtmlKggI0sSOEVHBgqKCqKCCilhzbzTNmKbEX8o1zXhveu7NzVWTGEv01SQmMRYUFaXZRUXA3gtCVKyosLSd8v6xhhBE3YVlZ2f3+fy1DLMzX0bk2XPmzDlUkiQCAABgqxi5CwAAAJATghAAAGwaghAAAGwaghAAAGwaghAAAGwaghAAAGwaghAAAGwaghAAAGwaghAAAGwaghAAAGyamYJw1apV/fv379+/f1ZW1sP2qa6uPnr0qOHHFATBFKXZEFwxY4miKHcJCoMrZixcMWOJomjymUGpGeYavX79ekRERG5urlarjYiIOHbsmEajeXC3goKCqKioixcvGnjY0tJSZ2dnk1Zq5XDFjFVWVubk5CR3FUqi1WodHBwopXIXohjl5eUajYZh0DlnqIqKCjs7O5ZlTXhMc1z99PT0oUOHOjs7e3l5BQcH79+/3wwnBQAAMIQ5gvDatWseHh761y1btrx69aoZTgoAAGAIcwShk5NTZWWl/nVFRQV65wAAwHI0KgivXr36008/vfHGG//+979rby8tLX3nnXdGjBjxxhtv3L59u2vXrgcPHtR/Ky8vLzg4uDEnBQAAMCGuMW9OSUlZv369TqerqKh48803a7ZPnjxZFMX/+7//W7p0aWJi4rZt29Rq9YwZM+7cudOvXz9/f/8Gn7FUR1QM0Txwl/RGJfGoZ/wNAADAY5hg1OjixYtXrFixa9cu/ZcXL14MCgq6evWqq6trZWWlp6fnjh07unTpkp2dbWdn17t374eNKCsoKAgPD//ggw9qb3RxcZkwYULNlxsKpcWnpd8GMRr2zzGQS05L24rJTwMxUO0xMGrUWBg1aiyMGjUWRo0ay9hRowzDPPYXslEtwnrl5uZ26tTJ1dWVEKLRaMLCwg4ePNi9e/cBAwY89r06nS4nJ6f2Fnd399GjR9d8GdOSXNcyCeniz/0Fobq6qqpq+XlmfRFd2V+oqjL5j2Jtqqurq3CZjFFVVaVSqeSuQkmqqqpYlkUQGq6qqopSiiA0XFVVlSRJhgehRqORIQivXbvm5uZW82WLFi2Ki4sNfK+zs/P333//6H2mdSF2duLf97LLe9v/XKTZeEVMjuEe7CyFBwmC4ODgIHcVSiKKIq6YUSRJQovQWGgRGoVSavLnCE0fhI6OjjVjRAkh5eXlJu9ceqYDIxEyOEvdQiNuieXUSEEAAGgo038M8fPzKygoqLn1WFBQ4OfnZ/KzCCKp4Mnea9J/jog6TFEEAAANZfoW4YABAyRJ2rx584gRI7Kzsy9fvhwbG2vaU3x/Wky+JO6NqfqhyOGLI8Lys+KCCDa6FXpjAMB8MjIyFi1a1MiDCIJgyGgOIITExMQ899xzTXHkRgXhli1bnnzyyaqqqurqajc3t/j4+OXLl6tUqm+//XbKlCmdOnU6efLkvHnzTNs1+v1pcW2BuGYwpysn/xfMuqrp1yeEZ3cIA73pl31YPEQBAOZx4MABnucnTZokdyE2Ye/evenp6ZYYhNHR0efPn6/50s7OTv8iISEhOjr67Nmzbdq00Q8fNZXkS2LyJXHtYE7NEh0h5I/7hZuLpFaO5InVuo/C2KmdGHy4AgAzCAoKGjdunNxV2ARRFNetW9dEB29UEKpUqoflnIuLS2hoaGMOXq/BrZhhvkyd0THPdmBGBZAWajKxDfPSHuHHc+KCfuwTrkhDAAB4PIWN2XVSkXrHiLZQE0JItxZ0Tzz3fEcmejOflC1oeTNXBwAAyqOwIHwshpLJ7ZlDY1QlVSRkLb/l9yZfbREAABTN9KNGLYGXPVkRyW67Kr20W2jfjCzox/o6oqcUAADqYW0twtqivGl+AhfqTkOT+bnHRAGNQwCQ28/nxRVnH/Psc/pl6atjeD7afKw5CAkh9hyZ04PdHcdtLBJ7JvM5NxCGACCnMa2Z1Rel5Q/PwvTL0uxc4al2Vv7H2aLYxLVu34ymDeNe68LEp/HTdgv3dHIXBAC2Ss2S36LZNQ/JQn0Kpgzl3Bv0SHRAQMCtW7caWyIhhBBRFDMyMgRBaORxqqqq3N3dHz3df0hISO0n8czPJoKQEEIJmdyeOTlWpWHJE6v5x3ZNAAA0kYdlYSNTkPyxMoMJSiSkurp6yJAhjV+vhuO42bNnc9yjxqNUV1eLopx/k61zsMzDuKrJ3HB2XKD04h5h1QXx675sa2cMogEAc9Nn4bhMgRBxSnuGmCIF65Ak6ddffz106FDbtm0nT56sVqv12zMyMnbv3u3p6Tlp0iT9g+Dbtm3z8PDo0qULIeT333/PyclJSEhITk4mhCxdutTOzm7kyJFeXl6GnLS8vHzhwoU3b94MDw+Pi4vTb9Ro7v9IycnJPXr02L59+9mzZ6OjoyMjI03zozaarbQIa4vwovkJ3BAfptd6fk6eUI3GIQCYXe12oclTkBAyffr0RYsWtWvXbv369SNHjtRv/Oijj1599VVPT8/Dhw9369bt9u3bhJBly5ZlZWXpdzhx4sSnn37asDNWVFT079//1KlT/v7+b7311meffUYI0el006ZN0+l0hJDPP/984sSJZ86cadasWWJi4vbt2xv/Y5qEbbUIa6gYktSFGRlAX9krhCXz3/Vj+7ZE0xAAGmvEVt6oz9aiRJKyRV4kYR500jbjJgH5LZprblf/ty5durRy5cpLly45OztPmTIlMDBw9+7dISEhH3/8cV5eXlBQECEkPj7+22+/fe+99+o9gn5F9GeffdbwJTlXrVrl6uq6cOFCQkhERESfPn1mzpxZZ59hw4a9//77hJDbt28nJydbSKPQRoNQL9CZbo7hUgrFJ7cJA7zof/uwJvw4BgA26PVgVjTmJl3uTemKVtRwpI8nHdzKuC46h4f//T569GhQUJCzszMhRKVShYaGHjt2zMnJydHRUZ+ChJD+/fvn5eUZdUZCyPXr119//XVCCMdxy5Ytq/2tI0eO9OnTR//6iSeeoJQWFBQEBATU3qdbt276Fz4+Pnv37jX27E3EpoNQL96fGdSK+TBf6Iw5uwGgcQYZsx5c+mUp+ZK4K55z5EhCBh/UXNLfL2w8rVZbe4H0iooKBwcHBweHyspKSZL0qz6Vl5frW3ssy9aMDq2oqHj0kZ2cnBITEwkhDFO3VK1WW/NaFMXq6uoHW5OmXVneVGzxHuGDHDnyWU82bRj3/WkxciN/4g4eNwSAplX7vqA9R9YP4R72TEXDHD58+MSJE4SQoqKiffv29evXLzAw0M3Nbe3atYSQioqK3377Td8z2bp1a33TUJKkNWvW6N+uVqvt7e31NxFrc3BwGD169OjRo2vuO9a2YcOG8vJyQsiaNWv8/Px8fHxM9eM0KbQI/9StBc0eyS05LQ5I4Z/vxMzpwWos8bMLACjeg6NjHhxH2khdu3Z95plnWrZseeDAgXfffbdt27aEkGXLlk2ePHnJkiVnzpwJDw//29/+Rgh57rnnwsPDBw4cWFZW1rlzZ/3bKaWvvvpqaGion5/fkiVLaro0H619+/Z9+/b19fU9cODAL7/88mCr0TIhCP+CoeSFTky8P/OPA0KXNfw3fdkYX3SUAoApPWyMqGmz0N3dfePGjYcOHfL19W3VqpV+Y3R09Llz506ePOnh4eHn56ff6Ovre+rUqRMnTvj6+rq7u9f0qX7xxRcfffSRVqvV32s0RERExEsvvXTmzJmgoCAXFxdCiEajuXbtmv4JirS0tJpHKZ5//vmnn366kT+jqSAI6+HtQFZEsllXpJf2CB0wZzcAmE7WFWl2rrAxhmuhrue7apb8OohNyOA5Sv7W6FnW7OzsevXqVWejvb19jx496mx0dHTs2bOn/rVKpap9hJoV1w3k6urau3fv2ls8PT31L2oHqlqtrnm0UXbKaLfKYlAreghzdgOASbV1IQ9LQT17jqwbzPVwb8iHbz8/P/1olKCgoAkTJjS4yIaJiYkJDw9vwBt9fHxqp6/5oUX4KPo5u59sy7y8R1hxTvyuH9vTA01DAGi4AKfH/w2x50hQ84b8qcnJydG/6Nq1a9euXRtwhMZISEho2BszMjJMW4mx0CJ8vA7NaPpwLukJZlQ6n5QtlGLObgAAK4IgNIh+zu7jiapKgQRhzm4AACuCIDSCq5osjGB/jmK/OCLGp/GXynDbEABA8RCERuvvRfMTuMGtmJ7JmLMbAEDxMFimIfRzdscH0Ff2CD2T+e8i2HBPDKIBsC2BgYFvvPHG1q1bG3OQmgnP4NFKSkpqJjI1OQRhw7VxpqmxXEqhODFLGOhFvwpnHzEkGgCszKRJkzp27NjIhXArKyvt7OyUMgOLvHx9fZvoyAjCxor3ZwZ4MbNzha5r+E97Mk+1x5zdALbiwSfTjVVeXq7RaBCE8sLVN4FmdmRuOLsxhv3mhBi1iT+JObsBAJQDQWgy3VvQ7JHck22Z/in8rByhUpC7IAAAMACC0JT0c3YfTVRd0ZLgNXzaZTQNAQAsHYLQ9PRzdv8vnJ22WxifKVx7zDqXAAAgJwRhUxnhR08kcp1dSde1urnHRBGNQwAAi4QgbEL6Obt3xXEbCsVe6/mDNxGGAAAWB0HY5Do0oxnDuRlPMPFbMWc3AIDFQRCaw/05u8fen7N79UVMywYAYCkQhObjpiYLI9iVUewHuZizGwDAUiAIzW2AFz00hotoyfRYx8/JE3RoHAIAyApBKAMVQ/4RwhwYxe27LvVM5rOvo2kIACAbBKFs2rrQLbHch2HMhExh8nbhVpXcBQEA2CQEoczi/ZmjiZyrmnRdg4XvAQBkgCCUn37O7pSh7NcnxKhN/CnM2Q0AYEYIQkvRw53ujedGBzADN/Fz8oQqzNkNAGAWCEILwjEkqQuTn8BduEe6rOHTMWc3AEDTQxBanFYOdEUk+1UfduouYXymcB1zdgMANCUEoYWK86cnxnKdXUkw5uwGAGhKCELL5cCROT3Y9GHcqovigI38sRKEIQCA6SEILV1XN7o7nnuhEzNkM5+ULZRhzm4AAJNCECpAnTm712DObgAA00EQKoZ+zu4fI9n3c8X4NL4Qc3YDAJgCglBhBnrTw2O4iJZM93X854dFAWkIANA4CELlqZmzO+uKGLqO34c5uwEAGgFBqFRtXejWYdy73ZmEdH7yduE25uwGAGgQBKGyjQtkTo5TuapJ17WYsxsAoCEQhIrX3I7MDWfXD2HnHxcHbeJP30VPKQCAERCEViLUnWaP5EYFMP1TMGc3AIAREITWo2bO7uMlJHgtn3FZIoRIhNTbQsScbQAAeghCa+PjSH+LZr/szUzdLUzeLvxwTpy+t+5DFjk3pDEZaDMCABCCILRW8f7M8USujQt5c79w5q70wu4/szDvpjRtt/C/PvinBwAgBEFoxWrm7C7Tka2/i+MzBYmQvJvS87uEtYPZ1s5U7gIBACwCJ3cB0LT0c3YvPCW+sV/odcvOXiVsGMr6OSIFAQDuQxBaP4aSl4KYQGc6LpO3q5LKebkLAgCwJOgatQl5N6V3coR9MVX9WtLu6/gfz+HRewCA+8wUhKdOnVq+fPnSpUvNczqorea+oL+jtGEo91Q7Zka2MHm7oEXTEADAbEG4dOnS/Pz8//3vf+Y5HdTIvSlN3S2sH/rn6JiFEezk9syBG1LPZKx6DwBgrnuEX3zxRWFhYVZWlnlOBzU87UnykLqjY/7Xh90bKJ27J0Vu5D/tyU7thB5yALBdGCxj5R42QLRvS9q3Je3pQcdnCruKpQURrCN+FwDAJpn4j9+QIUNu375de8umTZu8vLxMexYwlaDm9MAoLmmfEJbM/zqI7eqGxyoAwOYYEYQ3btwoLCxs06aNq6trzcaqqqrMzMzKyspBgwY1b948PT29CYqEJmTPkUUR7IqzYvRm/r1ubFIXdJMCgG0x9K9eWFiYv79/3759MzIyajZqtdrw8PDPP/985cqVQUFB58+ff9jbL1y4cPTo0YqKitzc3Fu3bjW2ajC1ye2ZXXHc0jPiU9uFMp3c1QAAmJGhQfjDDz+UlpZ26NCh9sbly5c7ODhs27Zt9erVY8eO/eyzzx729oyMjA0bNgwaNGjRokVnzpxpVMnQNDo1p/tHcW5qEprMH76N0aQAYCsM7RoNCgp6cOOGDRvGjh3LMAwhZPz48YmJiYsXL6737S+88IIhZ7l58yalf7lNFRgYeOTIkXp31mq1dXaGRzPkin0cTHo2Z4Zu5l4PEl7uaOsrVGi1WrlLUJjy8nJRFPEf03Dl5eU8z+v/ioIhKioq7OzsWJY1cH8HB4fHXt5GDZa5fPmyr6+v/rWfn9/NmzcrKys1Gk2DD+ju7l5aWmrgzpIkOTk5NfhcNsjAK/b3INKrlTQ+k9lfYvf9ALa5nRlKs1z4HTMKpdTBwQFBaDiGYTQaDYLQcCzLGhWEhmjU1ed5vqYalmUlSeJ5zFZiDTo0o/tGcb6OpPs6fv91dJMCgDVrVBB6eXnduHFD//r69esuLi74+Gw1NCyZG85+2ZsZlc7PPYa5SQHAajUqCAcMGFAziDQjI2PAgAGmKAksyJjWzL6R3M8XxIR0oaRK7moAAJqAofcIlyxZcu7cueLi4p9++ik3N3fatGmBgYHTpk0LCQl58803vby8Pv3005SUlCatFWTR2pnuGMG9dUDovo7/ZRDbxxO3fwDAqhgahN7e3qIofvzxx/ov9SNiWrVqlZOTs2zZsuLi4vT09J49ezZVmSArNUvmhrNRrcSRafzrwexbIQzCEACsBpUkSxkKUVBQEBUVdfHiRQP3Ly0tdXZ2btKSrEzjr9ilMmlSluBhT5cNYN3UpqrLcpWVleG2t1G0Wi1GjRqlvLwco0aNYuzjE4bA1QcjBDjRnXFcUHPSfR2/95qlfIQCAGgMBCEYh2PIZz3Z+eFMQgY/J08QkYYAoHAIQmiIkQHMwdFc+mVpVDp/C6NJAUDJEITQQH6OdMcILtSd9ljH7y5GwxAAlApBCA3HMWROD/abvuzYTHSTAoBSIQihseL86cHRXOYVaWgqX1whdzUAAEZCEIIJ+DrSbcO5CC/aY50u4zIahgCgJAhCMA19N+kPkdzTOwV0kwKAgiAIwZSiW9Hc0dyea9LgzfzVcrmrAQAwAIIQTKylPdkSyw3wpj3W6dLQTQoAFg9BCKbHUjKnB7tyEPfMDmFWjiAgDQHAgiEIoalEedP8BC7/pjR4M3+lHGEIABYKQQhNyNOebBnGxfoyoev4Lb8jCwHAEiEIoWlRQv4Rwvw8iHt+l5CULeiw1j0AWBgEIZhDpDfNT+DO3JUGb+Yva9E0BAALgiAEM/HQkM2x3JjWTM/1/OYiZCEAWAoEIZgPJSSpC/PLIG7abnSTAoClQBCCuQ3wovkJ3Ll7UkQKf7EUTUMAkBmCEGTgriEbY7gn2zK91/NrLqJhCAByQhCCPPTdpBtjuLcOiEnZQjXSEABkgiAEOfXyoHkJ3JVy0m8DfwHdpAAgBwQhyKyZHVkVzf69HdN7Pf8bukkBwOwQhCA/fTdpaiw364A4bTe6SQHArBCEYCnC3Gn+GK6kivTdwJ+/h25SADATBCFYEBcVWRXNvtCJ6bOB/+U8GoYAYA4IQrA4L3RitsRy7+eKk7cLFbzc1QCAtUMQgiUKdae5CZxOIv1S+LN30U0KAE0IQQgWykVFfo5iX+vCRGzkf0Y3KQA0GQQhWLTJ7ZktsdwHeeLk7UI5ukkBoAkgCMHSdW9Bc0dzvER6JvPHStBNCgAmhiAEBXBWkZVR7D9CmOjN/I/n0E0KAKaEIATFmNye2TaC+/ywOHm7oEU3KQCYCIIQlKRzc5o9kpMICUvmj95GNykAmACCEBTGSUV+iGTfDmEGbebnHkM3KQA0FoIQFGlye2ZXHPf9GXHydqFMJ3c1AKBkCEJQqk7N6f6RnKuahCXzR9BNCgANhSAEBbPnyNxw9t3uTDS6SQGgoRCEoHhPtWN2x3NLz4hjM4W71XJXAwBKgyAEa9CxGd0/ivNxIL3W84duoZsUAIyAIAQroWHJ3HD2kzAmdgu6SQHACAhCsCqJgczOOO7/nRXHZAh30E0KAAZAEIK16dCM7hvJ+TmS7uv4fdfRTQoAj4EgBCukZsnccPa/vZnR6fzcYyLCEAAeAUEIViuhNbN/FPfLBTEhXSipkrsaALBUCEKwZgFOdPsILsCJdF/HZ6ObFADqgyAEK6fvJp0bzoxO5+fkCegnBYA6EIRgE0YFMDmjuLTfpdHpwm10kwJALQhCsBX+TnRnHNfDnfRYx++5hoYhANyHIAQbwjFkTg92fl9mTAa6SQHgPgQh2Jx4f+bgaC7jshSzhb9WIXc1ACA3BCHYIj9Hun0E168l7bWe31WMhiGATUMQgo3Sd5N+35+dtE1ANymALUMQgk0b7EP3jWSzrkhDUvlidJMC2CQEIdg6X0eaNZzr70V7rNOlX0bDEMDmIAgB7neT/hjJPbNTmLJDmJ0r1NmhsEx6clvdjQBgHRCEAPcNakXzRnNXy6Wlp8VX9vwZe0VaKSFDeK0L/rMAWCf83wb4k6c9SY3lnuvE/L+z4vhMgRBSpJVGpwsL+rG9PKjc1QFAk+DkLgDAsrCU/LMH29+LGZ3GD7yrkhjhO6QggFVDixCgHoNb0W0juFOlTGGZJGEADYBVQxAC1KNIK724R9gUpevlQYem8pO3Czcq5a4JAJqGmYJQEITTp0+fP3/ePKcDaIya+4JhbuKmGG58G+bMXemJ1bq5x0Q8dw9gfcxxj7Cqqqpbt27BwcGlpaWiKKakpNjZ2ZnhvAANUHt0TFkZoYQs6s++skcIcqWrLoo/nhO/7cf2xC1DACtijhahSqU6ePDgqlWrUlNTCSE7duwww0kBGqacJwsj/jI6hhLydV82xofuiuemP8GMSueTsoV7OhlrBABTMkcQMgzj6OhICBEEobi42MfHxwwnBWiYjs1omHvdBh9DycS2DEPI5PbM8UQVIaTzan7FWVGOAgHAxEzcNbpw4cIbN27U3jJ16tSWLVvqX7/zzjuRkZGdO3c27UkBzMlVTeaGs5PbSy/tEZadEb/px3Zujp5SAAUzcRC2a9fO09Oz9haNRqN/8fHHH1++fHn58uWmPSOALELd6d547psT4oAU/m/tmI/DWCeV3DUBQIMYGoTz5s1LT08/ffr0nDlznnzyyZrt33777b/+9a+Kiorhw4cvWbIkOjq63rd/9dVX+fn5v/zyC8uyJqgawAJwDEnqwoxvw/zjgBC0mv+qDzM2EM8jASiPof9vdTrd+PHjHR0d7927V7Px0KFD7777bnp6+tWrV0tKSj788MN633v79u0PP/zw3Llzffr0CQsLS0lJMUHhAJbB24GsiGR/iGQ/yBXj0/iCUjxgAaAw1KhpM6KioiZMmPDiiy/qv0xKSqqoqFi0aBEhZNu2bRMmTLh+/XqDSykoKAgNDY2MjKy90dPT88svv6x3/7KyMicnpwafzgbhihnLqCumE8miM8wXx5kXOohvdBbVNtn3odVqHRwcKMVNU0OVl5drNBqGQV+CoSoqKuzs7AzvXDTk8jbqHuGZM2diY2P1r4ODg2/cuFFSUuLq6trgA6rV6vHjx9fe4uTkpFar6925urr6Yd+CeuGKGUun0xl+xdSE/F8IGd9Oev0A03sLmdeHDmnVpNVZIp7n1Wo1gtBwgiCo1WoEoeFEUTQqCA35bWxUEJaUlDg7O+tf61/cvn27kUE4YcIEA3dmWRZ3HI2CK2asBlyxABeyejBJKRRf2it2cSUL+rG+jjaUCvorhiA0nP6KIQgNx/7BhMds1NVv0aJFzS3Du3fvEkI8PDxMUBSAwsX7MycSuVB3GprMzz0mCrhvCGDBGhWEHTt2PHLkiP71kSNHvL29XVxcTFEVgOLZc2ROD3Z3HLe5SAxL5rOvIwwBLJShQXjhwoXc3NzS0tLCwsLc3Nw7d+4QQp599tm1a9dmZ2ffvn37o48+ev7555uyVADlad+Mbh3G/SuUmZglYAkLAMtkaBAuXrx42rRphJC0tLRp06bl5+cTQrp06TJ//vynn366c+fOnTp1euedd5qwUgDFivdnTozlWjkSLGEBYIGMe3yiSRUUFERFRV28eNHA/UtLS2uG6oAhcMWMZfIHTg7fll7eI1QLZEEE++CMplYAj08YC49PGMvYxycMgasPYD4hbnR3PDf9CWZkGpawALAUCEIAs6JYwgLAwiAIAWSgX8IieQj79QkxahN/4o6l3KEAsEEIQgDZhLnTvfHc6ABmQAqflC2UoacUQA4IQgA56ZewOJqoKqkiIWv5TUVoGgKYG4IQQH76JSy+H8C+tV/AEhYAZoYgBLAUkd700BhucCum13p+Tp5QJchdEIBtQBACWBAVQ5K6MHkJ3PESEryWT7uMpiFAk0MQAlgcX0f6WzT7ZW9m2m4hPo3/XYs4BGhCCEIAC4UlLADMA0EIYLmwhAWAGSAIASxdnSUsbmIJCwCTQhACKEO8P3NkDOeqJp2xhAWASSEIARSjmR2ZG86mD+d+vSD2Xs8fvIkwBDABBCGAwoS40T0jsYQFgMkgCAGUB0tYAJgQghBAqWovYTFoE38SS1gANAiCEEDZ9EtYjApg+qfwSdmClpe7IAClQRACKF7tJSy6rsESFgDGQRACWAksYQHQMAhCAKuCJSwAjIUgBLA2WMICwCgIQgDrhCUsAAyEIASwZljCAuCxEIQAVg5LWAA8GoIQwCZgCQuAh0EQAtgQLGEB8CAEIYBt0S9hkTYMS1gA3IcgBLBF3VpgCQuA+xCEADZKv4TFMSxhATYPQQhg09z+WMJi/nEsYQE2CkEIACTMnWaPxBIWYKMQhABACJawABuGIASAP2EJC7BBCEIAqAtLWIBNQRACQD1qL2HRdS2fjiUswHohCAHgofRLWPynN/MClrAA64UgBIDHiPdnjmMJC7BeCMj/j6YAABkOSURBVEIAeDwHLGEB1gtBCACGwhIWYJUQhABgnDpLWOy7Ls0/Xnd6tns6Mn2vgNUtQBE4uQsAAOXRL2HxTAfm5T3CbwWkmYqKkpjU5f4Hay1PEjP4l4MYhspbJoBBEIQA0ED6JSx+OCvOyhHO36OVAnm1HdHyJCGDfzmISWiNDidQBgQhADScfgmLOH9mTp7wr3wx5xp3VxCQgqAs+GUFgMZyU5N54exv0WzK78y1cgkpCMqC31cAMAEtT746JnzTm/9dS0dswZxsoCQIQgBoLC1PRqfzLwcxk1qLx8eyu69JCenIQlAMBCEANEpNCup7RL3syZFENv2yNCkLWQjKgCAEgEa5cE96tfNfRscEONGcUWzq79KqC3WfLwSwQBg1CgCNEuxGg93qPjAY5Ep3xrFDU3lnFR3mh8cJwaKhRQgATaKrG00ewj29k99VjAlmwKIhCAGgqfTxpCujuLGZfN5NZCFYLgQhADSh6FZ0YQQbl8afvIMsBAuFe4QA0LRGBzClOhK7Rdgxgm3tjPuFYHEQhADQ5J5qx9yrJkNShV3xnJe93NUA/BW6RgHAHF7pzDzVnhmayt+ukrsUgL9CEAKAmczuzsT60uFb+TKd3KUA1IIgBADz+bwXG+JGR6XzlZh2BiwGghAAzIcSsqAf62lPJ2YJPKadAcuAIAQAs2IoWTGQ5UXpmZ2CiEcqwAIgCAHA3FQM+S2aK9JK07PRQwryM0cQ6nS68ePHh4WFDRgw4NdffzXDGQHAwtlzZMNQ7sB16f1cZCHIzBzPEbIs+69//atTp05Xrlzp3bt3RESEj4+PGc4LAJbMRUVSY7mBG3lnlfhWV/ROgWzM8cvHMEynTp0IIa6urk5OTgyD33gAIIQQdw1JG8YuPCkuPIWRMyAbE7cIL126VGdLQEAAIaSqqqpfv37nz5+fNWuWt7e3aU8KAMrl40jThrGRm4RmKjKxLT4lgwxMGYQ6nW7mzJl1Nq5du5YQolarDx48eOPGjejo6GHDhnXt2tWE5wUARWvrQlNj2SGbeWc7OgKLF4LZGRqER48ezc3NvXbt2lNPPdWqVaua7fv371+xYgXDME8//XRoaKg+9h7Gw8MjJCTkwoULCEIAqK2LK10/lItP438dxEV6IwvBrAzqiKiurh4xYkRycvLs2bOLiopqtufm5g4ZMqRjx46BgYGDBg06duxYvW/Pz89//fXXFy9e/Pbbb+fm5kZFRZmmdgCwIr086C+DuAlZ/EEsXgjmZVCL0M7OrrCwkBDSrFmz2tu/+uqrV155ZcaMGYSQwsLCefPmLVq06MG3BwUFDR069NKlS7169Zo9e7a9/UMnn9fpdBkZGbW32Nvb9+vXz5AiAUDporzpkv5s/FY+Yzj3hCvahWAmjbpHuHfv3u+++07/Oioq6t133613N41GExMTY8gBy8vLP/nkk9pbvL29Q0JC6t1Zq9VSiv8qRsAVM5ZWq5W7BIUpLy8XRbExv2ZRbuQ/PZjYVJIaXd3a0fqbhuXl5TzPYyy94SoqKuzs7FiWNXB/BweHx17eRgVhcXGxu7u7/rWHh8eVK1caczRCSLNmzbKysgzcWZIkJyenRp7RpuCKNQCumFEopQ4ODo38vPW3IFLNiqN3qHfGsa0crPyjG8MwGo0GQWg4lmWNCkJDNOrqq9Xq6upq/evq6upH9HkCABjumQ7Mq52ZqE3C9Qq5SwEb0Kgg9PHxqRk7U1RUhPliAMBUXuvCJLamMVv4O9VylwLWrlFBOGbMmJUrVxJCJElauXJlQkKCiaoCACCf9GSjvOmIrbyWl7sUsGqGBmFCQkJYWJhWq3366afDwsL0twOnT59++vTpgQMH9u/f/+rVqy+99FJTlgoANufLPmynZnRMOl+FqbmhyRg6WGbevHk6na7mS09PT0KIh4dHfn7+3r17GYbp27evSqVqkhoBwFZRQhb1ZydmCU9uE1ZFs6yVD50BeRgahH5+fvVuV6vVeEAeAJoOS8lPUezINP75XcLSAYhCMD2M2QUAS2fHkLWDuXP3pJn70EMKpocgBAAFcODIxqHczqvSx4ewYBOYGIIQAJShmR3ZEsv9cFb871FkIZgSghAAFMPTnqQPZ+efEL8/jSwEkzHxwrwAAE3Kz5GmxbKRmwQXOzIuEB/lwQQQhACgMO2b0dRYdmgq76yisb4YRgqNhc9TAKA8Xd1o8hBuyg5+d7H1r1ABTQ1BCACK1MeT/hTJJWbyeVjIFxoHQQgASjXYhy6MYOPS+FN3kIXQcLhHCAAKNjqAuVdNYrYIO0awrZ1xvxAaAkEIAMo2uT1zT0eGbhF2xnFeWBQVjIeuUQBQvFc7M0+2pTGp/O0quUsBBUIQAoA1mNODHepDR2zly3SP3xmgNgQhAFiJL3qzwW50NBYvBCMhCAHASlBCFvRj3TV0YpbAYwo2MBiCEACsB0vJD5FstSg9u1MQ8UgFGAZBCABWRcWQ1dHcpTJpRjZ6SMEgCEIAsDb2HEmJ4fZdlz7IRRbC4yEIAcAKuajIllhu9UXp30dwtxAeA0EIANbJXUPShrHfnRQXnUIWwqNgZhkAsFo+jjRtGBu5SXBRkYlt8bkf6ocgBABr1taFpsayQzbzznZ0hB8mI4V64CMSAFi5Lq50/VDu2Z38jqt4ogLqgSAEAOvXy4P+Mogbn8UfxOKF8AAEIQDYhChvujiCjd/Kn8DihfBXuEcIALZiZABTqiOxqcKOODYQixfCHxCEAGBD/taOuacjQ1KFnXFsKwdkIRCCrlEAsDUvBTHPdWRiUoVbWLwQCCEIQgCwQW+HMHH+dPgWvhSLFwKCEABs06c92b4t6ag0vhLTkdo8BCEA2Kj/9mFbO9NRaXw1pmCzbQhCALBRlJDF/VlnO/rkNkHAIxU2DEEIALaLpWRlFFumk6buQhTaLgQhANg0O4asHcyduSvN3Ie7hTYKQQgAts6BIxtjuB1XpU8P426hLUIQAgCQ5nZkayy3/Iz41TFkoc1BEAIAEEKIpz1JH87OOy4uPYMstC2YYg0A4D4/R5oWe38h37GBaCfYCgQhAMCf2jejm2PZmFTeWUVjfDEZqU3ARx4AgL8IcaPrhnCTd/C7i/FIhU1AEAIA1BXuSX+M5BIz+fxbyELrhyAEAKjHEB/6XT92xFb+FBbytXa4RwgAUL+E1sw9HYnZIuyMYwOccL/QaiEIAQAeakp75l41GZIq7IrjWtrLXQ00DXSNAgA8yvQnmElt6NBUvgQL+VopBCEAwGP8M5Qd4kOHb+XLsJCvNUIQAgA83r97s11caUIGX4Wpua0OghAA4PEoId9FsG5qOjFL4DEFm3VBEAIAGISl5IdItkqUntsliHikwoogCAEADGXHkNXRXEGplJSNHlLrgSAEADCCA0dSYrjs69KcPGShlUAQAgAYx0VFUmO5VRek/xzF3UJrgCAEADCah4akD2O/PSEuPoUsVDzMLAMA0BA+jjRtGBu5SXCxIxPaoFGhYAhCAIAGaudCU2PYoam8s4oO98NkpEqFTzEAAA0X7EaTh3DP7OR3YvFCxUIQAgA0Sm9P+nMUNy6Tz72JLFQkBCEAQGMNakUXR7DxafwJLF6oQLhHCABgAiMDmHs6Epsq7IhjA51xv1BJEIQAAKbx93Z/Ll7o7SB3NWAw83WNSpKUnJxcUFBgtjMCAJjZy52ZZzswQ1P5W1i8UDnMF4TLli178cUXd+zYYbYzAgCY3zvdmBH+dPgWvhSLFyqEmYLwypUrq1evTkxMNM/pAABk9GlPtnsLOiqNr8R0pEpgpiCcOXPmf/7zH5ZlzXM6AAAZUUK+7cd6OdDxmYIOU7BZPFMOlsnLy3vttddqbwkODv7mm29+/PHHLl26dO7c2YTnAgCwZAwlyweyYzL4Z3YKKwayDIaRWjAqSYY+9XLnzp2SkpLAwMDaG8vKynJzcz08PDp37iyKYkVFRe3vsiyr0Wiio6Nv3rypUqmKioocHBzmzZsXHx//4PELCgqioqIuXrxoYD2lpaXOzs4G7gwEV8x4ZWVlTk5OclehJFqt1sHBgVL81b+vgifDt/JBzem3/ervDysvL9doNAyDR7oNVVFRYWdnZ9r+RYNahIcPHx47duzFixcppTrdn/d/8/Lyhg8fHhwcfPbs2ejo6O+//97R0fHBt2dmZupfzJgxIzQ0tN4UBACwPvYc2RjDDU3lZ+4TvuqDe0MWyqCPId7e3j///PP+/fvrbH/77benT5+enp6en5+/efPmvXv3Pvo43bp1a926dcMKBQBQIkeObIrhtl+VPjuMu4UWyoiu0aNHj/bo0aOmRXj37l1XV9fCwkJfX19CyAsvvODs7Pzll182uJSCgoKuXbvWafC2bt16165d9e6Pbitj4YoZS6vV1tvJAQ9TXl5ub2+PrtEH3aikMZmq59oJr3T8y0BSdI0ay9iuUQcHh8de3oYPlrl8+TLLsj4+PvovAwICDh8+3OCj6bm5ueXn59fewjDMw/52S5KEP+tGwRVrAFwxo1BKcY+wXk5OJGOENGAj4+msfqbDn3+XGYZBEBqFZVl57hHWq7KyUqVS1fzGq9XqOiNlGoBS6urq2siDAABYIH8nmjaMjdokuKhIYiCSz4I0/B+jZcuWFRUV5eXl+i9v3rzp5eVloqoAAKxQh2Z0cyz7yl4h5VI996TKefNXBIQ0Jgi9vb0DAgJqpkzbuXNnnz59TFQVAIB1CnGjn/dkEzP59Mt/ycLcm9LgzTyP8TRyMKhrtLy8fP78+cXFxaIofv75587Ozi+//DLDMDNnzpwxY8Ynn3xy4MCBoqKiiRMnNnW5AABKN6UDU1RO4rbymSO4Hs6EEJJ/S5q6S1gzmOXQYyoHg4JQFMWSkhK1Wv3mm2+WlJQIwv1RTzNmzHBzc0tOTvb09Ny9ezfG1wEAGOK9bkwFT6I38dtiqNpOmrpbXDMYqxjKxojHJ5oaZpZparhixsIDJ8bCzDJGeXmvsPS06OtI02LZNi64aAZpipll0A4HAJDHcx2Y5irpWoUUup4fmyl8f1q8Um4pLRObghXqAQBkoL8vmDWUv1Cp+eywODKA7imW3jsoOHA0zp/G+zMDvakKTRWzQBACAJhb/i3puZ3CmsFsS5Z0cieUsv85Iqwfyi3ox+bfkjIuS3PyhCO3pYHeNN6fGeFHfRzRcdqEEIQAAGZ14o70/C4heQjr70T1T2KP8KNVAjMqjU8fxoW601B3+o8Q5kYl2X5VTLkkzcoRvO1pfAAd3ArNxCaBIAQAMKu2zjRlKNvK4S+NvDGtme4taO3HJzw0ZFwgMy6QCBJ76JaUUijOyhEulUlR3sxgHxrnT+scARoMQQgAYFZqltSbYQ97fIKlJNSdhrqzc3qQ6xVky+/ixkLp7Rwh0JkO9qFxfkzflhQL/zYGghAAQDE87cnk9szk9n82E1/b92czcWQA42Uvd4kKhCAEAFCe2s3EaxVk6+/ixkJpVo6uDZqJxkMQAgAoW8s/mom8yO67Lm0sEl/bJxSWSZHejP5JDFe13CVaNgQhAICV4BgS4UUjvFjSk1wsldIvSxsLpaRsXVuX+88m9nDHrD/1QBACAFihQGf6Qif6QidSKbC7i6WMK+Lk7cLNKinGh4kPoEN8mOZ2cpdoMRCEAADWTMOSwT50sA/7WU9yoVTKuCz9dkGaukvXDs3EPyAIAQBsRZs/mokVPLvnmpRxRfz7dqGkShrqw8QH0BhfxkUld4lyQBACANgce65uM3HFWfGZHUK3FjTenxnsQ0PdbaiViCAEALBpfzQTmXKe7L0mpRSKiRkiQ8kQHzrYh8b6Ms7W3kxEEAIAACGEOPzRTJwbTi6USimXpEWnxKm7hJ4edHArJj6Adm5unc1EBCEAANTVxpkmdaFJXf5sJsamiipGn5R0mC/jZEXNRAQhAAA81KObiaMCaCflNxMRhAAAYJCaZqKWJ1lXxI2F0pBU0Y4h+tUwhvowalbuEhsEQQgAAMZx5Ei8PxPvTwghx0ukjYXSvOPiU9uFnh40zo8Z3ZoGOCmpmYggBACAhnvClT7hSv8RwtyuIplXxIzL0hcbRA1L9E/rD/Cmdha/kjCCEAAATMBNfX8lYfJHM/HzI8LYTEnfTBwTSP0cLbSZiCAEAAATq2km3qoiWVfEjMvSZ8mCA0cts5mIIAQAgKbS4o9m4oJ+bP4tKeOyNCdPOHJbGuhN4/2Z4X7U1wKaiQhCAABocsz9lYTpP0KYG5Vk+1Ux5ZI0K0fwtqfxAXRwK2agN1XJ1Ey0pNYpAADYAA8NGRfIrIhkb/xdtSKS1bBkVo7g/ZNufKaw6JR4pVyqs//pu9KH+WKdjZUCeWWvUCWYoB4EIQAAyIOlJNSdzunBHhzNnRirGteG7i6WgtfwYcn8rBwh47LEi4QQ0qEZLa6Q3jrwZ+hVi2R8phDsSk3y5CK6RgEAQH6e9vfvJgoSe+iWlFIozsoRLpVJUd7MYB/6Xnf2o3zhrQPCP4NJtUgmbBOG+9EXg0zTlkMQAgCABWHv301k5/QgRVoptUhKLZLeOqDr3JwKhDxbxpTppBH+jKlSkCAIAQDAYvk53l9JuFpkdxVLKYXi96eYQd6SCVOQ4B4hAABYPjuG9PeiF+6Rj3uIPk609v3CxkMQAgCApasWydgMYbgfndpenNeHanliwixEEAIAgEWrSUF9jygl5Ou+rAmzEEEIAAAW7bJWmtDmL2NE9VnY0p6a5DlCDJYBAACLFuhMA53rzsRGCXk92DRtOQW3CJOSki5cuCB3FUry5ptvnjx5Uu4qlOS9997Lz8+Xuwol+fDDD/ft2yd3FUryxRdfbN++Xe4qlGTu3Llbt2417TEVHIQ5OTllZWVyV6EkeXl5d+/elbsKJTl8+PDt27flrkJJjh49euPGDbmrUJITJ05cu3ZN7iqU5NSpU1euXDHtMRUchAAAAI2HIAQAAJtGJanuPN9yuXXrVlxc3PXr1w3cv6SkxNnZmeMw3sdQuGLGunPnjqOjo0qlkrsQxbh79669vb2dnZ3chSjGvXv31Gq1Wq2WuxDFKC0tValUGo3GwP03btwYFBT06H0sKAgJITdu3CgtLZW7CgAAsBK+vr6P/WRmWUEIAABgZrhHCAAANg1BCAAANg1BCAAANg1BCAAANk2RQfjpp5/GxMS0a9du48aNcteiDK+++mrr1q1VKlWHDh1+/PFHuctRgKlTp7Zs2ZJhmMDAwK+//lruchTj5s2bISEhzzzzjNyFKMDs2bPb1iJ3OcqwYsWK9u3bU0oDAwP3799vqsMqMggJIc8++yzDMJhizUAtWrRIS0urrKycO3fu1KlTDx8+LHdFlm7KlCknT54URfGHH36YNWtWdna23BUpQ1JSkpOTU3FxsdyFKMDNmzcTEhLS/yB3OQqwbt26t99+e9myZTzPZ2Vl+fr6murIigzCt99+e8KECQ4ODnIXohj//Oc/O3TowLLssGHDOnbsiCB8rIiICDc3N/2Ltm3bXrp0Se6KFGDjxo13795NTEyUuxDFcHV1bfMHuWtRgE8++eSDDz6IiIhgWTYwMNDHx8dUR1ZkEEKDXbx48ezZs7169ZK7EAXIy8v79ddfZ86cqdFo4uLi5C7H0t29e/ett95asGCB3IUoyfz58729vfv167dhwwa5a7F0oigeO3bs3r17YWFh7du3nzVrVnV1takOjtm2bIhWq50wYcJrr73WqVMnuWtRgLy8vC1bthw5ciQ6Ohrz0j3WzJkzX331VT8/P7kLUYynnnoqKSmpefPmW7ZsmTBhQmZmZt++feUuynKVlJRUVlauXr16/fr1giCMHDnSycnpvffeM8nB0SK0FRUVFSNHjgwODv7oo4/krkUZnn/++dWrVx8/fjwnJ2fhwoVyl2PR9u7du3Pnzp49e+bm5v7+++937949dOiQ3EVZuvDw8I4dO7Zs2XLKlCnjxo1bu3at3BVZtObNm3Mcl5SU5OPj4+/vP2PGDBMOlsTnXJtQXV09btw4d3f3RYsWUVp3oWd4BJVK1bVr16KiIrkLsWhVVVWBgYHvvPMOIaSwsPDWrVuzZ89Gd5/h8L/ysViW7dChQxMdXJFBePbs2Xv37pWXl1+8eDE3N7dDhw7Ozs5yF2W5JEkaM2bM1atXFyxYoP+c7uPj4+XlJXddlqu0tHTt2rVRUVH29vZ79uxZs2YNPq0/WlRUVFRUlP71f//73/T0dKTgYy1ZsiQ6OtrR0TE9Pf23337DwNHHevHFF+fNmzdo0CBBEL7++uvx48eb6siKnHR7+vTptYezf/fdd2FhYTLWY+F0Ol14eHjtLdOnT58yZYpc9Vg+rVY7ZcqUgwcPVlZWtm3b9rXXXhs3bpzcRSnGTz/9tG/fvvnz58tdiKWbOHHi/v37Kyoq2rVr9+abb44aNUruiiydKIrvvffejz/+qFarJ02a9P7775tqiTRFBiEAAICpYLAMAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYNAQhAADYtP8Pckkh4uFFFj4AAAAASUVORK5CYII=", "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" ], "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" ] }, "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.2" }, "kernelspec": { "name": "julia-1.6", "display_name": "Julia 1.6.2", "language": "julia" } }, "nbformat": 4 }