{ "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 a bulk silicon lattice,\n", "see Input and output formats for more details." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using ASEconvert\n", "\n", "system = pyconvert(AbstractSystem, ase.build.bulk(\"Si\"))\n", "model = model_LDA(attach_psp(system; Si=\"hgh/pbe/si-q4.hgh\"))\n", "basis = PlaneWaveBasis(model; Ecut=5, kgrid=[3, 3, 3]);" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "DFTK already defines a few callback functions for standard\n", "tasks. One example is the usual convergence table,\n", "which is defined in the callback `ScfDefaultCallback`.\n", "Another example is `ScfPlotTrace`, which records the total\n", "energy at each iteration and uses it to plot the convergence\n", "of the SCF graphically once it is converged.\n", "For details and other callbacks\n", "see [`src/scf/scf_callbacks.jl`](https://dftk.org/blob/master/src/scf/scf_callbacks.jl).\n", "\n", "!!! note \"Callbacks are not exported\"\n", " Callbacks are not exported from the DFTK namespace as of now,\n", " so you will need to use them, e.g., as `DFTK.ScfDefaultCallback`\n", " and `DFTK.ScfPlotTrace`." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In this example we define a custom callback, which plots\n", "the change in density at each SCF iteration after the SCF\n", "has finished. For this we first define the empty plot canvas\n", "and an empty container for all the density differences:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots\n", "p = plot(yaxis=:log)\n", "density_differences = Float64[];" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "The callback function itself gets passed a named tuple\n", "similar to the one returned by `self_consistent_field`,\n", "which contains the input and output density of the SCF step\n", "as `ρin` and `ρout`. Since the callback gets called\n", "both during the SCF iterations as well as after convergence\n", "just before `self_consistent_field` finishes we can both\n", "collect the data and initiate the plotting in one function." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using LinearAlgebra\n", "\n", "function plot_callback(info)\n", " if info.stage == :finalize\n", " plot!(p, density_differences, label=\"|ρout - ρin|\", markershape=:x)\n", " else\n", " push!(density_differences, norm(info.ρout - info.ρin))\n", " end\n", " info\n", "end\n", "callback = DFTK.ScfDefaultCallback() ∘ plot_callback;" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Notice that for constructing the `callback` function we chained the `plot_callback`\n", "(which does the plotting) with the `ScfDefaultCallback`, such that when using\n", "the `plot_callback` function with `self_consistent_field` we still get the usual\n", "convergence table printed. We run the SCF with this callback …" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) α Diag\n", "--- --------------- --------- --------- ---- ----\n", " 1 -7.774222897082 -0.70 0.80 4.8\n", " 2 -7.779029743647 -2.32 -1.52 0.80 1.0\n", " 3 -7.779318616829 -3.54 -2.58 0.80 1.5\n", " 4 -7.779350231472 -4.50 -2.92 0.80 2.5\n", " 5 -7.779350723614 -6.31 -3.26 0.80 1.0\n", " 6 -7.779350850199 -6.90 -4.30 0.80 1.0\n", " 7 -7.779350856066 -8.23 -5.13 0.80 2.0\n" ] } ], "cell_type": "code", "source": [ "scfres = self_consistent_field(basis; tol=1e-5, 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+gvaeTAAAgAElEQVR4nO3dd0AUZ/4/8OeZGXZhAaWJoBTBhkFFUVQUUIqgAoI9MVETS0w5412Siya/+ybeJUaTfC8xF89vzqixxFw0FhQElSI2VBCNsYuKSrErSmen/P5YQwy2BZad3Z3366/dcXfmM4D73qfMPFSSJAIAAKBUjNwFAAAAyAlBCAAAioYgBAAARUMQAgCAoiEIAQBA0RCEAACgaAhCAABQNAQhAAAoGoIQAAAUDUEIAACKZkJBWFFR8eGHH+r/eq1W23LFmD6cvtwlyInneblLkBNOX8m3xhQEweCnb0JBeOvWrTVr1uj/+pqampYrxvTh9OUuQU41NTVK/iisra0VRVHuKmRTV1en5N9+XV2dwX/7JhSEAAAAxocgBAAARUMQAgCAohkvCIuLi69fv260wwEAAOiDM85h3njjjaKiovLy8ujo6A8++KDJ+ynXEiuGWLMNt9+sIW2sm1UhAAAokzFahCdOnDhy5MjWrVszMjJWrlx548aNJu8qq1Qck8HXCH/Y+N0Z8a0DwhPeAQAA8DTGCMLc3NywsDBKKcdxQUFBR48ebfKuEryZsT7MqHS++rfriJafFZMui9+HPdJIBAAA0IMxgvDu3bt2dna6x61atbpz505z9vZKF2a8LzM6g6/myfKz4qZL4sYo7tHOUgAAAH0YY4zQxcWluLhY9/j27dtt2rRp5g5f6cJIhPTfqfa1F1NjODVSEAAAmqpZLcKcnJwpU6YEBARMmDDh4e3Hjx8PCgrSaDQBAQG5ubmDBw/OyMjQarXl5eX5+fn9+vVrXs2EECKIRBDI/utSWrFy77AAAADN16wWYWVlZWBgoJOT0+HDhx/ePnHixJdeeungwYPLly8fP378hQsXJk2aFBoaqtVqP/roo1atWjWv5gfjgvkjar+5oHkxi3+5C/NVMKvCJZEAIKusrKxvv/22pY8iCALDMJTSlj6QqYmJiZk2bVpL7Jk2/55133333erVq/fu3at7mpubGx0dfePGDZVKJUmSl5fXsmXLYmJieJ5nGIZhnphXp0+fDgoKiomJqd8yYcKEuLi4Bi9bdYHZUkR/DBX46go7O7ulBczC44y7DVkTwvvaN/NUzElFRUX9yKsC4fRtbW0V+FGoU1lZaWNj85QPE7n885//PHLkyMSJE+UuxALl5OQUFxevWrWqurpapVKxrL5DYiqViuOe0eQz/BhhQUFB165dVSoVIYRS6u/vf+7cuZiYmGeWolarbW1tx48fX7+lR48earX64ddsuSJtK5U2RTFqlisX6tRq9azuxE4t/ecMiUi3WtSfTvBVykdDXV1dgx+Ooij89LVarVqtVmwQ8jyvVqtNMAhZlvXz8xs3bpzchVggURQ3bdqkVqtFUWxUEOrzd2L4IHx4jihpzDRRhmE0Gk2D4cYGoj1IrBfRzY5hWVb3s5juR0b5kMvl0vO7hO0l5NsQVmOk+wTIqf70lQmnz7KsYoNQd/omGIQmWJIloZSyDzHgng3/a3N2di4vL69/WlZW1vxpovXsrMhj54g6q0mgC81P5ERCgpL4E3cxgwYAAPRi+CDs2rXrmTNnamtrCSGiKJ44caJr164GP8pj2VuRH4awcwKY8G381yeUu1wZAADor1lBeO/evfz8/CtXrlRUVOTn51+8eJEQEhgY2KlTp/nz51dVVS1atEij0YSHhxuoWr1M7szsjeNWnBMnZwsVil7GHAAAnq1ZQfjrr7/OnDkzLS2NZdmZM2fWzxtet27d/v37PT09N2/evHnzZuP3m/s50EMJnKOa9Enif7mNblIAMGmVPKl7ch+WREhZnRGrUZ5mzSoJDQ1tcAWhTufOnTMzM5uz5+azZsnXwezGQnH4dn5uAPtWd0ah8woAwORNzhYEiayPfMz10BIhsw8I+65JhxM5fIq1EAuf4zTGhzmYwP10URyVLtyplbsaAIDH+SiQ2X9dHJXO1/5xHR2JkD8fEJaeET/uyzYtBX19fZuz4M/DRFHMyMgQhOYu9aPVal1cXKqqqp7ymt69excUFDTzQPqz8CAkhHjb0T1xXC9nEriZ338d3aQAYHJ6OtHMEVzuTWl0xu9ZqEvB/5wRN0ZxsZ5NbAzW1tY2/64pOjzPDx06tLq6upn7YVn2ww8/1F1r/iR1dXWiaLwJjwq44I4QK4bMC2T7uIhjMvjXujEf9m7idysAgBaiy8LIVH50Br8pilOxBkjBR61fv/7IkSO+vr5TpkypvyVFZmbmvn37XFxcXnjhBScnJ0JIdna2s7Nzjx49CCElJSWHDh0aPXp0UlISIeT7779Xq9Xx8fHu7u76HLGqqmrp0qU3b94cMGBAfHy8bqO19YOF1Ldu3dqjR4+9e/cWFBREREQYeWZlPctvEdaL92IOJ3JZpVJ0Gn/1aY1yAAAZPNwunJVj+BScPXv2kiVLOnXqlJycXH/3yk8//fSNN95o06bNiRMnevfuffv2bULIqlWr6ud5nDlz5tNPP23aEWtqagYPHnzq1CkvL6+5c+fOnz+fEMLz/MyZM3WX2H3++eejRo06e/asg4PDuHHjsrKyDHCejaeIFmE9D1u6K5b7+KjQJ0n7fRgX44GGIQC0oIwS6bNfGzeo1rEV3VEiSZLU05kuOiEsOtGI907uzEzq9PjmTVFR0erVqy9fvtyqVaspU6Z07Nhxz549gYGBn3zySV5enr+/PyEkMTHx3//+94cffvjYPSQmJhJCXnnlFf1v87thwwZ7e/ulS5cSQsLCwoKCgt55550G1xEkJibOmzePEHL37t2kpKSIiAh9z9ZwlBWEhBCWknmB7GB3ZlK2MKYD/d/+rJWCWsUAYFR9XOicno24GZhEyJJTIkMkK5bYsvSdHgzXmK/r3Rye+E/Hjx/38/PTLf5jZWXVt2/f48ePOzg4qNVqXQoSQkJCQnJzcxtxPEIIITdv3nz77bcJIQzDrFq16uF/+vXXXwcMGPCgtm7dOI4rLCzs2LHjw6/p1auX7oGHh8fu3bsbe3SDUFwQ6oS706OjuCm7+ahU/sdwtr0tmoYAYHiOahLVXt+PF93smLRicfNQztOWRKby35yUNkUZZu3xysrKmpqa+qfV1dUajUaj0ehm0+juW6vbSAhhWbZ+dujTp3cSQmxtbceMGUMIefTmt5WVlfWTYkRRrK2t1e3/YaZw02DltobaWJNtMdzoDkzQFj7lCmaTAoCcGswRfew80mY6fvz4iRMnCCHFxcX79+8PCQnp0KGDq6vrhg0bCCHV1dXr168fMmQIIcTb2/vIkSOEEEmSNm7cqHu7SqXSaDSPLqKg0WgSExMTExMTEhIePejWrVsrKysJIZs3b27Xrp2np6dhTsagFNoi1KGEzO7OBLelL2QJ6SX0i/5Y3RcAZPDYKyUazCNtfruwR48e06dPd3V1zcvLe//99zt37kwIWbFixaRJk1asWFFQUNCvX79JkyYRQqZOnRocHDx48OCKior6jlNCyKxZs4KCgjw9Pb/77rvevXvrc9CuXbsOGjTIw8MjNzf3xx9/NM0FOhQdhDr92tAjo7hX9wkDt/I/RbCdWqGbFACM6u2DwndnxeRobugf+1Hrs3BClrApqrnXfTk7O6empv7yyy8eHh7t2rXTbYyIiDh//vzp06fbtGlT31xr3779mTNnTp486eHh0aZNm/re0YULF3788ccVFRX29vougz5w4MA333zz3Llz3bp1041QqlSq69ev29raEkLS0tLqL6V45ZVXdDFsfAhCQghprSLrItjVBWJIMv/VAPaFjqb4nQUALFWwKx3h2TAFdXRZuPuqZJCrn1UqVb9+/RpstLGxCQwMbLBRo9EEBQXpHusCTMfKysrR0bFRB3V0dOzfv//DW1xdXXUPHg5UtVot12rbCMLfTe7MdHekz+8S0ookhazuCwCmYLzv075893SiPZ2aGIMeHh662Sh+fn7PP/9803bSZEOHDn04RPXXrl07Kysrg9fzJPiw/wPd6r6v7xf6JvHrItgeTf3jAwAwBYcOHdI96NGjh+5OMcaku/SwCdLT0w1bydOhD7Ah3eq+cwOYiFSs7gsAYPkQhI9Xv7rvuEwBK4EBAFgwBOET6Vb3bach/bfwR7G6LwCAhUIQPo1udd8FQcyI7fzXJ0SEIQCA5cFkmWcb3YHp40JfyBJ2XZVWhLFO8szvBQBTx3HcypUrW3qiR/0d0RTl7t279bctNTgEoV50q/t+8osQuJn/YQgb4qa4v0IAeKaZM2caYUW9mpoalUplmrdoaVEeHh4ttGcEob44hswLZAe1lZ7fJUzvSrG6LwA0YG9v36dPn5Y+SlVVlbW1tQKDsOXgR9k4Q9vTgyPZXaXSUKzuCwBgERCEjeZhS7NiuVA32idJu6MYE2gAAMwbgrApdKv7rg3npu0VZh8QtLjsHgDAbCEIm063um/BfSk0hS8sR9MQAMAsIQibRbe67wu+zKBkrO4LAGCWEITNpVvdd8tQ7s8HhdkHhDp0kwIAmBUEoWEEtaH5idy1ajJwK3/+PpqGAABmA0FoMLrVfd/yZ0KS+R8voGEIAGAeEIQGNrkzkzmC+/QXcXK2UMnLXQ0AADwLgtDw/B3pgZGcREhQEn/8DrpJAQBMGoKwRdhbkTVY3RcAwBwgCFuQbnXf78+JY7G6LwCAqUIQtiw/B3owgWuvIb038wdvoJsUAMDkIAhbnG513y/7M6PSsbovAIDJQRAayagOzMEEbt1FMTFduFMrdzUAAPAbBKHx6Fb37e1MAjfz+66hZQgAYBIQhEalW913WSj7/C5h3hEB/aQAALJDEMogqj3NT+RyrmN1XwAA+SEI5dHWhqQN40LdaOBm7Xas7gsAIB8EoWx0q/v+N4KbgdV9AQDkgyCU2RB3emQUd/6+FJKM1X0BAGSAIJRfG2uSEsNN7Mj038L/XIiGIQCAUSEITYJudd9tMdz7eeLsA0KtIHdBAACKgSA0IUFt6JFR3PVqMjAZq/sCABgJgtC0tLIiP0Wws/2ZkGR+7Xl0kwIAtDgEoSnSre674BhW9wUAaHEIQhPl70jzEjhrjvRN4n/F6r4AAC0GQWi6bDiyNIR9P4CJxOq+AAAtBkFo6iZ3ZvbFcysLxLGZwukyqbTqMa3DvJtoMgIANBGC0Ax0bU0PjOTaa0hkqhCVKjTIwv93WFhxDu1FAIAmQhCaB93qvv8eyFyrlvptEYorH2Th/+QLN6rJvwey8pYHAGC+EITmZFQH5pdRnJOa9NjIn7pH/ydfuFZF/hPCMlTuygAAzBaC0Mx42dH8RG64JxOyQ/3LLWlpKFIQAKBZEITmx4ohPvYkwFHaUSylFWF0EACgWRCE5kc3LpgZVbuwH5uQLqy/iCwEAGg6Tu4CoHHqxwUrK8jbPRhrjryULfASmdgR32kAAJoCQWhOki6LN6vJ0lC2fljwjW4MS8nU3YKGI4neyEIAgEZDEJqTBG8mwZs0mBwz04/p40xHpvNakYzzQRYCADQOgtCcPGl+aN82NGMEF5MmVGjJK12QhQAAjYAgtBDPOdCMEezQVKFCS2b5IwsBAPSFILQcXVvTvfFsVKqgFcnbPZCFAAB6MVIQnj179tChQzzPT5061ThHVCZvO7onjotK5W/USAuDcN81AIBnM1K7YdmyZfn5+YsWLTLO4ZTMXUOyYrm0ImluniB3LQAAZsBIQfjFF1+88847xjkWtLUh2bHcrlLpzRwB6zMBADwdRpIsk6Oa7BzO/XJbmrlPEBGGAABPZuAxwqFDh965c+fhLdu2bXNzczPsUUAfrVVk53AuMZ1/KVtYPZjl8J0HAOBxGhGEt27dunz5sq+vr6OjY/3G2trarKys2traIUOGODg4pKent0CR0ES2HEmJ5sZnCWMyhfURrBqzZwAAHqFvMyEoKMjDw2PgwIEZGRn1GysrKwcOHLhgwYI1a9Z069btwoULT3r7xYsXjx8/Xl1dnZ+ff/v27eZWDXpTs+TnSJajZHQGX83LXQ0AgOnRNwhXrVpVXl7epUuXBhutra2zs7M3btw4duzYhQsXPuntGRkZW7dujYiIWLp06blz55pVMjSSiiHrI1kXazp8B1+ulbsaAAATo2/X6HPPPffoxq1bt44bN45hGELI+PHjx4wZ89133z327a+++uozD1FVVXX58mVKf7+P2Lx5854y17SysvLhFytNY0//X4HkrTwuJlXYOJi358x+/gx++5IkKfYnUFVVJQiC7pNHgaqqqnieV+zpV1dXq1QqltV3pMfa2prjnpF0zZosU1JS4uHhoXvs4eFx69atmpoaa2vrpu1No9F4e3sXFhbq+XpJkuzs7Jp2LAvQhNNfGUHePijEZ7M7hnPO6haqy0gU/tsnhNja2io2CCmlNjY2ik0ChmGsra0Ve/osyzYqCPXRrB8lz/P11XAcJ0kSz2MYynRRQr4awI7wpGHJ/NUquasBADANzQpCNze3mzdv6h7fuHGjVatWCv+Sbhb+0Yed3JkJ38YXV5p9BykAQPM1KwjDwsLqJ5FmZGSEhYUZoiRocXMCmGldmdAU4cJ9ZCEAKJ2+Y4TLli07f/78tWvX1q5dm5+fP3PmTB8fn1dffTUgIOC9995zc3NbsGDB1q1bW7RWMKC/9mTsrUhkqpA+nO3cWqFDTQAARP8gdHd3F0Vx/vz5uqe6GTHt27fPy8v7/vvvS0tL09PTg4KCWqpMaAGvdWPsrUh4qpAWw/ZwQhYCgELpG4SxsbGP3e7j4/OPf/zDcPWAUb3YieEYEpXGJ0dz/dogCwFAiRQ6ARfqTfBllodyI3fyB25gvBAAlAhBCCTOi64czI3cyWeWIgsBQHEQhEAIIcM86MYo7sVd/I5iZCEAKAuCEB4Ic6MpMdzLu/lNl0S5awEAMB4Dr0cIZq2vC00dxo3YzlfyZFInfEkCAEVAEMIf9HamWbFcTJqgFcnULshCALB8CEJoqJsD3R3HRqUK5XVkdndkIQBYOHzMwWP42NOsWHbxKfGToxgvBAALhyCEx/O2o3vjuXUXxbl5gty1AAC0IAQhPJGbDckcwW0vkt7LRRYCgMVCEMLTuNqQ7DhuzzXp9f2CiCsMAcASIQjhGRxUJGM4d7ZMmrJb4DFiCAAWB0EIz2ZnRVJiuOvV0qTdghZZCACWBUEIetFwJDmaq+HJ6Ay+BiOGAGBBEISgLzVL1kWy1iwdlc5X83JXAwBgIAhCaAQVQ36KYNva0GHb+XKt3NUAABgCghAah6Xk+8FsTycamcrfqZW7GgCAZkMQQqNRQv41kA1pS4em8bdq5K4GAKB5EITQFJSQLwew8V40LIUvrcIFhgBgxhCE0HTzAtkpnZmIbUJxJbIQAMwVghCaZU4AM8ufCU0RLtxHFgKAWcIyTNBcbz7HWDEkLEXYOZz1d6RylwMA0DgIQjCAV/0YOysSnSakDWN7OiELAcCcIAjBMCZ2ZDhKhqbxW4dy/V2RhQBgNjBGCAYz3pdZEcYlpPP7r2O8EADMBoIQDCnWk/4UwY3J4DNKkIUAYB4QhGBgQ9zphkhu4i5+y2UsVAEAZgBjhGB4IW40bRgXv5PnRTLGB1+2AMCkIQihRfRxoakx3PDtfAVPpnRGFgKA6UIQQkvp5Ux3xXIx2wWtSKZ3RRYCgIlCEEIL8nOgu2PZqDShQkv+3B1ZCACmCJ9N0LI62NOsEeyS0+I/jmLuDACYIgQhtDgvO7onjttQKM7NE+SuBQCgIQQhGIObDcmO5bJKpT/lCLjAEABMCoIQjMRJTXYO547ckl7fJ4gIQwAwGQhCMB4HFUkfwRXclyZlCzxGDAHANCAIwahsOZIczd2qkV7MFrTIQgAwAQhCMDYNR5JjOK1IRqXzNZg9AwByQxCCDFQMWR/BOqrpiO18hVbuagBA2RCEIA+OISsHsx3s6Ygd/H1kIQDIB0EIsmEpWR7G9namkdv427VyVwMASoUgBDlRQhYFs4Pd6dBU/maN3NUAgCIhCEFmlJD/7c9O8GXCUviSSlxgCADGhiAEkzAngHm5MxOaIhSWIwsBwKgQhGAq5gQwf+3JRKQK5+8jCwHAeLAME5iQ17sxHCWDU4Qdw9nujlTucgBAERCEYFpm+DH2ViQmTUgdxgY4IQsBoMWhaxRMzvMdmW9DmGFp/MEb6CMFgBaHIARTFO/FrAjj4nfyWaXIQgBoWQhCMFHDPemGSG7iLn5nCbIQAFoQghBM12B3ujGKe2kXn3QZC1UAQEvBZBkwaYPa0u3DuLidPC+SsT743gYAhocgBFMX6EIzRnAxaUKFlrzcBVkIAAaGjxUwA8850IwR7Hu54vNZDRcwvFlD/pSDVQ0BoOkQhGAeuram++PZlCvSmIzfY+9mDYnfyY/qgD9jAGg6fIKA2ejcmp4ay6WXSNFpPPktBef3ZSPb4bp7AGg6jBGCOfGyIyfHcv4b+Lhd6lqCFAQAA0CLEMyMpy3ZP5LLvU1PlZEtl8XMUkmLaysAoBkQhGBmbtaQGXv5H0PqBrShVyvJR/mC8xpt/E5+dYFYVid3cQBghhCEYE7qxwWj3MSUGJaXyAsdmYsTrMb5MClXpA4/aUOS+a9PiEVY4BcA9IYxQjAbD8+OKS8nKoasi2QnZAqEiG8+x0zuTKp5NqNUTLkiLUwSnNQ03pvGeTKD3ChGEQHgKYzRIuR5ftKkSb169erTp88///lPIxwRLNLlCmlB0B9mx+iykP9tjNCGI/FezH9C2JKJVquHsNYseXWf0OEnfuY+IfmKWIehRAB4HCO1CF955ZWIiIjKysoBAwYMHjy4b9++xjkuWJK+Lo9p2qkYMrt7w+9zDCV9XGgfF3ZeILlYLiVflj47Jr64SxjsTsf5MAneTGuVUSoGAHNgjCDkOC4iIoIQYmtr2759+8rKSiMcFEDH157O7k5nd2du1ZDUIjHlivTWAW13RzrOhxntQz1t0W8KoHQGDsKTJ09WV1c/vMXf39/Gxkb3eOvWrZWVlaGhoYY9KIA+XKzJ5M4YSgSAhgwchBs2bCgtLX14y9/+9jdPT09CyO7du//+97+npaUxDKaqgpx0Q4nxXuT/BrFHb0vJV8RX9wlVPInxoHFeNMaDUeEvFEBJ9A3CI0eOHDhw4NSpUxMmTAgLC6vf/ssvv3z99dd3794dNWrUlClTPvroo8e+/cCBA3/5y19SUlJcXV0NUDWAITx2KPHl3UJkOybOi2IoEUAh9P3q+8knnxw4cCAlJeXUqVP1G69evRoeHu7v7z99+vR58+atWLHise8tLy8fNmyYu7v73//+95kzZx44cMAAhQMYlK89nd2d2RfPnR1nFedFU65I3rgqEUAZ9G0Rbtq0iRASHh7+8MYVK1aEhYW9++67hJDKysqPP/546tSpj77XxsYmKyur/mmHDh0eewhRFKuqqtatW1e/pWfPnl26dHlSSYIgCIJy19/B6bfQ6TtakRd9yYu+pJpnMq9K24rEhcckJzWN8yKxHnRgW2IKQ4m606dKHdbUnb4kKfQLCk6/Uf/3GYZ55v+UZo0RHj58OCQkRPc4JCTk5MmTVVVVGo2m4TE4rk+fPs/cW21tbVVV1fr16+u3UEq9vb2f9Pq6urra2tomFW4JcPotffoMIUNdyVBX8s9A8utdmlpCZ+xjagQa5S4Oby9GuUsyDiXW1tZyHKfYIKytrWUYRrGzDWpraymlSj59SZJYltXz9SqViuOekXTNCsLr1687OTnpHjs7OxNCrl275uvr27S92djYuLi4bNy4Uc/XC4LwaOgqB07fmKc/yJYM8iDz+z8YSlx0Rpx5UJJxKFEURY1Go9gglCTJxsZGsUlACLG2tlbs6VNKVSqV/kGoj2YFoUajqamp0T3WXTVhZ2dngKIATBWuSgSwPM0KQk9Pz0uXLukeX7p0SdekM0BRACYPVyUCWIxmNa7Hjx+/YcOG+/fvE0JWrFgxduxYxbbWQbHqb3BaPNHqPyEsIWTGPsEHNzgFMB9Uz6lH06ZN27x5c3l5uUqlUqvVa9asiY2NFUXx5Zdfzs7Odnd3v337dkZGxpNmhOrj0qVL4eHhhYWFer6+vLzc3t6+yYczdzh9Uz593VDiz4XimTIpogWGEisqKmxtbRU7RlhZWankMcKqqioljxFWV1cbfIxQ3yCsrKysq/t92VM7OzsrKyvd44sXL967d69Hjx7PnJnzdAjCRsHpm8Xp36whaUViyhVpZ4moG0oc40M9mj2UiCBEECr29FsiCPWNLltbW1tb28f+U5OniQJYvDaPDCUuSBKcMZQIYEqwMC+AMdTf4HTJIPbAdSmlSJyxT6jGDU4BTACCEMCoWEpC3GiIG7sw6PcbnL6yW2iJoUQA0Ae+hQLIpv4Gp6cfucFp8SM3OP31jrS9uOFGrUi+OiEq9F5bAAaCFiGA/OqHEqt4NvOhocRxvjTeiwl0oZQQLzv6Zg6vFZl4rwffXwWJTNktDHDFQCNAsyAIAUyI5pGhxJeyfx9K3BjFjcngCSHhTkSQyOTdwgBX+pY/+nUAmgVBCGCKHjuUeKZMCHVj3j0kzn2OSb8hIgUBDAJBCGDq6m9wWlolJV+RyrXS9INW7TVSSFt6pULyskPPKECz4OskgNlop6HTuzKuNvT97oKtlbTxktgnie+9mf+ffCHvplKXpwNoNgQhgNkQJDIpWxjgSuf68znxXJ1Iloeyy0JZlpLX9wttf9BOzhZ+LhTLtXIXCmBWEIQA5qE+BXXjgg4qkhzNfXFcLK2S5gWyhxO5QwlciBtdXSB6/KgdmsZ/fUIseeQaDAB4FIIQwDwcvCEF/3F2jIOKbB3KpRU96BT1saev+jHJ0dyl561e9WPyb0ndN/L+G/i5ecK+a+g4BXgiTJYBMA+D2tJBbRvOi3FUkyWDGt592FFNxvkw43yIIP1+O7cqngzzoHFeNLo9ozbk/YoBzB6CEMBiPXoNxr9OijIkfMQAABtCSURBVC/tEsLcabwXM9KbcbORu0QAE4AgBFCE+mswbteSrFIx+bI0N0/ra0/jvGi8F9PHBddggHIhCAGUxfm3jlNeZA/ekH4uFEelCxxDhrbHOhigUAhCAIXimAcdp18Hk5N3pZQr0mfHxJd3C5HtmDgvOtKbccA6GKAMCEIAIP6O1N+RzglgbtaQtCIx5Yr01gFtd0ca78UkeFM/B3ScgiVDEALA7+rXwajm2f3XpeQrYlSaqGaIbihxiDvl0HEKFgdBCACPYcORqPY0qv2DjtOfC8W5ecLlCincnYnzookdmFZWcpcIYCAIQgB4Bn9H6u/IzgsklyukHcXSz4XiG/uFXs50nA8z2od62qLjFMwbghAA9OVtR1/1o6/6MVU80S0gvDBJcFLTeG8a58kMcsMSwWCWEIQA0GgPLyD8y20p+Yr454PClQppmAcT702HezB26DgF84EgBICmYynp40L7uLDzAh/cvGbpGXH6HqGfK43zZMb60PboOAWThyAEAMOov3nN3VqSUSomX5bmHRHaaeg4XxrvxQS6oOMUTBSCEAAM7NG7fr+ULVTzJAZ3/QaThCAEgJaCu36DWUAQAoAx1Hec3qohu66KyZelObnajq1w12+QH4IQAIzKxfpBx2mNwO67JiVfEUelC1a/3bxmsDu1ws1rwLjwFwcA8rBmSVR7+nUwe+UFbms0205D5x0R3Ndqx2cKqwvEsro/vPjQDenPBwXpj3u4WkXGZAhig60AjYQgBAD56W75vS+eOzXWKs6LplyROvykDUnmPzsmnr0nEUL6u1INS17b93sWXq8mCen8LH+GQa8qNA+6RgHAhLjaNLzrd8Q20Zp90HGaUSK+tk/4Z29yvZokZPCf92OHuCMGobkQhABgiurv+v3VAJJ7U9p6WXz7kHC9WnKzoSMyrapE8X/7IwXBMBCEAGDSGEoGuNIBruynQaSwXFp7Xpp/VGptLQW1QQqCYWCMEADMhoajW6+IG4fwba2J7zrt7Vq5CwKLgCAEAPNwvZrE7+Q/78cObiseSWR97Kjfz9rr1XKXBeYPQQgAZqA+BXXjgpSQgwlcD0fq97P2SiWun4BmQRACgBmoFaSvgxvOjsmK5cb4MINThHP3kIXQdAhCADADXnY02PUxs2OWhbKf9GEiUoVfbiMLoYkQhABg3l7sxPzfIGb4dn7fNWQhNAWCEADMXrwX898Iblwmv70YWQiNhiAEAEswxJ0mx3BT9/DrL4py1wJmBhfUA4CF6OtCM0Zww7cL9+rIDD98ywd94W8FACzHcw40O5b9/Ffxi1/RLgR9IQgBwKL42NO98dya8+LcPEHuWsA8IAgBwNK42ZDdsdzuq9Ib+7FaITwbghAALJCjmqQP5wruS5OyBS16SeGpEIQAYJnsrEhKNFcjkNEZfDUvdzVgwhCEAGCx1CxZH8m6WNPhO/j7WrmrAVOFIAQAS8ZSsiKM7e1MI7fxt2rkrgZMEoIQACwcJeSrAexYHyYshS/GUhXwCAQhACjCnADmtW5MWIpw/j6yEP4Ad5YBAKV4y59xUJHwbcK2GLan02PWsgBlQosQABRkcmfmm4FMTBp/4AbahfAAghAAlCXRm/khnBuVzqeXIAuBEAQhAChQZDv6cyT3Uja/6RIutgcEIQAoUqgbTYvh/pQjfH8OWah0mCwDAAoV6EJ3xXLRaUJZHflLd7QKlAu/ewBQrq6t6d449tvTWKpC0RCEAKBoXnZ0bxy3o1ialSNg8owyIQgBQOlcbciuWO7obWlKtsBjxFB5jBSEJSUlu3btysvLE0X8lQGAyXFQkZ3DuZs10thMoQa9pApjjCCsqamZMWPGtm3bvvzyy8GDB2u1uAk8AJgcDUe2RHNqlsTu4CvwKaUkxpg1am1tnZqaqnscFBR0/vz5bt26GeG4AACNomLIj+Hsa/uEyFQ+dRjnrJa7IDAKI3WNiqKYnp6+ePFiBweHjh07GuegAACNxVKyNJQNc6ODU/jSKsyeUQRDtgh5nn/33XcbbFy0aBEhRBTFvLy806dPOzs7G/CIAAAGRwn5oj/rYi2GJgvpI1hfe9ye28IZIAirq6sppdbW1izLjhkz5vGH4bgPPviAEPL888+npaUlJCQ0/7gAAC1nTgDTSkXCUoTtw9jujshCS6ZXEAqCMGPGjPz8/KKiovT09D59+ui2a7XaqVOnJicnS5I0fvz4b7/9NjQ09NG3379/38rKysbGprq6+vz5823btjXkGQAAtIzXuzGtVSQqld8azfVrgyy0WPqOEXbv3v2bb77RarU8z9dv/O67706dOlVaWlpcXHzw4ME1a9Y89r1FRUVhYWF9+/YNDg6eMGHCgAEDDFA4AEDLm9iRWRbKJezk913DeKHFopLUiN9u69atd+7c2b9/f93T4ODg6dOnT5s2jRCyePHizZs3Z2ZmNrmUU6dO9e3btz4mKaUvvfTSuHHjnvT6iooKOzu7Jh/O3OH0lXz6lZWVGo2GUoW2UaqqqqytrRnGePcD2XeDeTmHW9KPj24n/5XQxj99k1JdXa1SqViW1fP11tbWHPeMvs9mjRFeuHDBz89P99jPz+/ChQvN2ZuNjU3r1q11Q4n1+3zKh50kSUr+KMTpK/n0CSG2traKDUJKqY2NjTGTYJgdSWkljdxJFw1gx/vKnEAMwyg5CFmWbVQQ6qNZQXj//n2NRqN7bGdnV1ZW1py96WbcREVFNWcnAAAtoa8LzRzBDdsu3KsjM/wUGkKWqlm/zjZt2ty7d0/3uKyszNXV1RAlAQCYom4OdHcs+8Vx8fNf5e8gBQNqVhD6+/vn5+frHufn5/v7+xuiJAAAE9XBnu6J49aex7JNFkXfrtGcnJyqqiqe5/Py8srLywcNGmRjYzNz5sxZs2YNGTKE5/nFixevXLmyJUsFAJCfmw3JjuVid/Cv7xf+PZBlFDpQa1H0nTX68ssvl5SU1D/94YcfdJcDLl68eNmyZQzDvP766zNmzGhOKZcuXQoPDy8sLNTz9eXl5fb29s05olnD6Sv59CsqKpQ8WaaystLIk2UeUwNPRqXzLtZ01WDWyriFYNaowSfLNO7yiRaFIGwUnL6STx9BKHsQEkJqBTJxl1ArSD9HcjbGWL/gAQShwYNQoT9KAIBmUrNkfSTrakOHbefvY9kmc4YgBABoIpaS5WFsoAuN2MbfrJG7GmgqBCEAQNNRQr4awI7zYQan8MWVpjLSBI2CIAQAaK45Aczr3ZiwFKHgHrLQ/BhxhBcAwHLN8mcc1CQ8VdgWwwY4KXQek5lCixAAwDAmdWL+PZAZlsbnXEe70JwgCAEADCbBm1kbzo3O4HeWIAvNBoIQAMCQItrRrdHcy7v5DYW4Jal5wBghAICB9WtD04dzw7cL97Vkahe0N0wdghAAwPD8HWl2HBudJpTVkrd7IAtNGn49AAAtwtee7oljvz+HpSpMHYIQAKCltNPQzBHczmLpTzmCiNkzpgpBCADQglxtyK5Y7tgdacpugcfsGZOEIAQAaFmtVWTHMO52rTQmU6hBL6npQRACALQ4DUeShnLWLBmxnS/HUhUmBkEIAGAMKob8GM52ak0jU/nbtXJXAw9BEAIAGAlLyX9C2CHudHAKX1qFyTOmAkEIAGA8lJDP+7GTOjGhycKF+8hCk4AgBAAwtjkBzF97MoO3CcfvIAvlhyAEAJDBa92YL/oxQ9P4QzeQhTJDEAIAyOOFjsyyUC5+J59ZiiyUE4IQAEA2cV7050juhSw+6TIutpcNghAAQE6D3WnaMO71fcKqAmShPLD6BACAzPq40KxYLiZNuFdH3vJH+8TY8BMHAJBfNwe6J45dfEr87BjahcaGIAQAMAkd7OmeOO6/F8S5eQImzxgTghAAwFS42ZBdsdzea9Lr+7Bsk/EgCAEATIijmuwczhWWSy9mC1r0khoFghAAwLTYciQ5htOKZFQ6X83LXY0CIAgBAEyOiiHrIti2NnTYdv5endzVWDoEIQCAKWIpWRbG9nGhEan8zRq5q7FoCEIAABNFCflyADvehwlJ5rOvPWbA8PAtzKgxAAQhAIBJmxPAvNCRxqQKqUV/yMJlZ8V5+QKPCTXNhjvLAACYunmBrB1HEtOFtUNIrBshhKw4J266JG6K4jg0Z5oNQQgAYAbe7cnaceTFbOH/+jHEStp4SdoUxVmzcpdlERCEAADm4bXnWDVHZ+wlfq3F/FFWaqSggaBRDQBgNgSJ9HAQC+6T6XuFOowOGgiCEADAPOjGBTOjhZx4JuWK2HsTf+E+Zo0aAIIQAMAMrDgnbigUdeOCvZ1pbgJXriV9k/h1F9EwbC4EIQCAqfvltrSpUNz80OyYzq1p5gjW155+kCfO3CfUCrLWZ+YQhAAApq6XM90SzTWYHdO5Nc1N5I6O5srqyMBk/jy6SZsKQQgAYAZY+viNrazIugh2ph8TvJX/6QK6SZsCQQgAYPZe9WN2DOP+J1+cnC1gwYrGQhACAFiCQBeaP4qrE8kgdJM2EoIQAMBCtLIiP0Wwf+7ODEpGN2kjIAgBACzK5M5MWsyDbtIqdJPqAUEIAGBpAl3okVGcViIhyXzBPXSTPgOCEADAAtlbkf+Gs3/uzoSk8P9FN+lTIQgBACzW5M7M9mHcR0fQTfo0CEIAAEvW25nmJ3K8RIKS+JN30U36GAhCAAALZ29Ffgxn5wQwEan82vPoJm0IQQgAoAiTOzO7YrkFx9BN2hCCEABAKZ5zoAdGcoJE+ibxJ9BN+hsEIQCAgthbkbXh7NwAZkgKv+wsukkJQRACACjQ5M7M3nju6xPi5GyhUvHdpAhCAAAl6uZAD4zkREKCFN9NiiAEAFAoOyvyw5AH3aTfnVFuNymCEABA0XTdpP86qdxuUgQhAIDSdXOguQmcNUf6JvHH7yiumxRBCAAAxIYjS0PY9wOYiFT+6xPK6iZFEAIAwAOTOzN747jl58TJ2UKFVu5qjMV4QajVaj/66KMjR44Y7YgAANBYfg700EjOUU36JvG/KqOb1HhBuGDBgnXr1h0/ftxoRwQAgCaw4cjXwewHvZhIZXSTGikIT58+feHChejoaOMcDgAAmknXTbrinDjJ0rtJjRGEgiC89dZbCxYsMMKxAADAUPwc6KEEzklN+ibxxyy3m5Qz4L7y8vJee+21h7f06tVr+fLlX3755QsvvNCuXTsDHgsAAIzAmiVfB7NrzotRqfzferGzu1vgFEt9g7Cqqur48eMVFRWRkZEPbz979mxWVlabNm1GjhwZFBSUn5//6HuPHTt2+vTpJUuWFBUVJScnOzk5xcfHG6B2AAAwikmdmH5t6PhMYe91aXko21old0EGpVe25+TkODo6jh8/ftiwYQ9vT0tLCw4OPnPmzDfffDN06FBRfPyY6g8//JCfn3/48OEJEybMmzcPKQgAYHa6tqaHErj2GtJvC//LbYvqJtWrRRgQEHD79u3CwsLAwMCHt3/00Ueff/759OnT6+rq/P3909LSYmNjn7KfkSNHurm5PeUFoijevXu3/qmdnZ2VlZU+FQIAQEvTdZNuLBSHbeffD7CcblIqSfoG+/HjxwMDA7XaB5OHbt261aZNm5s3b7q4uBBC3nrrLVEUFy9e3ORSTp061bNnT3t7+/ot77///htvvPGk11dUVNjZ2TX5cOYOp6/k06+srNRoNJRSuQuRR1VVlbW1NcNYyKdwY5nC6Z8vp5P3W/nYSUv6862tjNo6rK6uVqlULMvq+Xpra2uOe0aTr+mTZUpLS1UqlS4FCSHt2rXLzc1t8t4IIRqNxtPTs7CwUM/XS5Kk5I9CnL6ST58QYmtrq9ggpJTa2NgoNggZhpE9CHvZkdxRZE6uEJHBrotgezkb70+RZdlGBaE+mv6jlCTp4f+HDMM8aYwQAAAsjK6bdEEQM2y72V903/QgdHNzq62tLSsr0z29du0aLpAAAFCU0R2YPXHcygJxdIZQVid3NU3V9CBs27Ztjx49UlNTCSGCIOzYsaPBlRUAAGDxurSmB0dynrak92b+0A2znE2q1xhheXn5u+++e+fOHVEUZ86c6eDg8NlnnxFCPvzww9dff72goCA/P9/a2johIaGFqwUAAJOjZsnXwewQdzEh3Sxnk+oVhCqVKioqihAyfvx4QoiNjY1u+9ixY728vHbu3BkfHz9x4sRnzswBAABLNaoD09uZPr9LyL4qrQhjHdVyF6Q3vaJLrVaPGzfusf/Ur1+/fv36GbQkAAAwSx3s6e5Y7r1cofdm/qcIdoCreUxsNrMGLAAAmDJdN+miYGbkTv6zY6JZjBkiCAEAwMASvZm8RC7psjgqXbhbK3c1z4IgBAAAw/O2o9mxnLcd6b2ZP2Das0kRhAAA0CJ03aRfBzMJpt1NiiAEAIAWlODNHE7kki6LienCHZPsJkUQAgBAy/Kyo3vjuG4OpPdmPue6ybUMEYQAANDiOIYsDGL/FcyMyuDnHRFMqp8UQQgAAEaS4M3kJXA7iyWT6iZFEAIAgPF42dE9cdxzjqT3Zn6/aXSTIggBAMCodN2k3wQzo02jmxRBCAAAMhjpzRxO5NJLpIR0/ras3aQIQgAAkIenLd0dy/VxoYGydpMiCAEAQDYcQ+YFsosHytlNiiAEAACZxXsxhxO5jBIpZjt/vdrYR0cQAgCA/DxtaXYsN6gt7beF33fNqA1DBCEAAJgEXTfp8lB2QpZgzG5SBCEAAJiQqPb0UAKbWSpFp/HXjNJNiiAEAADT4mFLd43gQtxo/y383pbvJkUQAgCAydF1k64IY1/YJcw6wP/jiNjgBTUCeTNHqBUMcCwEIQAAmKjIdvTgSPbYbbL0rPBmzu+hVyeS8ZlCD0eqZg1wFAQhAACYLg9bmjWCm9aVWVkgjs8SCCF1IhmbIYzwpK91M0yEcQbZCwAAQAvhGPL3QDbUjUlM58PLWRdrKc6LMVQKErQIAQDALES1o8dGcSfu0ru1kgFTkCAIAQDALNSJ5C8HxU8DxQBn+l6uISbJ/AZBCAAApq5+XHBGZ/FfA2glTwyYhQhCAAAwaQ1mx1BCFg9kDZiFCEIAADBpJZXSBN8/zBHVZWFbG6r06wiXLl1aVlYmdxWyWb58+a1bt+SuQjarVq26du2a3FXIZu3atUVFRXJXIZv169dfvHhR7ipks2nTprNnz8pdhVH52NMXOz1Iq+Tk5JMnTxJCKCHv9GCUfh3hypUrr1y5IncVslm7du2FCxfkrkI269atU9pnwcM2btyo+yxQpqSkpF9//VXuKmSTnJx89OhRuauQTWpqam5urmH3acZBCAAA0HwIQgAAUDQEIQAAKBqVJKMuBPwUBQUF3bt39/Dw0PP1JSUlrq6uVlZWLVqVySotLXVxcVGpVHIXIo+rV686OTmp1Wq5C5HHtWvXWrdubWNjI3ch8rh+/bq9vb1Go5G7EHncuHHD1tbW1tZW7kLkcevWLWtrazs7Oz1fP3HixI8//vjprzGhICSENGomWG1trWI/BwlOH6ev4NOvq6uzsrKilMpdiDwUfvparZZlWYbRtzvT3d39mV8ZTSsIAQAAjAxjhAAAoGgIQgAAUDQEIQAAKBqCEAAAFM0sV6i/fv364cOHS0pKIiMjO3bsKHc5xnbq1Km0tLTS0lJfX99Jkya1atVK7oqMKisra//+/Xfv3vX09Jw0aZKLi4vcFclAkqQ1a9a0bds2JiZG7lqMavfu3fW31mNZdtq0afLWY3x37txZuXJlUVGRl5fX5MmTnZ2d5a7IeNauXVtZWVn/tGPHjpGRkQbZs1m2CMPCwj799NM5c+bk5eXJXYsMoqOjL1265OXllZaWFhgYeO/ePbkrMqr169cLguDr63vw4MGAgICbN2/KXZEMvv/++1mzZi1ZskTuQoxt9erVP/3008WLFy9evFhYWCh3OcZWVFQUEBBw+PDhDh06XLp0SWl3Hr98+fLF3/ztb3/bv3+/ofZslpdPiKLIMEyvXr3mzp37/PPPy12OsdXU1FhbWxNCBEHo1KnT559/Pm7cOLmLkkfXrl0//vjj8ePHy12IUV29ejUqKmrkyJGnTp3asmWL3OUY1bRp0/z8/P7617/KXYg8JkyY4Orq+s0338hdiMyKi4t9fX0LCgq8vb0NskOzbBHqfymlRdKlICGEUlpbW2tvby9vPXI5f/78rVu3nnvuObkLMbY33nhj/vz5Dg4Ochcij9zc3M8///znn3+uq6uTuxZjS0tLGzNmzKpVq5YsWaLABnG95cuXR0REGCoFiZkGIejMnz+/Xbt2UVFRchdibB988IGHh4e/v/+CBQu6d+8udzlGtXbtWo7jEhMT5S5EHp6enq6urmVlZQsXLgwKCnp4xMji3bp1q7y8/C9/+cvZs2fPnTvXp0+fw4cPy12UDCRJWr169dSpUw28UzMVEBDw3//+V+4qZLNq1SoPD48LFy7IXYgMKisrr169umHDBmdn55ycHLnLMZ6bN2926tSpuLhYkqSFCxeOHDlS7opko9Vqe/bsuWjRIrkLMZ7bt28TQr766ivd03feeWfs2LHyliSL9PR0Z2fnmpoaA+7TLGeNwk8//fT+++9nZmb6+vrKXYsMNBqNRqMZM2ZMcnLypk2bgoOD5a7ISFJTU+/cuZOQkEAIuXbtWmVlZURERFZWltx1yYDjuP79+yuqe9DR0dHGxqZ+LMDf33/Pnj3yliSLFStWvPTSS4a91y6C0Pxs2rTp7bff3rlzp5+fn9y1GBvP86Io6tbc0Gq1x44dM3APiWmLjY319/fXPV61atWxY8e++uoreUsysurqat0NlCsqKnbt2vXee+/JXZHxUEoTExMPHjwYHR1NCDl48KACB8jLysqSkpJycnIMu1uznDU6a9asAwcOnDp1ys3NzcnJ6dtvv+3bt6/cRRlJXV2dnZ2di4tLu3btdFtmz549adIkeasymmvXrgUEBAQHB9vb2+/bt69Dhw6pqanKXI3os88+y8nJUdqsUVdX1wEDBtjb22dnZwcEBGzZskVRC7GdOXMmIiIiPDy8qqrq6NGju3bt8vHxkbsoo1q8ePGqVasMfuGcWQZhQUHB/fv365926dJFOTMnJUk6cuTIw1s8PDzatm0rVz3GV1RUdPTo0Zqamk6dOgUGBspdjmx0XaNKu6HE5cuXjx49Wltb26VLl969e8tdjgzKysoyMzNtbW1DQkL0X5PPYhQWFnIc5+npadjdmmUQAgAAGAounwAAAEVDEAIAgKIhCAEAQNEQhAAAoGgIQgAAUDQEIQAAKBqCEAAAFA1BCAAAioYgBAAARUMQAgCAoiEIAQBA0f4/xj/vDaNSeMgAAAAASUVORK5CYII=", "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.8.4" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.4", "language": "julia" } }, "nbformat": 4 }