{ "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.844955072242 NaN 1.97e-01 4.0 \n", " 2 -7.850336839205 -5.38e-03 2.93e-02 1.0 \n", " 3 -7.850617210242 -2.80e-04 2.75e-03 1.3 \n", " 4 -7.850646808108 -2.96e-05 1.25e-03 2.5 \n", " 5 -7.850647346859 -5.39e-07 6.00e-04 1.3 \n", " 6 -7.850647505598 -1.59e-07 5.12e-05 1.0 \n", " 7 -7.850647511555 -5.96e-09 4.92e-06 2.3 \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+gvaeTAAAgAElEQVR4nO3deUBU5d4H8N/znDMzLIKCiLKIK5aJaCoWiAIKKAqaJqbdrNzvvS1Wr2XXlmu9t9fsduvaXlY3uy6VlQqooIAKbuWGu7nggiioKCD7nOX9Y4zIREFn5szy/fwFZ2bO+Xkc5jvPc57nOUxVVQIAAHBWXOsCAAAAtIQgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp4YgBAAAp9ZoEBqNxnnz5k2aNGnNmjU3ef3x48crKyubfjxFUZpRnUPDqainKAqW+jNRVRWnwgSnop6qqvi4qGeJU9FoED7//PNVVVVz5sx5/fXXf/7558aeNmXKlJ07dzb9eFVVVXhzm9TU1ODNbVJbWyvLstZV2IS6ujpJkrSuwiZIklRXV6d1FTZBluXa2lqtq7AJiqLU1NSYfbdiYw988803J0+edHV1feKJJxYvXty/f3+zHxsAAEBzN24RlpaWGgwGV1dXIurUqdOpU6esWhQAAIC13DgIBUGo77WTZVkQBCuWBAAAYD2/dY2ePHkyJSXFxcUlOTnZ29tbUZSrV696eHgcP368c+fOGpYIAABgOddahHv37r333nvz8/NzcnL69OlTUlIyefLk2bNn5+bmfvDBB48//rimRQIAAFjKtRbh/Pnz//KXv8ybN4+IRowY8fnnn7/22mv/+c9/Vq9e/fHHH/fs2bOx1yuKUlxcfObMmYYbPT09W7VqVf/rxvPqshPKJ5ECa/CcY2Xqs9vl1KFiw40AAABWxkyTGVq3bp2amhoREUFEn3zyyfLly7Oyspry+t69e587d85gMDTcOHny5Oeff77hln/sF85Xs/fDpOqqSjc3t/wK/uhW3cL7jfe0dN6pFFVVVS4uLpxjTQOqrq7W6XSi2OgYZudRW1vLOdfpdFoXoj2j0agoynWfLc5JkiSj0WgavejkTNMn3Nzcmv4SNze3W37MikRUU1Nz+fJlPz8/0yY/P79z58418RgtW7ZcsGBBVFTUzZ/2Zji9ukt+Lk/3Tm86L7s9vl1ZMlgI8dI38SgOiXOOIDQRBAFBaKLT6RCEJgjCegjCeoqiiKLYrCBsCk5EjDHGWP08d1VVGTN/h+XrfYV2rjR5qzguW1kcLYR4oU8UAAC0x4nIYDB4eXkVFRWZNhUVFdW3Ds3rsWC+9ZJw8qpaiwVVAADANlzrjBo6dGhaWprpGmFaWlp8fLzZj3SsTB2XLa+MrnvrsCEyVXq+J3+1jyCgWQgA9m/9+vUnTpyw0M4VRZFl2ck7zH18fMaOHWuhnYvjxo177733Zs+eHR0dXVlZefHixYMHDy5atMi8hzGl4H+jhY46dfkQ4ZntysJflPWF6n+jhS6eCEMAsG+ffvppeXk5plxbSElJyfHjxy0YhMnJyS1atGjXrl1eXl5KSkrPnj0//PBDLy8vMx7jWJn6ULa8OFro4cUqKoiIFoQLHjracE4NT5H+0U+YfjcGjACAfZs+fbrlPqmd3J49eyZPnmy5/YvJycmmnzp06PDUU09Z4hjuOlocI9zT6nctv3/0E9YXqn5u9MgGeX2h+kmk0BqjwwAAwOqs0RTzd2PXpaBJXAAL8WLbR4ldPCn0B2ltgfPOKQQAAK1o3yfpItCbYcKiKGHGZnnGZrkK92IDAAAr0j4ITWID2L4HxQojha2U8krQNAQAACuxlSAkolZ6WhIjzO3LE9Kl+XsVBWkIAACWZ0NBaJLciW8fJa4pUOLWSgWVCEMAALAsmwtCIurQgm0YIY4M4n1XSIuPYxEaAHBYR0rV0rpbPKesjg6XolVgQbYYhETEGc0M4VnDxbf2KeOy5Mu1WhcEAGAB+VcpaZ10kywsq6PEdVL+VSvW5HxsNAhNenqzn0aKfm7UZ4W08Ty+EAGAoxnenr3UWxieLt3w635ZHSWtk2b15CPa3+YKXO+9997FixfvqMQ7duTIkW+//faGD504ceKrr76ybjk3YNNBSESuIi0IFz4bKEzcKM/cJtfKWhcEAGBWwwLZq32ExIzrs9CUgv/Tk4/qcPsf1AsXLrx06dKdlDd+/Pjdu3ffyR5OnDixevXqGz508uTJb7755k52bha2HoQm8QFs92jxTAX1WyntvYymIQA4lD9moVlS8DpVVVXl5eXXbSwvL6+/B98NHTp06OrV5vXMqqpaUlIiSddmhY8YMeLrr7+uf/Tq1au1tbZ1ucs+gpCI2rjQijjh+VAeuwaTKwDA0TTMQrOnoKqqL7zwQufOnXv16jVw4EDTTfd++umn7t279+vXz8/P79133yWi5cuXP/TQQ/Wv8vHxeemll44fPz59+vR+/fotW7asKcc6cOBAjx49IiMj/fz8Pv74YyJasWLFAw88sG/fvl69ej355JP33nuvr6/vvHnzzPJPMws7uyf4o8F8UDs2cZOceU75apAQ4I47VwCAjdpwXq1u5lJZ8YF8QKpEROM6Mx1na5qz8KSrSDF+N/5IXLZs2YYNG/Lz811dXWfMmPHCCy8sXLjwoYceeu211x577LFTp0717ds3LCysrq6uwnRjBCIiunLlynPPPZeamvr+++9HRUU1pQZVVZOTk1999dUJEyYcP368b9++MTExRqOxqqpKkqR9+/bNnj37gw8+OH78eGho6PTp05v+r7MoOwtCIurowTYMF/+1X7l3hfTvcOHhLnbTqAUAp/L1MaW4unmdV5JC5ypVlWhrMe242LwxEb4uLMZPuOFDK1eunDJlipubGxE9/fTT0dHRBw4cqK6ufvTRR4moY8eOY8eOXb16dUhISBOPVVVVdfToUSLS6XQ9evSo337kyJHLly+PHz+eiLp27Tps2LD09HR/f3/To97e3g8//LDpoYCAgPz8/Gb9Ay3H/oKQiEROs3vxwf7skY3yylPqZwOFVnqtawIA+L3/DLpxLDXG1CP6dbRgENjru+W0oaK3me7JU1xcXH9zPS8vr9LS0pKSEi8vL8autSBbt25dUlLS8CWqqt7k2uHp06dNdyvy8fFZsWJF/fYLFy60atWqfrfe3t5XrlypD0JPT8/6ZxoMhrq6W82gtBY7bk6FtWF5o0U/N7p3hZRThGuGAGDHGl4XbGwc6W1TVXX//v2mn/ft29elS5fg4OAzZ86UlpaaNubl5QUHB3t5edXH4bFjx0xBqNPpZPn6tmn37t1zc3Nzc3MbpiARKYpSUFBQVlZm+nXv3r1du3Y1z7/BkuyyRVjPNLlieHt1QrY8thP7532C3o6THQCc1B9HxwwLZERCYoZkrnbhkiVLgoODfX19n3322WeeeaZTp04jRoyYOHHiM888s2XLll27dv33v/8losOHD7///vudOnX67LPPTA27nj17fvzxx0VFRWFhYcHBwbc8kF6vnzp16hNPPJGVlVVYWDh27NjU1FQz/AMsyRFyY2gg2zNaPFVBA1KkX8rQNAQAe9LYGFHztgtfe+21Q4cOLVmy5JVXXvnLX/5CREuXLh0yZMiXX35ZWVm5c+fO1q1bt27dOj09fc+ePatWrZo/f/7s2bNdXV0XLFgQExOzf//+4uLiphyoY8eOU6dOXbhw4dWrV7ds2eLq6tq9e/fx48e3bdv2z3/+c/3TpkyZEhgYaIZ/mDnYd4uwnq8rrYoTvj6mRKZKL/cWng7hGE4KALbv5jMlzNgubNmy5fz58xtu0el0zzzzzHVPCw8PDw8PN/1cP8Phr3/9a7OONXTo0KFDh9b/2rNnz549exLR7Nmz6zc+++yzRHTs2LFm7dlCHKFFWO/RYP7TKHH5SWXYWulcFZqGAGDrvj+pPB96s/mCwwLZnN7Cj6dw+wELcpAWYb3OHmzjCPFf+5WwlfKHA/gD5luUAQDA7KbcdevPqMQgRnSbnVzvvPNOYGBgTk7O7b28uWJiYvLy8pr+/NDQ0Jdfftly9TSRowUh/Tq5ItqPPbJR/vGk+tEAoYVO65oAALQQFxendQk34+vr6+vrq3UVjtU12tB9vmz3aNFVpNAfpc2YXAEAAI1w2CAkIg8dfRopvH0fH5slvbhDNqKPHQAA/sCRg9BkTEe+d4zu4BU1MlU6iskVAADwe44fhETU1pVS4sUpd/EBqdKCA7hxBQAA/MYpgpCIGNH0u3lOovj1cWV4unS+SuuCAADANjjgqNGb6N6K/TRS/Eee3Hel8dNIISnIWb4HAIClnT179uDBg1pX4ZgsfZ8K5wpCIhI5ze0jDA3kEzfKy/PVjyMFd6c7BwBgZhMnTpwzZ87ChQsttH9VVetv6eCcunfvbrmdO2kIhPuy3aPF53+SQ3+Q/hstRLR16ncYANyhUaNGjRo1ykI7lyTJaDS6urpaaP/gvH2Dnjr6NFKY35+PyZTm7pZlDKEBAHBKzhuEJmM78bwxup8vqpGp0vFyhCEAgNNx9iAkonautHqoOKkbj0yVPjuCWfcAAM4FQUj06+SK7BHix4eVBzPlSzVaFwQAANaCIPzNPa3Y9pFicEsK/dG4ugDdpAAATgFB+DsGgd4MExZHi3/dIs/YLFdJWhcEAAAWhiC8gcH+bN8YsUqifiul3ZfQNAQAcGQIwhtrqaf/Rguv9eUJGZhcAQDgyBCEN5Pcie8YJW46rw5MlU5gcgUAgCNCEN5CUAuWPUJ8qDMPT8HkCgAAB4QgvDVGNDOEZ48QPzqkJGfJl2u1LggAAMwHQdhUIV5s+yixiyf1/EFai8kVAACOAkHYDC4CvRkmLIoSpm/G5AoAAAeBIGy22AC2/0Gxwkj9V0l5JWgaAgDYNwTh7WilpyUxwt/78IR0af5eRUEaAgDYLQTh7UvuxLePEtcUKHFrpbOVCEMAALuEILwjHVqwDSPEkUG8zwpp8XFMrgAAsD8IwjvFGc0M4WuGif/Yo4zLkq9gcgUAgF1BEJpHPx+2Z7To50b3rpA2nkc3KQCA3UAQmo2rSAvChc8GChM3yjO3ybWy1gUBAEATIAjNLD6A7R4tnqmgfiulfZfRNAQAsHUIQvNr40Ir4oTnQ/mQNZhcAQBg6xCElvJoMP95lJhWoAxNlwoxuQIAwFYhCC2okwfbMFyM9ef3rpCWncDkCgAAW4QgtCyR0+xePG2oOHe3Mi5LPlZOO290y/u1BSrajAAAmkAQWkP/NixvtOjnRkPWSI9tkrdf+F3q/XOfsuiYIqPFCACgBQShlVybXBEpXKmlBzOl3KJrWfjPfcquS+riaEHEfwUAgBbw6WtVwwJZ3mgxxJsnpEs/nOFv71eRggAA2hK1LsDp+LpS+jDhzTyatIVCvNRdo0WBaV0TAIATQ0tEA4xI5CzSl46UqU9txQo0AABaukWLsKysrEWLFoIgWKcaJ2G6LpgSY8yvdYlMlQsq6cdYQYfvJAAAWmj00/f777/v1q1bQEDA3r17rVmQw2s4OibUi/Y9qNtSrESmSuVGrSsDAHBKjQZhQkLCkSNH7rvvPmtW4/D2lKj7L6tLYn4bHdPZg/JGi8XV6oAUqQAL0AAAWF2jQeju7s45euvM7N7W7Oto4brRMUEt2Knxuql38YgUOa8EWQgAYFUiEW3evHn79u0Nt4aHhw8YMECjkpzUzBAe6E7xa6VFUWJCewwkBQCwEpGIPD09AwICGm719PTUqB6n9mAnHuDOxmRKc/sI0+9GcxwAwBrEgwcPhoaGhoaGal0JEBHd78tyE8XhGXL+VXVeGGYYAgBYHB83bpx6owWft2/fHhcXd/z48aeeemrWrFnWr8xpdfFkW0eKW4rVSZvkOixACgBgYczb2zs1NTUiIuL2Xt+/f39vb+927do13BgfHz9mzJgbPr+iosLd3Z0xNHWoqqrKxcWlsRFJNTJN3cpLamnZIKWV3sqlWVt1dbVOpxNFrHNEtbW1nHOdTqd1IdozGo2KohgMBq0L0Z4kSUaj0dXVVetCtKcoSk1NjZubW9NfotfrbznwU7znnnt++eWX2w5CQRB69uzZrVu3hhs7d+7c2F+yTqfT6XQIQvr1VDT2P6TT0TeD6aVd6pB1LC2OBbVw5DMmSRKC0ERRFARhPUVRcCqIyPSBiVNBRIqiyLLcrFPRlLgRW7ZsWVJScttl6fX6xMTEqKioJj5fEARBEBCE9OupuPlXlfn9yd9NGbhGWRXH+/o47EkTfqV1IdozvSVwKohIURTGGE4FEamqqigKTgURmd4SZj8VvKyszMfHx7w7BTOaGcI/iOAJ6dLqAkwxBAAwP3748OHrOjbB1jzQgacNFafnyh8fxuAZAAAz476+vuHh4VqXAbfQvw3bnCS8d1CZuU1W0DIEADAf/t133+GKnV3o5MG2JIl5JepD2XIN7t0EAGAmPCQkROsaoKm8DbQuQdRzGrJGulSjdTUAAA4B63jZGYNAi2OEuAAWniIdK0MnKQDAnUIQ2h9GNLePMLsXH5QmbS1GFgIA3BEEob2aehdfFC0+mCmlnMZQUgCA24cgtGPxASwjQXxqm/LBIWQhAMBtQhDat1BvtjVJ+PwIplUAANwmBKHdC3BnGxPFfZfVsVlytaR1NQAA9gZB6Aha6SkjQXQXafAa6UK11tUAANgVBKGD0HP6OloYGsgiUqVfMK0CAKDJEISOwzSt4qXePDpNyi1CFgIANAmC0NFM6sb/Gy2OzZK+OYGhpAAAt4YgdECxASwzQZy9Q5m7G2uSAgDcAoLQMfX0ZttGCimn1embZQktQwCAxiEIHZa/G8tJFM9WqmOz5CpMqwAAaASC0JG10FFKnNjOlaJXS8WYVgEAcCMIQgcncvokUvhTFx6RIh0uxVBSAIDrIQidwswQ/vc+PGa1tOk8shAA4HcQhM7i0WC+NEZ8KFtachyDZwAAfoMgdCKD/VnWcPHlXZhWAQDwGwShc+nhxbYmiWln1Ck5shEtQwAABKET8nOjTYnihRp1RIZUbtS6GgAArSEInZG7SCvjxOCWbGCqdLYSw2cAwKkhCJ2UwOjDCGFyNx6RIueVIAsBwHkhCJ3azBD+7v08fq2UfhZZCABOCkHo7B7sxFfFi5M2SQuPYPAMADgjBCFQuC/LTRLf3q/M3S2jYQgAzgZBCEREXT3Z1pFiZqE6GdMqAMDJIAjhmtYGyhwuVkmUkC6V1WldDQCAtSAI4TcuAi2NEe7xYpGp0pkK9JICgFNAEMLvCIzeCxem3sUHpMq7LyELAcDxIQjhBmaG8PfC+bB0aXUBshAAHByCEG5sdEeeGi9Oy5U+OYzBMwDgyBCE0Kj7fNnmJPHfB5SZ22QFLUMAcFAIQriZzh5s60hxT4k6Pluuwb2bAMARIQjhFrwNtD5BFDkNWSNdqtG6GgAAc0MQwq0ZBFoSI8QFsKg06dRVdJICgENBEEKTMKK5fYRne/KBafJOTKsAAAeCIIRmmHoX/3ygkJQhpZ7BUFIAcBAIQmieoYEsJV7882blw0PIQgBwBAhCaLawNmxLkvDhIUyrAABHgCCE29HRg21JEvddVpOz5GpJ62oAAO4AghBuk5eBMhJEN5EGr5EuYloFANgtBCHcPj2nr6OFoYEsPEU6WoZOUgCwSwhCuCOmaRVzevOoNGlzEbIQAOwPghDMYHI3/nW0+GCW9G0+hpICgJ1BEIJ5xAWwzATxhZ+VubuxJikA2BMEIZhNT2+2baSw6rQ6c5sso5cUAOwEghDMyd+N5SaKx8rVsZlyFaZVAIA9QBCCmbXQUUqc6OtKMaul4mqtqwEAuBUEIZifyOnTSGFMRx6RIh0pRScpANg0BCFYyuxe/NU+PHq1lINpFQBgwxCEYEGPBfOlMWJylrT0BKZVAICNQhCCZQ32Z9nDxTk7MK0CAGwUghAsrocX2zZSTDujTs2VJbQMAcDGIAjBGvzcaOMIsahKHZEhlRu1rgYAoAEEIVhJCx2tihe7eLKBqdLZSgyfAQBbgSAE6xEYfTRAmNyNR6TIey8jCwHAJiAIwdpmhvB37udxa6WHN1y/EFuNTNNycadfALAqBCFoYGwnvjJWTD2jDFkt1WdhrUzjsuR+bZirqGVtAOBsEISgjYi2bM9oMe+yGrZSUolqZUrOkkcEsRl34z0JAFaFDx3QTFdPdmisrrBSvTdFHLeBkIIAoIlGO6F+/vnnb7/9tri4OCIiYvr06aKI7iowv3audGCsrtMy48kK1UOvugrKsEDu66p1WQDgTBqNt7y8vNjY2KCgoFdffbWsrOxvf/ubNcsCJ1Er06RN8v/1VX4p57+Uq5mFNHOb0d+NJXVgsf482o+JaCICgIU1GoTTp083/fDYY48tW7bMWvWAE6m/LvhoB0WnE17Zw67U0sWJur0lauoZ5cUd8ukKNcaPxwawkR14OzQTAcAybvF9W1GUjz76aNy4cdapBpzHH0fHzAsTvAz05Ba5jw+b20fY+YB44EFdYhDLLFS7Lzf2Wym9uEPeXKRi+iEAmJdIREuXLn3nnXcabk1OTp49e7aqqk899VRoaOjo0aM1Kg8c1sUa9eEubHyX330VmxcmfHhIqZXJRSAiautKjwbzR4NJUoTtF9S0AuWZ7XJBpTo0gCd1YHEBvJVem+IBwJEwRVEYYzd87Nlnn5Uk6b333mvsCUQUFRX1+uuvR0VFNfF4FRUV7u7uN9mh86iqqnJxceEcF8Gourpap9M1cUBW/lU1s1DNLFTXFSohXiwpiMcGsD4+DvKWqq2t5ZzrdDqtC9Ge0WhUFMVgMGhdiPYkSTIaja6uuDxAiqLU1NS4ubmZd7d8xowZN3zgf//3fzMyMgYOHPj9999v2rTJvEcFuG2dPdj0u/l3Q4TzD+vm9hGu1Kl/2iB3/laasVleflK5ihW9AaCZhMOHD0+YMMHLy+u6By5cuNChQ4fS0tIrV65wznv27HnD1//nP/8JDQ1t2bLl1QaIqLHvcXV1dXq9Hi1CIjIajaIo4lQQkSRJgiA0t3Gs49TZk8UG8Cd78JEd2JVaWpav/M92Ofu8er6KvF2ojYv9nVtZlhljgiBoXYj2FEVRVRUTt4hIURRFUdBPQESqqkqSZPZTwaKioiZMmNBYu/CWevfuffbsWb3+d9dqJk+ePHv27Bs+v7Ky0s3NDZ/+RFRdXW0wGNA1Ss3sGr25KpltKmbp54T15wWRqTHtlGH+ckxbxcVOkgVdo/XQNVoPXaP1FEWpra1t1qlwc3O75TdL0d/f/9y5c7ddVsuWLRcsWND0a4SMMVwjNBEEAdcITURRNFcQehCNa0Xj7iIiOnhFTTujLsxXpm1X+/uyxPb8gY6sQwubfu/p9XoEoQmCsB6CsJ6pZWz2a4SiqqqIJXBIPbxYDy82uxe/XEtZ55TMQvWtFMVFoMQglhTEB/kxPb6EAACRWFxc3PT2HIA98jZQciee3IkUlfaUqJmF6tzd8r7LapQfSwriw9uzQHd8FwRwXuL27dsXLlyodRkA1sAZ9fVhfX3Y7F78Ug1tOK+knlZf3CH7uV5b1C3Kj+nQTARwMuL48eO7dOmidRkA1ubjcq2ZKKtC3q+Lup26qg7257EBLCmI+5n5MgQA2Cjxiy++0LoGAC0J15qJwtw+dKGa0s8qaWfUF3cYO3uw2ACW2J5HtGUcXacAjgvz2AB+4/vrom6yKmwrvrao25kKNdqPm4bYeGEMI4DDwWRVgBsQGEW2Y5HtBAqjU1fVdYVq2hl15jZjF09mSkSHWdQNABCEALfQ0YNNv5tNv5tqZGFzkZp5Tpm4Ub5cq8YH8KQOLD6At8Ta3wD2DEEI0FQuAsUGsNgA4c2wa2t/f31MmbRJ7t362trffX3QSgSwPwhCgNvR+VozkVdJtLVYTT2jjMlUBEZxASw2gA0L5B5YHAbATiAIAe6Im3itmbggnPKvqqmn1c+OKFNz5P6+LNafj+rA7m6FZiKATUMQAphNZw82M4TNDOGVEmWfU9LOqLFrFQOn2ACWGMTiA7ihkbV/JYXEP0zkv+FGADA7/J0BmJ+7SElB/NNI4ewEMSVe6OzB3juo+C42xq2VFhxQCirVhk8+X0XRq6XLtXTdxoFpUmmdVcsGcE5oEQJYVv3a3yW1lH1OySxU56+SXQVWv/a3nxu92kcYkSGlDRVbEBFRcTWNWi+91V9ohfGoAJaHIASwkta/rv39sSr8ce3vp3rwxAzpxxgmqzRmg/RWfyHaDxcXAawBQQhgbQ3X/i6qpjUFytoCNeuc0sbA+q5iLXT0RZQwsB1SEMBKEIQAWmrnSpO78cndSFKEVafVGblSQRWNXC/F+PG4ABYfwLp4IhEBLAtBCGATSmpp/j55WbRqVNmre1hiEPvpgjovT9Fx0/QMNsSfe2OlUwALQBACaK+4mpLWSW/1F8K9Zc6ZKAp/3yWnDRU/jby2hM3yfHXG5ms3xIj15wPbscZmYgBAcyEIATR2sYYSM6R/hwsD2rLaWiKi+ACmqkJShpSeIP66hA1JirD3sppZqM7fJ4/NUsPasFh/HhvAsPw3wB1CEAJozEtPX0UJPbx+F2dDA5mfm9CiwR+oyH8bYlM/E2NsllItqYPacdOcfX83ZCJAsyEIATQmcrouBU1CvRtNtfqZGETX+k4zC9XZP8vehmvTEyPbMRf0nQI0DYIQwL7V953KqpBXgr5TgGZDEAI4CKHB9MQKI22/oGaeU2ZsVgoq1ah2PDaAJbRn7d2RiQDXQxACOKAWumv3xKCw3/pO/7ZDbqVnpskYQwO5J24UBUBECEIAh/fHvtPPjiiTNsl3t7p2QfHe1oyjoQhODEEI4Cwa9p2a7ieceU6ZsVk+XaHG+PHYADY0kHVogUgEp4MgBHBG9fcTpjAqqqbcIiWzUH1tt+IiXFvIJj6At8S9L8A5IAgBnF07198mYxy8oqad+a3v1LSQTZQf0+HWpeC4EIQA8Jv6uydWS7SlWM08p7y4Qz5apt7nyxLb85EdWCcP9J2Co0EQAsANuDboO71QTZuKlMxC9e39iv7XRcBj/bkXFqEf+lIAABlgSURBVAEHh4AgBIBb8HW9fiGb5fnq9FxjF89riTjIj+nRdwp2C0EIAM1w3SLgqWd+6zs1LWTT1wd9p2BnEIQAcDt+XQRcmNuHLtXQhvNKZqH6YKZSK6sD2/HYAJYUxP3ctK4SoAkQhABwp3xcbrAI+Is7cANFsA8IQgAwJywCDnYHF7gBwCJMC9nM7sXXJ4hnJ+hmhwrnqtTkLLndEuO4LPmzI8rZSvW6l5TW0aj1coXxdxtrZRqbJZ+vsl7l4GzQIgQAi6tfBHxBeMO+U9nr10XAhwVyDx210tPErixpnZQSL7oQEVGdQuOy5WGBDJcbwXIQhABgVTdcBHxarlzfd/rX7nzkOumHGNIxeiRLHhbI/tIdfVdgQQhCANBGw0XAy4204ZyyvlB9eINSWqd282Thq6m9G3+wM1IQLA5BCADa89TRqA58VAciotMV6toCdc4u+dRVeiccA2vA4vBVCwBsi58bW3tWfb03Temmhq2Uss9dP6YGwLzQIgQAG1KnUHKWPCyQTe1KiqJ6GoRh6dK3g8XRHdE0BEtBixAAbEV9CtZfF3wzTHitDx+fLf3nqKJtbeDA0CIEAFuhqDTjbj68/e8af3/rLbRxZa/uUkpqaVZPfHcH88O7CgBshYtA16WgydS7+OYk4YtflJnbZFwwBLNDEAKAHejQgm1JEnddUh/fJBvRSwpmhSAEAPvgbaB1CeKlGvXBTLlK0roacCAIQgCwG24irYoTfV1p8BrpUo3W1YCjQBACgD0ROS0cKAwLZFFpUsEflu0GuA0IQgCwM4xobh9h+t08IkXefxlZCHcKQQgAdmlmCH8zjMetlbYUIwvhjiAIAcBe/akrXxIjjs2UVhcgC+H2IQgBwI4N8Wcp8eK0XOnzXzCpAm4TghAA7FtYG5abKM7fq8zdLWtdC9glBCEA2L0uniwnUVx1Wp25TVbQSwrNhCAEAEfg50YbRoh5JerEjVh6BpoHQQgADqKVntYliHUKJaRL5UatqwH7gSAEAMdhEOibwUJXTzZktXShWutqwE4gCAHAoQiMPokUxnbiEanS8XJcMIRbQxACgAOa3Ys/G8Kj0uS8EmQh3AKCEAAc0xP38A8H8Pi10vpCZCHcDIIQABzWAx34d0PERzZKy09iICk0CkEIAI4s2o9lDRdn/aR8chhZCDeGIAQABxfixTaOEN49oLy4A0vPwA3cLAjPnTt3+vRpScKtoAHAvnXyYDmJ4vpCdUqOLKFlCL8nNvbArFmzjhw5otfrjxw5smLFirvuusuaZQEAmFdbV9o0QhybJSVnyUtjBNdGP/zA6TTaInz77bfT0tJ+/PHHKVOmLFq0yJo1AQBYQgsdpcSLLiINz5DK6rSuBmzGzbpGL1y4sGvXruzs7MjISKsVBABgOXpOS2OEvj4sMlUqrMS0CiAydY1mZGQcOnSo4da4uLiQkJBVq1ZlZmaWlpaiXxQAHAYjevs+oY2LMjBNTh8mdGvJtK4INCYSkYeHR5s2bRpudXV1JaJp06ZNmzbthx9+eOmll7755httCgQAsIDZvXg7NxqUJq2ME+/3RRY6NbGwsDAiIiIiIuK6B2praw0GAxGVl5ebfgAAcCSPBXNvA41cJ30dLQ4LRBY6L3HatGlr1qz54wOjRo0yGo2yLFdVVaE5CAAOKSmIp8az0ZnSv+4TJnTBvGonxXQ63YkTJ9q3b//Hx0pLSwVB8PDwuMnr+/fv36ZNm3bt2jXcGBsbO3r06Bs+v6Kiwt3dnTF8+aKqqioXFxfO8bdH1dXVOp1OFDGenWpraznnOp1O60K0ZzQaFUWxTnfU4TIalS389W71me62OMdQkiSj0Wi6YuXkFEWpqalxc3Nr+kv0ev0tP2ZFPz+/o0eP3jAIW7VqdctjCIJw9913BwcHN9zYsWPHxv6SdTqdTqdDENKvpwJBSESSJCEITRRFQRDWUxTFOqci1Ie2JKrD19G5avGd+2zu48n0gYl3BREpiiLLcrNORVPiRmzZsmVJScltl6XX60eOHBkVFdXE5wuCIAgCgpB+PRUIQvr1VAiCoHUh2jO9JXAqiEhRFMaY1U5Few/KSaSkddKUzfT5QEFnS3+XqqoqioJ3BRGZ3hJmPxW8rKzMx8fHvDsFALA7XgZalyCW1KhjMqUqrCzpTPj58+e7du2qdRkAANpzE2lVvNjWlQ1eI12q0boasBYeFRUVFBSkdRkAADZBYLRwoDAskA1Kkwqw9Ixz4J9//rnWNQAA2BBGNLeP8OfuPCJF3ncZWej4eIcOHbSuAQDA5jzdg78ZxuPWSluKkYUOzpaGRgEA2JI/deVLY8SxmdLqAmShI0MQAgA0aog/S4kXp+VKn/9ii3PtwSwQhAAANxPWhuUmivP3KnN3y1rXAhaBIAQAuIUuniwnUVx1Wn16m6ygl9ThIAgBAG7Nz402jBD3lqgTN8p16CV1LAhCAIAmaaWndQmiUaGEdKncqHU1YD4IQgCApjIItGywEOzJhqyWLlRrXQ2YCYIQAKAZBEafRApjO/HwFOl4OS4YOgIEIQBAs83uxf+nJ49Kk/NKkIV2D0EIAHA7/noP/3AAj18rrS9EFto3BCEAwG16oAP/boj4yEbpu3wMJLVjCEIAgNsX7ceyhovP/6x8chhZaK8QhAAAdyTEi20cIbx7QHlxB5aesUsIQgCAO9XJg+UkipmF6uQcWULL0N4gCAEAzKCtK20cIZ6rUpOz5GpJ62qgORCEAADm0UJHKfGii0gJGVJZndbVQJMhCAEAzEbPaWmM0M+HRaZKhZWYVmEfEIQAAObEiN6+T5h6Fx+YJh8tQxbaAVHrAgAAHNDMEO5loEFp0oo4MdyXaV0O3AxahAAAFvFoMF84UBi1Tko/i3ahTUMQAgBYSlIQT40XJ+dIy05gUoXtQhACAFjQfb4sc7j44g7l7f3IQhuFIAQAsKx7WrFtI4Wvjykzt8noJLVBCEIAAIvzd2ObRoi7LqmPb5KNaBnaGAQhAIA1eBloXYJYUqOOyZSqsPSMLUEQAgBYiZtIq+LFdq4sZrV0qUbrauBXCEIAAOsRGH02UEhozwalSWcqcMXQJiAIAQCsihHN7SP8uTsfkCrvu4ws1B6CEABAA0/34G+G8bi10uYiZKHGEIQAANr4U1e+NEZMzpJWFyALtYQgBADQzBB/ljpUnJYrLTyCSRWaQRACAGipnw/LTRTf2qfM3S1rXYuTQhACAGisiyfLTRJXnVaf3iYr6CW1OgQhAID22rnShhHi3hJ1XLZci5ahdSEIAQBsQis9rUsQOdHwDKncqHU1zgRBCABgKwwCLRssdGvJhqyWjpTe4AkltYSuU7NDEAIA2BCB0ccDhFEdeO8fjT+c/N1Q0rOV6tC1Un45otDMRK0LAACA6718Lxc5PbxR/lyiCZ2IiM5Wqg+slz+MELp4Mq2rczQIQgAAW/RiL97OjabkyiU1LNGfxufKH0YI9/kiBc0PQQgAYKMeD+buIj2cLb/tIv4YL/RvgxS0CFwjBACwXeG+rJMHXaihJcex9IyloEUIAGCjTNcFvxrIBJKGrKNTV+n7WEGH9ou54YwCANii+tEx/dtQaCt172hdTrESmSqV1WldmcNBEAIA2JxqiUavlz8a8NvomC6etHe0ePKqet8q3NHXzBCEAAA2x1Wk9ATxutExQS3YkWTdn7vziFR51yVkodkgCAEAbFFrww02ehvomRD+fjhPSJdSz2D4jHlgsAwAgJ0Z3ZH7u7ExmfLp3vTkPWjP3CmcQQAA+3OfL9ucJHx0SJmJOzfdMQQhAIBd6uTBtiSJ+y+rD2bKVZLW1dgzBCEAgL3yMlB6guipp+jVUnG11tXYLQQhAIAd03NaFCX8qQsPWyntvYxO0tuBIAQAsHszQ/jb9/G4NdLaAmRhsyEIAQAcwbjOfGWcODlH+uQwplU0D4IQAMBBRLRluUniuwcwlLR5EIQAAI6jqyfbOlLcU6KOz5arMZS0aRCEAAAOpbWB1ieIOk5D1kgXa7Suxh4gCAEAHI1BoMUxQnwgC0+RfilDJ+ktIAgBABwQI5rbR5jTm0elSTlFyMKbQRACADisyd344mgxOUvCDe5vAkEIAODIYgNY1nDx5V3K3N2y1rXYqFsE4fLlyxcuXGidUgAAwBJCvNjWJHH1GXVyjmxEy/APbhaEhw8fnjdv3rJly6xWDQAAWIKfG21MFEtqKSFdKq3Tuhob02gQKory/PPPv/baa9asBgAALMRdpB9jhR5eLDJVOl2B4TO/aTQI//3vf48ZM8bf39+a1QAAgOUIjBaEC0/34JGp8o6LyMJrRCJ64IEHzp4923Drv/71rw0bNqSkpOzevVujwgAAwCKm382DWrCR66RPIoVRHTBkkkQiWrly5XVbV69eXVZWFh8fX15efuzYsVmzZr399ttalAcAAOY3LJClJ4hJGfKeEnVuH0HrcjTGf/755z9uHTFiRE5Ozvr16z/66KPevXsjBQEAHEwvb7ZtpLDqtDp9syw591BS/sYbb9zk4YCAgGnTpt3kCaqqlpSUnP+9q1evmrtOAAAwswB3lpMoFlaqieukcqPW1WiHGQyG6upqxtjtvb53795nzpzR6/UNN06ZMuXFF1+84fMrKyvd3Nxu+3COpLq62mAwcI4OeqqurtbpdKIoal2I9mpraznnOp1O60K0ZzQaFUUxGAxaF6I9SZKMRqOrq6ul9q/QC3t02y/x7wbWBrpZ6CDmoShKbW1ts06Fm5ubINyi71esra0tKSnx8fG5vbJatmy5YsWKqKioJj6fMebu7o4gJCJBEFxcXBCERCSKIoLQRK/XIwhNEIT1LB2ERLQwmhYcUOKzeUq8cG9r2/18VhRFp9O5uZk5rjkRqSoG0QIAOLWZIfzf9/Oha6XVBU6XCFyv1992cxAAABzGg514Srw4LVf66JBzDZ7hQ4YMQUclAAAQ0f2+bHOS+P4hZeY2WXGaliF/6aWXtK4BAABsRWcPtiVJ3HtZTc6SqyStq7EKPmDAAK1rAAAAG+JtoHUJoptIg9dIF6q1rsbyMGQRAACup+f0dbQwLJCFp0iHSx28kxRBCAAAN8CI5vYRXrmXx6yWNp535CxEEAIAQKMe78aXxogPZUtfH3PYoaQIQgAAuJnB/mxzovhGnvLiDtkhG4YIQgAAuIXglmzrSHFzkTohW66Rta7G3BCEAABwa60NlDlcZIxi10iXarSuxqwQhAAA0CQuAi2NEWIDWHiKdLTMcXpJEYQAANBUpqGks3vxqDQpt8hBshBBCAAAzTP1Lr4oWhybJS074QhDSRGEAADQbPEBLDNB/NsOZe5uux88gyAEAIDb0dObbR0ppJ5Rp+bKRntuGSIIAQDgNvm7sU0jxOJqdUSGVFandTW3C0EIAAC3r4WOVsaJ3VqyyFTpTIVdDp9BEAIAwB0RGH0QIUy9i0ekyrsu2V8WIggBAMAMZobw98N5QrqUesbOLhiKWhcAAAAOYnRHHuDORq+XT/Wip3rYTUPLbgoFAADb178N25wkfHJYmblNVuyklxRBCAAA5tTJg20ZKe6/rD6YKVdJWlfTBAhCAAAws1Z6WpcgtnWliBTpbKWtNwwRhAAAYH4ip08ihUndeESKnFdi01mIIAQAAEuZGcL/dT+PXyutLbDdLEQQAgCABSV34qvixck50ieHbXRaBYIQAAAsK9yXbU4S/33ARoeSIggBAMDiuniyrSPFvBL1oWy52saGkiIIAQDAGrwNtC5B1HMavEa6UK11NQ0gCAEAwEoMAi2OEYYGsohU6ZcyW+kkRRACAID1MKK5fYQ5vXlUmrTpvE1kIYIQAACsbXI3vjhaHJctLT6u/VBSBCEAAGggNoBlDxdf2aXM3S1rWwmCEAAAtNHDi20bKa4+o07KkY3atQwRhAAAoJl2rrQxUbxSSwnpUmmdNjUgCAEAQEvuIv0QK/RrwwakSKeuajB8BkEIAAAaExi9GSbMDOHhKdKWYmtnIYIQAABswvS7+VdR4uj10rf5Vr1giCAEAABbMTSQZQ4XX/jZqkNJEYQAAGBDQr3Z1iQh5bQ6fbMsWaVliCAEAADbEuDOchLFc5XqiAyp3GjxwyEIAQDA5rTQ0ap4sWtLFpkqpZ25QcPwlzLVXBmJIAQAAFskMPowQpgUzJOzlJd3/u6S4cEr6oRs+WK1ecaXimbZCwAAgCU825P7u9NjG+XLNfRBBCOiI6XqIxvlpTFCF09mlkMgCAEAwKY91Jn7u7P4NVJJLftbDzZpm7w0RujeyjwpSAhCAACwfQPbst0PiP1WyVuLhewRQnBLs6Ug4RohAADYBYWoiwd18aAN5r6LIVqEAABg60zXBZfFsAC99Pg2QVKUv95jtoYcWoQAAGDTjpSqEzZcuy6o5/TdEGFdofrRIbNNtkcQAgCA7ZIU+vMW+bvBv42O0XP6drCQflbdU4LpEwAA4OhETlnDReH3g2MMAq2IEwQzjZixdouwoKCgoqLCyge1TYWFheXl5VpXYROKioquXLmidRU24eLFi5cuXdK6Cptw6dKlixcval2FTbhy5UpRUZHWVWipPvDKy8sLCwuv23jnrB2ETz/99ObNm618UNs0Z86ctWvXal2FTXjjjTeWL1+udRU2YcGCBV9++aXWVdiEL7/8csGCBVpXYROWL1/+xhtvaF2FTVi7du2cOXPMvltcIwQAAKeGIAQAAKd2p4NlJEnau3dv059/5cqV/fv3u7m53eFxHcDFixcPHTq0adMmrQvRXlFR0bFjx3AqiKigoKC0tBSngohOnTp19epVnAoiOnbsWFFREU4FER06dOjixYvNOhX9+vVzd3e/+XOYqt7R8NNly5Z98MEHotjUQK2oqHBxcWn68x1YZWWlXq/X6XRaF6K96upqQRD0er3WhWivpqaGMWYwGLQuRHu1tbWqqrq4uGhdiPbq6upkWXZ1ddW6EO0Zjca6urpbBltDX3zxRdeuXW/+nDsNQgAAALuGa4QAAODUEIQAAODUEIQAAODUEIQAAODUrDp688qVK/n5+YGBgW3btrXmcW2Noig7d+48efKkn59fZGQk5877daS8vHznzp0XLlxo3bp1VFQUBo4SUX5+fmVlZc+ePbUuRDOXLl06ffp0/a/du3d38glXhw8f3rt3b6tWrQYMGODh4aF1Odo4depUSUlJ/a+CIPTu3dtcO7deEI4cOTIjI0MUxX/84x/PPvus1Y5rgyIiIqqqqkJCQvbt2+fi4rJhwwanfXM/8cQT58+f9/f3/+WXX0pKSnJzc/38/LQuSksFBQVhYWEtWrRomATO5scff3zllVdCQ0NNv3788ce3HP7uwGbNmrVkyZJBgwaVlZXt3Lnz5Zdf1roibXz//fcZGRmmn0+cOOHq6nrw4EGz7V21lqNHj9bW1iYkJLzzzjtWO6htOnTokOkHo9HYo0eP999/X9t6bIGiKAMHDnzzzTe1LkRjCQkJzz33XFBQkNaFaOnTTz8dO3as1lXYhBUrVgQGBl64cEHrQmzLoEGD5s+fb8YdWq9TLjg4GB1fJt27dzf9IIqiv79/TU2NtvXYCEmSvLy8tK5CS1999ZWvr298fLzWhWivvLw8Ozv7wIEDimK2m6/aoyVLlkyfPr2uri43N7e0tFTrcmxCfn7+tm3bJk6caMZ9Ou/VKVuwbdu2HTt2JCcna12IltLS0saMGXPPPff07t170qRJWpejmaKiov/7v//75z//qXUhNqGgoGD+/PnDhg0bMGBAwytDzubEiRM5OTmJiYlvvfVWcHDw+vXrta5Ie59//vnw4cPNew0FQaiZ48ePjxs37qOPPurQoYPWtWgpNDR02rRpjz/++I8//vjTTz9pXY5mnnjiiblz57Zp00brQrQ3adKkQ4cOZWRknDhxokWLFn//+9+1rkgzVVVVVVVVO3fuTE1NfeONN5544gmtK9KYJEmLFi2aPHmyeXeLNT+1cfr06djY2FdffXXChAla16KxoKCgoKCghISEioqKd999NzIyUuuKNHDkyJH09HQfH59NmzadPXv28uXLM2bMmDdvnre3t9alaaB+AV6DwfDQQw998cUX2tajIT8/v/DwcEEQiGjIkCEzZsyorq525kVH09PTZVlOSEgw724RhBo4e/bskCFDZs2aNW3aNK1rsSGlpaXNWkvXkfj6+r777rumn11dXbdv3963b19cUyeivXv3BgQEaF2FZqKjo48cOWL6OT8/39vb25lTkIi+/PLLxx57zOz3KrDeotvLly/ftWvX8uXLO3bsGBYWNn78eDPOArEvffr0KS8vHzt2rOnX+++//4EHHtC2JK3Ex8dHRET4+Pjk5eV99913GzZs6Nu3r9ZFaSwjI2P69OnOPH1iypQpbdu29fPz27Vr1w8//LBx40anfVcUFxf36dPnkUce6dSp01tvvfXkk08+99xzWhelmQsXLrRv3z4vL69+vKG5WC8IN27cePTo0fpfBw8e7LRzgxYvXlxVVVX/6z333OOc/YFEtH79+m3btpWWlgYFBY0bN87f31/rirR39uzZnJychx9+WOtCNLN58+bs7OzS0tL27duPGzfOmVuERFRYWLho0aLy8vK4uLghQ4ZoXY6W8vPzt2/fbok/DdyGCQAAnBpGjQIAgFNDEAIAgFNDEAIAgFP7f0hW59ajczr/AAAAAElFTkSuQmCC", "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.4" }, "kernelspec": { "name": "julia-1.5", "display_name": "Julia 1.5.4", "language": "julia" } }, "nbformat": 4 }