{ "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.844661443943 NaN 1.97e-01 0.80 4.2\n", " 2 -7.850303427131 -5.64e-03 2.92e-02 0.80 1.0\n", " 3 -7.850646784427 -3.43e-04 2.99e-03 0.80 3.0\n", " 4 -7.850647500427 -7.16e-07 4.37e-04 0.80 2.5\n", " 5 -7.850647510464 -1.00e-08 1.57e-05 0.80 1.0\n", " 6 -7.850647511576 -1.11e-09 6.89e-06 0.80 4.0\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+gvaeTAAAgAElEQVR4nO3deUBU5f4/8Oc558ww7LKIyKICiqCihqAioCIooKAimmlftcyye7+V9uumZtnX9r37tbw3pbp+szLLDRWVTSVBwQVcUHMpRRFFUwGBGZaz/P6YIlLMGRjmzPJ+/TUcZ+Z8HIv3POf5nOehkiQRAAAAa8XIXQAAAICcEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVEIQAAGDVjBSEGzdujI6OHjVqVF5e3v2e09TUVFpaqvt7CoJggMqsCT4xfYmiKHcJZgafmL7wielLFEWDrwxKjbDW6I0bNyIjI4uLi+vq6kaOHHnq1CkbG5t7n1ZWVhYTE3Px4kUd37a2ttbR0dGglVo4fGL6qqurc3BwkLsKc1JfX29nZ0cplbsQs6FWq1UqFcPg4pyuNBqNUqlkWdaA72mMTz83N3fcuHFOTk5eXl4DBgwoKioywkkBAAB0YYwgrKys9PDw0D729PS8du2aEU4KAACgC2MEoYODQ0NDg/axWq12cnIywkkBAAB00aEgvH79+vr165ctW7Zq1arWx+vq6pYtWzZx4sRFixZVV1eHhIQUFxdr/+jYsWP9+/fvyEkBAAAMiOvIi7ds2ZKens7zfGNj49NPP91y/LHHHmtsbHz22WfXrFkzderU3NxcjuOef/756urqYcOG9ezZs91nrG0mCoao7pkl/bWBdFW1+10BAMB6GaBr9PPPP1+7dm1+fr72x4sXLwYHB1+9etXV1bWhoaFr1675+fn9+vXbv3+/QqEYMWLE/fqjysrKRowY8c9//rP1QQcHh4SEhJYft12WPj8rbRjDqNg/eiC/OCvtrSTfjkKj2gOga1Rf6BrVF7pG9YWuUX3p2zXKMMwD/4Ps0IiwTSUlJUFBQa6uroQQlUoVFhZ2+PDhwYMHx8TEPPC1jY2NP/zwQ+sjHh4erV8Y343cqGdScsTvogWhqamxsfGrX5it5XRdtNDYaPC/iqVpampqxMekj8bGRoVCIXcV5qSxsZFlWQSh7hobGymlCELdNTY2SpKkexCqVCoZgrCyslKbglru7u6VlZU6vtbJyWnTpk1//Zz5A4hSKf7XAfarYbbflasyrorp8dy9F0vhXoIg2NnZyV2FORFFEZ+YXiRJwohQXxgR6oVSavD7CA0fhPb29q2HHWq12t7e3rCneDyQkQiJ22PjphIzEzgbpCAAALSX4b+G+Pr6lpWVtUw9lpWV9ejRw+BnEUSi4cmB69KHJ8RmLFEEAADtZfgRYXR0tCAImZmZiYmJBw8evHLlSutuF4P48qyYfkk8EN/4dbnd+yeEr86Ln0WxsV64GgMAxpObm5uWltbBNxEEQZduDiCExMfHP/HEE53xzh0KwszMzJkzZzY2NjY1Nbm6uiYnJ3/11VdKpfJf//rXrFmzBgwYcPLkyRUrVhi27+7Ls+LmMnFTHNesJv8vhHWxoStPC3N/FEZ1px8NZ3ETBQAYx6FDh3ienzFjhtyFWIUDBw7k5OSYYhDGxsb+8ssvLT8qlUrtg9TU1DFjxpw/f97f39/d3b1DBf5Z+iUx/ZK4OY6zYUkzIeT3+cKd5ZKXPem/sfnNMPbJIAZfrgDACIKDg6dNmyZ3FVZBFMUtW7Z00pt3KAgVCoWLi0ubf+Ti4jJ06NCOvHmb4ryYRB/mru6YuYHMpJ7EzYbMCGCeLhC++Vn8LJLt74I0BACABzOznl0HBWmzR9TNhhBCBrnS/cncvL5M7E5+QaFQzxu5OgAAMD9mFoQPxFAyuw9zbIqiqpEM3MRnXun03RYBAMCsGb5r1BR42pK1o9m916S/FQh9nMlnkayPPa6UAgBAGyxtRNhaTHd6NIUb4k6HpPMrTooCBocAILfvfhHXnn/Avc85FdI/T+L+aOOx5CAkhNhyZHkoW5DEZZSL4en8oV8RhgAgpym9mI0Xpa/un4U5FdKrxcKs3hb+y9mkWMVn3ceZ5iRyLw1mJmbz8wuEO81yFwQA1sqGJRti2U33yUJtCm4fx7m365bonj173rp1q6MlEkIIEUUxNzdXEIQOvk9jY6O7u/tfL/c/aNCg1nfiGZ9VBKHWND/mp6kKFUv6beQfeGkCAKCT3C8LO5iC5PedGQxQIiFNTU1jx47t+H41HMe9+uqrHPdX/ShNTU2iKOfvZMtslrkfFxuyIoKd5if9bb/wwwVx5Qi2lyOaaADA2LRZOG23QIg4pw9DDJGCd5Ek6fvvvz927FhAQMDs2bNtbGy0x3NzcwsKCjw8PGbMmKG9EXzv3r1du3YdMGAAIeTKlSuHDx9OSUlJT08nhPznP/9RKpUTJ0709PTU5aRqtXr16tU3b96MiIhISkrSHlSpfvsrpaenh4aG5uXlnT9/PjY2dvTo0Yb5q3aYFY0IW0R50pIUbqw3M3Qrv7xEaMLgEACMrvW40OApSAh59tln09LSevfuvXXr1okTJ2oPvvnmm88884yHh8fx48cHDx58+/ZtQsiaNWv27NmjfcLp06ffeeed9p1Ro9FER0efOXOmR48eixYtevfddwkhzc3N8+fPb25uJoS89957jzzyyLlz55ydnVNTU/Py8jr+1zQI6xoRtlAwZMEAZmJP+swBIWQT/1kkOwZrdgNAh03I4vX6bi1KZEGhyIskrCudsVe/RUA2xHJdlG3/0aVLl9atW3fp0iVHR8c5c+b4+fkVFBQMGjTorbfeKikpCQ4OJoQkJyf/+9//fuWVV9p8h8mTJxNC5s6dq/uWnD/88IOLi8vq1asJIVFRUcOHD3/++efvek5iYuKyZcsIIbdv305PTzeRQaGVBqGWnyPdEc9tvyzO3SeM9KQfD2cN+HUMAKzQCyGsqM8kXfFN6Wq9qOLIcA8a56XfJTq7+//+Li0tDQ4OdnR0JIQoFIohQ4acPHnSwcHB3t5em4KEkOjo6JKSEr3OSAi5cePGCy+8QAjhOG7NmjWt/+jEiRPDhw/XPu7fvz+ltKysrGfPnq2fM3jwYO0Db2/vAwcO6Hv2TmLVQaiV3IMZ48W8cVTohzW7AaBj9Lq2lFMhpV8S85M5e46k5PLBXSTtfGHH1dfXNzQ0tPyo0Wjs7Ozs7OwaGhokSdLu+qRWq7WjPZZlW7pDNRrNX7+zg4NDamoqIYRh7i61vr6+5bEoik1NTfeOJg27s7yhWOMc4b3sOfJuOJudyP3nnDg6gz9djdsNAaBztZ4XtOXI1rHc/e6paJ/jx4+fPn2aEFJeXl5UVBQZGenn5+fq6rp582ZCiEaj2bBhg/bKZK9evbRDQ0mSNm3apH25jY2Nra2tdhKxNTs7u8mTJ0+ePLll3rG1bdu2qdVqQsimTZt8fX29vb0N9dfpVBgR/mGwGz2QzH1xVhy5nX+0N/N2OGuPjwcAOsG93TH39pF20MCBAx9//PFu3bodOnTo5ZdfDggIIISsWbNm9uzZX3zxxblz5yIiIh599FFCyBNPPBERETFq1Ki6urp+/fppX04pfeaZZ4YMGeLr6/vFF1+0XNL8a3369BkxYoSPj8+hQ4fWr19/76jRNOE3/Z8wlDwVxCT3YBYfEgZt5leOYBN8cKEUAAzpfj2ihs1Cd3f3jIyMY8eO+fj4eHl5aQ/Gxsb+/PPPP/30U9euXX19fbUHfXx8zpw5c/r0aR8fH3d395Zrqu+///6bb75ZX1+vnWvURVRU1N/+9rdz584FBwc7OTkRQlQq1fXr17V3UGRnZ7fcSjFv3rzHHnusg39HQ0EQtqG7HdbsBoBOseeq9GqxkBHPudm08ac2LPl+DJuSy3OUPNrhVdaUSuW9+8La2tqGhobeddDe3j48PFz7WKFQtH6Hlh3XdeTi4jJs2LDWRzw8PLQPWgeqjY1Ny62NsjOPcasssGY3ABhcgBO5Xwpq2XJkSxwX6t6eL9++vr7abpTg4ODp06e3u8j2iY+Pj4iIaMcLvb29W6ev8WFE+Fe0a3Y/GsD8/YCw9mdxVSQb3hVDQwBov54OD/4dYsuR4C7t+VVz+PBh7YOBAwcOHDiwHe/QESkpKe17YW5urmEr0RdGhA/Wx5lmJ3IL+jOTcvgFhUIt1uwGALAgCEKdUEJm92FOpSoIwZrdAAAWBUGoB+2a3d/FsO+fEJOz+bJaTBsCAJg9BKHeojzp0RQuzgtrdgMAWAI0y7SHds3u5J70mQNCWDq/KpId0Q1NNADWxc/P7x//+EdWVlZH3qRlwTP4a1VVVS0LmRocgrD9/B3pznhu+2Vxxl5hFNbsBrAyM2bM6Nu3bwc3wm1oaFAqleayAou8fHx8OumdEYQdhTW7AazWvXem60utVqtUKgShvPDpG0DLmt1fnsWa3QAAZgZBaDCD3WjhRO7R3szI7fySw0KDIHdBAACgAwShIWnX7C5NVVytJwM28VlXMDQEADB1mCM0PO2a3XuuSn/bLwRizW4AANOGEWFnGeNFj2HNbgAAk4cg7ETaNbvzk7jtl8WhW/nDvyIMAQBMDoKw0wU605zxWLMbAMBEIQiNoWXN7gaBBGPNbgAAU4IgNB4XG7I66o81uy/V4UopAID8EITGFv37mt3h6VizGwBAfghCGWjX7C6axB36VQpP5w9cx9AQAEA2CELZaNfsfjOMmbFXmJ0n3GyQuyAAAKuEIJRZcg/mxBTOxYb029icdkbE2BAAwMgQhPJzVpIVEWxWIvflWTFmB/8T1uwGADAiBKGpeMiNFk7kZgYw0VizGwDAiBCEJqT1mt0hWLMbAMAoEIQmR7tm9/9GsE/vFx7eLVzXyF0QAIBFQxCaqAm+9HQq18+FDNzcjDW7AQA6D4LQdLWs2b3tsjhsK3/kJsIQAMDwEISmLtCZ5o7nnuvPJGdhzW4AAMNDEJqB39bsnvrbmt0bLmJZNgAAg0EQmg3X39fsXl6MNbsBAAwGQWhmoj3psSlYsxsAwGAQhOanZc3ugzek8HS+8AaGhgAA7YcgNFf+jnRXAvdmGPPIHmF2nnCrUe6CAADME4LQvLWs2T1wE7/2PNbsBgDQG4LQ7GnX7N4+jl15Gmt2AwDoDUFoIULdaRHW7AYA0B+C0HJo1+w+kcpp1+zOrpAIIRIhbY4QcRUVAEALQWhpvOyods3u+QXCw7uFtefFv++/e6XSgzeklFyMGQEACEEQWqqWNbsXHRIu1EpP5gstQ8Dim9JzhcK/R+CfHgCAEAShBdOu2Z2byNU1k53lYmquIEqk+Kb09/3C5jjW257KXSAAgEng5C4AOleIK81P4r48Ky4sEsJ2Kh1thO3jkIIAAH/AiNDyMZQ8GcRsiuWuNzIXa6Vq7F8BANAKgtAqFN+U/qdEOJzQMLI7DdvC/28pligFAPiNkS6Nbt26de/evRzHffjhh8Y5I7RomRd0EsnWsdzf9wtvHBPyb0hfRrNdlHIXBwAgNyONCG/dujVgwIDs7GzjnA5aHP5Vml/wp+6Yf0eyT/RlztdID23hi7BgNwBYPSMF4dy5c8eNG2ecc0FrXnZk69i7u2PeH8p+HsV+PIyZlMO/dxz31gOAVTPwpVGNRiNJf/q9amtrSyl6FGVzvwbRYR6UEBrqTmfuFQ7ckP4zknWzMXJpAAAmwcBBOH369KqqqtZHNm3a5OHhYdizgKH0dKA/TuDePCYM2cJ/PZqN9sRXFgCwOroG4a1bt44cOXLp0qW4uDh/f/+W4+Xl5WvXrm1oaEhNTR08ePC2bds6p07oLBxDloeyUd2kGXuFeX3pqw+xDNIQAKyJrnOEYWFhy5YtW7RoUXFxccvBysrKsLCw69evq1SqUaNGHThw4H4vz83NXbdu3e3bt9PS0s6dO9fRqsHQ4rxp8WTuwHUpbid/VY1JQwCwIrqOCH/++WeWZUNCQlofTEtLi4yM/OSTT7Q/vv/+++np6fd7B1dX11dfffWvz8LzfOugJYTY2NgMGDBAxyKhI7rZkl0J3BtHhSFb+P+M5BJ9MTAEAKugaxCyLHvvwby8vOnTp2sfx8fHf/DBB/d7eVxcnC5nqa2tffLJJ1sf8fLyWr9+fZtPrq+vRxuOXnT5xP4RSIZ3Yeblk4k+wluDeYV1r7hQX18vdwlmRq1Wi6KI/zF1p1areZ5nGOv+P00fGo1GqVS2GUltsrOze+DH26FmmWvXrrU0wnTr1q2mpqa+vt7e3r7db+ji4lJSUqLjkyVJcnBwaPe5rJCOn1iCAznmRR77kSbs5daPYf0crfqXGv4b0wul1M7ODkGoO4ZhVCoVglB3LMvqFYS66NCnz3GcIPy2rR3P85RSwxYHcnFXke3x3MwAZthW/ocLWI8NACxZh4LQy8uroqJC+7iiosLNzU2lUhmiKpAfJWTBAGZnArf0iDg7T9DwchcEANA5OhSESUlJmzdv1t5Bv2nTpqSkJANVBaYizJ2WpHDNEonK4M/XoJsUACyQrkG4ePHisWPHlpWVvfHGG2PHji0tLSWEPPbYY1VVVePGjXv00UfXrVu3ZMmSziwV5OGkIN/FsAv6M9EZ/Dc/4zIpAFgaXZtl5s2bd+fOnZYf/fz8CCGOjo4HDx7Mzs5Wq9WffPKJm5tbp9QIJmB2HyasK52+W8i+In0WxdpjR2cAsBS6/j7r06dPm8dtbW0nTZpkuHrAdPXrQg9N4pYcFsLS+fVj2EGu6AwEAEuAnl3Qgy1HVkSwSwczcTv5FSdxmRQALAGCEPQ2qzezP5n7v/PilFyhqlHuagAAOgZBCO0R6EyLJnK+9uShLXwhdvcFAHOGIIR2smHJigh2RQQzOYdfXiJge18AMFMIQuiQST2ZI5O53AopPpOv1MhdDQCA/hCE0FG+9jRvAhfZjYZuac6pwMAQAMwMghAMQLu779ejubn7hOUlgoA0BADzgSAEg4n1osWTuaIbUtxOvqIeYQgA5gFBCIbkYUt2JXCTezLhW/md5chCADADCEIwMO22FVviuGcPCAsKhSbcdg8Apg1BCJ1imActSeGuqknkNv6XOxgaAoDpQhBCZ3FWkh9i2f/qzURs49f/goEhAJgoBCF0Iu1l0l0J3LJicXaeoMbuvgBgehCE0OmGuNOSFE6QSHg6f7IKl0kBwLQgCMEYHBXk2xh28SAmZge2rQAA04IgBOOZ3YfZl8T955w4O0+oa5a7GgAAQgiCEIwsuAs9OIlzsSFh6fyxW7hMCgDyQxCCsalYsiKCfSuMScjkV5zErhUAIDMEIcgj1Y8pmsitvyBOyRFuY3dfAJAPghBk08uR7kviBrmR0C38/usYGQKAPBCEICcFQ5aHsp+OYFJzsbsvAMgDQQjyS+7BHJnM7bkqjdvFX1PLXQ0AWBkEIZgEH3u6dwIX5UmHpDdnXcHAEACMB0EIpoKlZHko+20M90S+sKBQaMZt9wBgFAhCMC0x3enRFO78HSk6gy+rxdAQADodghBMTlcV2RHPzfBnRmznMy4jCwGgcyEIwRRpt63YOpZbWITdfQGgcyEIwXSFd6XFk7lKDRmxjf8Zu/sCQOdAEIJJc1aS78ewz/Vnorbz67C7LwB0AgQhmIHZfZjd47m3j4mz84R67O4LAAaFIATz0N+FFk3kJELC0/nS27hMCgAGgyAEs+GgIF+PZpcMYsbsxO6+AGAwCEIwM7P7MAXJ3Jpz4tTdQnWT3NUAgPlDEIL56etMiyZx3nbkoS180Q1cJgWADkEQglnS7u778TBmUg7/3nHsWgEA7YcgBDOW0os5PInbdlmcnCPcwu6+ANAuCEIwbz0c6I8TuIfcSOgWPr8SI0MA0BuCEMwex5DloeyX0eyMvcLyEkFAGgKAPhCEYCHivGnxZK7wujR2J39VjTAEAF0hCMFydLMluxK4eB9myBZ+VzmyEAB0giAEi8JQsngQs34M91QBdvcFAJ0gCMECjepOj6Zwv9yRorbzF7C7LwD8JQQhWCZ3Fdkez80MYIZv5X+4gIEhANwXghAslnZ3350J3NIj4uw8QYNtKwCgLQhCsHBh7rQkhWuWSOR2/nwNLpMCwN0QhGD5nBTkuxh24QAmOoP/5mdcJgWAP0EQgrWY3YfZM4F77zh29wWAP0EQghXp14UemsS52JCwdP44dvcFAEIIghCsjS1HVkSwbwxh4ndhd18AIARBCNZpqh+zL4n76rw4JVeowrYVANYNQQhWKtCZFk7kfO3JQ1v4A9dxmRTAeiEIwXrZsGRFBLsigknJ5ZeXCNjeF8A6IQjB2k3qyRyZzOVWSPGZfKVG7moAwOgQhADE157mTeAiu9HQLc0rT4tbL93dRKPhyZtH0VkDYJk4uQsAMAna3X1HejKzfhQUVBIkMqXXb18TNTxJyeVnBOBbI4Blwv/bAH8Y40WPTuZ6O9O5+4QvzwqkVQrO6YP/WQAsE0aEAH/iYUtyErn3j4t/3y+cDuJO3UEKAlg4BCHA3SghiwcxwV3o1FwS5CJN6okUBLBk+D8coA0anvz7J+GNh/jbDSTg++a159EpA2CxEIQAd2uZF/zvQOH0NEUvB/rGUSFmB38WuzgBWCJjBKEkSVu2bHnhhRcWL15cWlpqhDMCtNtd3TFOCrI3ietuRwOcaPR2fnmJ0CjIXSIAGJQxgrCpqWnPnj2JiYnDhw+fMGFCeXm5EU4K0D7FN6WZf+6OcVKQjHjORUmOpnAnq0jIZn7PVQwNASyHMZplbGxsPv30U+3jb7755tSpU76+vkY4L0A7RHnSKE9610EnBflgGEsI2RjLbr8sPr5PGOVJPx7OuqvkKBEADMqoc4SXLl0qLS0dPny4MU8KYFjJPZjTUzkve9JvY3PaGSxQCmD2DDwinDZt2oULF1ofWbduXd++fQkhN2/enDp16qpVq7p06WLYkwIYmT1H3g1np/sz8wuEdb+IqyLZoC53DyIBwFzoEYQVFRVnz54NDAz08fFpOdjY2JiZmVldXT127FgvL68NGza0+dqqqqrk5OTXXnttzJgxHS0ZwDQ85EYPJHP/Oi1GbefnBTGvhbI2rNw1AYD+dL00Gh0dHRQUlJycnJGR0XKwsbFx1KhRH3300Y8//hgSEnL48OE2X9vU1BQfHx8eHq5UKnNzc69evWqAwgFMAMeQBQOYE6ncL3dIyGY+twIXSgHMj64jwm+//dbb2zsuLq71wY0bN2o0moKCAo7j3nrrrddff3379u33vra5uVk7EMzNzSWEuLi4eHl5tXkWSZKqqqpaH2EYxtnZWcciAWThZUc3xLLbL4vz8oWRnvSj4WxXNNEAmA9dg7BHjx73HtyxY8ekSZM4jiOEpKamLl++nOd57Y+t2dvbv/vuu7qc5fbt2/7+/q2P9OzZMz8/v80n19XV6VQ6/A6fmL7q6+slSddB3mgXciiBvnuK7beBXRbCPxYgWOG0oVqtFgSBUiv8q7eTRqNpbm5mGKxtoiuNRqNUKllW13kIOzu7Bz65Q80yFRUVo0aN0j729vbmef769eve3t7tfkM3N7eLFy/q/nxHR8d2n8s64RPTC6XUwcFB9+c7EvJxJJkdJM0vYDZdIaui2GAra6JhGMbOzg5BqDuWZVUqFYJQdxzH6RWEuujQpy+KYsu/n7YsnucNUBSAORvsRgsnco/2ZqK380sOCw1YiQbAtHUoCD09PW/cuKF9fP36dUpp9+7dDVEVgHljKHkqiClNVVytJyGb+Bw00QCYsA4F4ejRo7OysrSPs7KyIiMjlUqlIaoCsATd7cja0ew/h7NP5gsP7xZ+bZC7IABoi65B+Pnnn8+fP//s2bPr1q2bP3/+sWPHCCGzZs0qLy9/7LHH3nnnnZdffnnp0qWdWSqAWUrqQU9P5fydSL+NzStOiliKBsDU6Nos4+/vL0nSkCFDtD9qV4dxcnI6fPjwV199VVNTk5WVFRYW1lllApgzO468G87OCGDmFwibysRVUWw/K2uiATBlVPfu8M5WVlYWExOje9dobW0teiD1gk9MX3V1dXp1jT6QKJEvzopLDwvzgpjloazK4laiqa+vR9eoXtRqNbpG9aLv7RO6wKcPYDytm2gGbOKz0UQDYAIQhADGpm2iWRHBzi8QHt4t3NDIXRCAdUMQAshjgi89lcr1cyEhm9FEAyAnBCGAbOw4sjyUzUnkvr8gjszgT1UhDAFkgCAEkNlAV7p/IvdUEBO7k19QKNRjdSYA40IQAsiPEjK7D3NsiqKqkQzcxGddwdAQwHgMvEM9ALSbpy1ZO5rdc1X6234h0JmsimS97XEfAkCnw4gQwLSM8aLHUrgh7jQ0nUcTDYARIAgBTI4tR5aHsrmJ3A8XxegM/iSaaAA6E4IQwESFuNKCZG5+EBO3k19QKNQ1y10QgIVCEAKYLm0TzfEpiqpGMmgzv6scQ0MAw0MQApi6brZk7Wj2y5Hs80VCcjZ/pR5xCGBICEIA8zC6Oz2awg1xp0PS+RUnRQFpCGAgCEIAs6FtoslP4rZdFodu5Y/cRBgCGACCEMDMBDrT3PHcgv5MchaaaAAMAEEIYH60TTSnpioaBBK8kd9cJspdEYAZQxACmCtXG7I6iv1mNPvyETE5my9HEw1AuyAIAczbqN+baMLQRAPQLghCALOnYsnyULYgidt+WQxP5w//ijAE0AOCEMBC9HGmOeO5hQOYidn8gkKhFk00ALpBEAJYjruaaDZeRBMNwIMhCAEsjbaJZl0M+2qxmJzNX67DlVKAv4IgBLBMIz3p8SlcVDfmoS38e8fRRANwXwhCAIulYMjiQczBSVzuVTEsnT+EJhqAtiAIASxcbyeak8gtHcxMzObnF6CJBuBuCEIAqzDNj/lpqoIQEryRX3seTTQAf0AQAlgLFxuyOor9LoZ977iYnM1fQhMNACEEQQhgbaI96bEpXFQ3JnQL/95xkcfgEKweghDA6mibaA5N4vZcFcPS+YM3MDQEq4YgBLBSAU40K5F7+SFmcg4/v0C4gyYasFYIQgCrNs2P+WmaQsWSfkJp1FYAABkUSURBVGiiAWuFIASwdl2UZEUEuz6Gff+EmJSFJhqwOghCACCEkChPejSFG+vNhKfzy0uEJgwOwWogCAHgNwqGLBjAFE3iim5I4el8EZpowDogCAHgT/wdaWYC98pDTEoOPztPuN0od0EAnQxBCABt0DbRuNiQgZvRRAMWDkEIAG3TNtF8P4b94IQ4IYsvq8WVUrBMCEIA+CuR3ejRFG6cNzN0K5powDIhCAHgAbjfm2gO3pDC0vlCNNGAZUEQAoBO/B3prgTurTBm+m5hdp5wC000YCkQhACgh+QeTGkq52JDBm76o4lGbGuI2OZBABOEIAQA/TgryYoIdvs49pNT4pgd/He/iI/sEe7axeJUlTRuF48sBLOAIASA9gh1p0UTuUk9mQWFwp1mafoeofn3LDxTLf1XnvDpCJahspYIoBsEIQC0k7aJ5mgKZ8/R/EopbqfQLJIz1dKMvcK6GDa4C2IQzAMndwEAYN687emmOPaHC+JTBWLIdqWrjbgxju3jjBQEs4ERIQAYwMP+TFYCW9MsXVGT7nZIQTAnCEIAMIAz1dLT+4W94/jeTiRwA6/h5S4IQGcIQgDoqJZ5wSAnKT+JdVGSoI18oyB3WQC6QRACQIf8VC3N3Cv8MOa37hgFQ45M5hQMGbAJt0+AeUAQAkCHdLOlP8T+qTvGliMlKRzLkLeOYWVSMAMIQgDoEFcb0tvp7u4YJwXJT+K+/Vn835PIQjB1CEIA6BRdVSRnPPu/J8X/O4csBJOGIASAzuJrT3fEsy8dFnaUY7YQTBeCEAA6UX8XunUcN3cfv68SWQgmCkEIAJ1raFf6XQw3bTdfchNZCKYIQQgAnW6MF02LYpOy+bM1yEIwOVhrFACMYVJPpqaJjNsl5CexPRywBhuYEIwIAcBIZvdhXghhxmdhd3swLQhCADCe5/ozST3o+Ey+rlnuUgB+hyAEAKN6J5wd7EYn5fANWIwUTIMxglCSpOzs7JUrV3799dfV1dVGOCMAmCxKyGeRbFcVfWSPwONWezABxgjC5ubm7OxsR0fHCxcuDB06VK1WG+GkAGCyGEq+Hs02idIT+QK6SEF2xugaVSqVH374ofbxzp07L126FBwcbITzAoDJUjBkYywXn8kvLBRWRLBylwNWzUi3T/A8/9FHH5WXlw8aNCgoKMg4JwUAU2bHkYxx3Ogd/HvHxcWD0K8AsjFkEDY3N0+YMOGug9nZ2YQQSqm/v78oiunp6Xfu3HF2djbgeQHATDkryY54Nnq74GJDngpCFoI89AhCnudra2sdHR057o9XiaJ44sQJhmFCQkIUCsV3333X5mtZlp02bRoh5OzZszk5OVOnTu1g3QBgGbzsaHYiO2qH0EVJHvZHFoIMdApCnuejoqKOHz/e0NBQVFQ0bNgw7fGampq4uLimpiZBELp06ZKZmenm5nbvyysrK6urqwMDAy9cuFBUVLRo0SJD/g0AwMwFONFd8ezYXbyjgib6YtEZMDadvn8xDLN8+fLLly87OTm1Pv7pp5+6uroePXr0+PHjHMelpaW1+fKGhoZly5YNHTp0wYIFH3zwQb9+/QxQOABYkBBXmj6Wm/MjX4BNKsDoqCTp8Z+ds7NzdnZ2y4hw0KBBS5YsmTFjBiFkzZo1aWlphYWF7S6lrKxs0KBBrq6urQ/6+vru2LGjzefX1dU5ODi0+3RWCJ+Yvurr6+3t7eWuwpyo1WpbW1tK2zmq21vJzDuo2DKyaaCLtcShRqOxsbFhGFwT1pVGo1EqlSyra6exnZ3dA5/coWaZ8vJyPz8/7eNevXqVl5d35N0IIc7Ozrt37259RKFQODo63u/5f/FH0CZ8YnqhlOKrg14YhrGzs2t3EE50JKsU4oz9zI9JrJ+jVVwjZVlWpVIhCHXHcZxeQajTe3bkxdpk1j5WqVT19fUdrIZlWX9//w6+CQCYrym9mBsaMnaXkJ/EdbeTuxqwDh36GuLp6Xn79m3t41u3bnXv3t0QJQGAVXs6mHk8kInP5G9jkwowig4FYXh4eEFBgfZxfn7+0KFDDVESAFi7lwcz8d50QhY2qQBj0PXSaFpaWlVVVWNj49dff52Xl/f00087Ozs/99xzSUlJAQEBPM+vWrXqruk9AIB2e38Y+2S+kJLLZ4zjbLAEG3QmXUeEd+7cqaqqWrhwoYODQ1VVlSiKhJCoqKgNGzZkZmbu2bNn27ZtYWFhnVkqAFgRSsjqKLaLks7cK2BlbuhU+t0+0anKyspiYmIuXryo4/O1y9x0akkWBp+YvnDDib7q6+s70jV6ryaRJGfx/k70s0jLHBWq1Wp0jepF39sndIFPHwBMl5IhG+O4I79Ky4qxjS90FgQhAJg0RwXZlcBtvih9WIptfKFTIAgBwNS5q0h2Ivvv0+KXZ5GFYHhG2o8QAKAjvO1pdiI7KkNwVpKpfvgGD4aEIAQA89Dbie5KYMft4p2UdJy3VSzABsaBL1YAYDYGutItY7k5efyRm6bS7g4WAEEIAOYkwoOmRbPJWfzpamQhGAaCEADMTHIP5sNhbMIuoawWWQgGgDlCADA/j/ZmaprI2F1CQTLXzVbuasDMYUQIAGbp7/2YmQF03C6+CptUQMcgCAHAXL02hI31ohOy+Hpe7lLAnCEIAcCMfTScDepCH9nD87jVHtoLQQgAZky7SYUkkcf2CSJaZ6BdEIQAYN4UDNkQy5XXSc8VYmFuaA8EIQCYPVuObI/nCm9Irx/FFVLQG4IQACyBk4JkJnDf/SL+8ySyEPSDIAQAC9FVRbIT2RUnxTXnkIWgB9xQDwCWw9eeZieyY3YKHrZ0gi8W5gadYEQIABYl0JluiWPn7uP3VaKLFHSCIAQASxPela4fw03bzZdgkwrQAYIQACxQTHeaFsUmZfNnsEkFPAjmCAHAMk3qydxpJvGZQn4S28MB84VwXxgRAoDFmtWb+UcIM3aXcEMjdylgwhCEAGDJnu3PTOlFk7L52ma5SwFThSAEAAv3djgb6kYnZfMNWIIN2oIgBAALRwn5dyTrYUsf2SNgkwq4F4IQACwfQ8nXo9kmUXoiX0AXKdwFQQgAVkHBkI2x3IVaaQE2qYA/QxACgLWw40jGOC6/Unr3OK6Qwh8QhABgRZyVJCuB+79z4uozyEL4DW6oBwDr4mFLdiawozKELkoy3R+DAUAQAoD18XekuxLYuJ28k4ImYpMKq4dvQwBgjQa40K1juTk/8gXYpMLqIQgBwEoN86DrYrjU3fyxW8hCq4YgBADrFedNV0Wy47P4czXIQuuFOUIAsGopvZgbDWR8lpCfxHW3k7sakANGhABg7eYHMXMDmXG7+NuNcpcCckAQAgCQpYOZRF86PouvwyYV1gdBCABACCHvDWVDXGhKLt+IJdisDIIQAIAQQighq6JYFyWdsVfAytxWBUEIAPAblpJvYlg1L/33fowKrQiCEADgD0qGbI7jSqukl48gC60FghAA4E/sOLJtHJdeJn1wAgtzWwUEIQDA3dxsSHYi+9lP4hdnkYWWDzfUAwC0wdueZieyo3cIzkoyzQ9jBkuGIAQAaFtvJ7oznh23i3dS0HgfbFJhsfA1BwDgvga60i1juVl5/IHruKPCYiEIAQD+SoQH/c9ILjWXP12NLLRMCEIAgAdI6kE/Gs4m7BLKapGFFghzhAAADzYzgKluJGN3CfnJnKet3NWAQWFECACgk7/3Yx7tTeN38VXYpMKyIAgBAHS1PJSN86YTsvh6Xu5SwHAQhAAAevhwGBvUhT6yh2/GrfaWAkEIAKAHSsjn0aySoY/9KIhonbEICEIAAP2wlHwzmq1QS88WYmFuS4AgBADQmy1Hto3jDt6QlpcgC80eghAAoD2cFGRXAvf9BenjUswWmjcEIQBAO3VVkZxE9tPT4ppzyEIzhhvqAQDaz8eeZif8tknFlF4YWpgl/LMBAHRIH2eaPpb9237hx2voIjVLRg3CkydP3rp1y5hnBAAwgvCudP0Y7uE9fPFNZKH5MV4Q7t69e9iwYRkZGUY7IwCA0cR0p59HscnZ/BlsUmFujDRHWF9f/8EHH8ycOdM4pwMAML6JPZmaZhKfKexLYns6YCNfs2GkEeGSJUtefPFFW1us2Q4AlmxWb+YfIczYXcJ1jdylgM4MOSI8d+7cypUrWx/x9/dfuHBhYWFhc3NzbGzs1q1bDXg6AAAT9Gx/plIjJWfzu8dzjgq5qwEdGCAI79y5o1KplEplt27dpk+f3vqPnJycCCGvv/66k5PT/Pnz9+/fX1paGhgYGBER0fHzAgCYprfC2HpemJTN70zgVKzc1cCD6BSEv/zyy+LFi0tKSurq6m7cuNFyvKKiIjU19cKFC01NTa+//vpzzz0XGRl578tXrFhRW1tLCKmpqQkODg4MDDRU9QAApunjYeyjecL0PcKmWJbDfWqmTad/H0ppQkLC22+/XVVV1fr4iy++OHDgwOvXrx86dGjZsmVnz55t8+WBgYFDhgwZMmRIcHBwSEiIm5ubAQoHADBhDCVrR7G8KM3dJ6CL1MRRSdL136i0tDQ0NLS5uVn7o1qtdnV1PXbsWFBQECFk5syZAQEBb7zxRrtLKSsrCw0NHTx4cOuD3bt3X716dZvPr6urc3BwaPfprBA+MX3V19fb29vLXYU5UavVtra2lKJh8jcankz+UTHQRfogtO2dfNVqtUqlYhiMGXWl0WiUSiXL6nrF2c7O7oEfb/vnCCsqKpqbm1uuc/bt2/fcuXPtfjctOzu7pUuXtj5ia2t7v9/dkiTh17pe8Im1Az4xvVBK7ezsEIQtHAjZmUhidvCf/KxcOriNX8cMwyAI9cKyrF5BqIv2B2FNTU3rfz8HB4fq6uoOVqNQKOLi4jr4JgAApsNZSTITuJEZvD1HFgxA4Jmi9v+rdO3aVa1WNzU1aX+srq7u2rWrgaoCALAcHrZkZwL7Yam4/hdsUmGK2h+EXl5eLi4uxcXF2h+PHDkSEhJioKoAACyKvyPdlcAuLBJ2lqN1xuToFIRNTU25ublFRUWSJOXm5ubn5xNCFArFvHnzli5dev78+W+//bawsHDWrFmdXC0AgLka4EK3juUe38fnVyILTYtOc4QajSYtLY0QMmXKlLS0NFdX1+joaELI66+//sorr0yePNnT03Pbtm0eHh6dWywAgDkb5kHXxXBTd/OZCdxDbmgpMhV63D7R2crKymJiYi5evKjj82trax0dHTu1JAuDT0xfuOFEX/X19egafaCPTwqvHBGPTOb6daEtt09cVUtP5gtbx3K4+/6v6Xv7hC7wkQMAGNX/G8BO7MEM2ypcrvttHFKpIam5wrKHsAaNPPCpAwAY2/oxbLQnHbhZ+LWBVmpISg7/z+HscA+MpOVhpP0IAQCgtZ3x7OgM/qEdnJ+D+K9IpKCcMCIEAJDHt2NYIpHS29L8AuEfB4XMK1J92wuxQefCiBAAQAaVGjI1V9g8mi+6rdx7jbjZ0I9Kham5UlAXGudN47yYUd2pAkMVo0AQAgAYW8u84ECHppE+lGVI8U1pVzzXJJID16Xcq+KSw8K5GmmYB43zYuK8aag7OnE7EYIQAMCoWnfHqNWEEPLiQOaDE+J/5QnfjGbjvGmcN0vCya8NJO+amFshTd0tNvBStCcT503H+1Ife2SigSEIAQCMiiHkkwg2vOuf8uzFgczOcoltdS20q4pM82Om+RFCyIVaKbdCyq2QXjosdFHSOG8a503HeTPOSuOWbqEQhAAARuVhSzxs2xjVjfe971DP35E+FUSfCiKCxB67JeVWSGlnxMd/FFomFEd2p0pMKLYXghAAwGywlAxxp0Pc6eJBjIYn+zGhaAgIQgAAs2TLkZYJxZsNZO81MbdCmrZbVPPSSE8mzpsm+lJfTCjqAEEIAGD23B80oTjWm+mCCcX7QBACAFiUlglFUWKPYkJRBwhCAADLxGBCUTcIQgAAy3e/CcV6XhrlycR50wQf2sPBSjMRQQgAYF3anFBcelhwttYJRQQhAID1euCEYrQntTHkJrimCEEIAABtTyguLxFOVknhXS18QhFBCAAAf/LHhCIhtc3k4A0p96o4v0Asr/9tQjHeh/a0oAlFBCEAANyXo+KPLptralJwXcytkF4rEVWs9jiN82JcbOSusmMQhAAAoJPudnd32Wy4ID2V3xzgZN4TighCAADQW0uXDS+yx29LuRXSeyeEqbvNckIRQQgAAO3HMX902dQ1k6LfJxQv10mju5vHhCKCEAAADMOhrQnF14+KNoxJTygiCAEAwPDMaEIRQQgAAJ3rgROKD7lRRr6rp1iBHAAAjEQ7obh4EJOTyF2ZoVg8kK1qkuYXCJ7fNj+8W0g7I5bVSve+6myN9MZR8a6DDQL57wNCo2CIqgzwHgAAAHpqPaFYqSH5lXdPKMZ6Ma42hBAS6EwrNeKiQ8L7Q3+7ltokkod3C+N9DXNxFUEIAAAy87RtY0JxfkGzv+NvE4ofDmP/cVBYdEh4LYQ0iWT6XmG8L3062DAXNRGEAABgQlomFJtEtvC6lHtVXFYsnKqSRnSjd5rJE3VMbbM0oQdjqBQkCEIAADBNSoaM6k5HdWffGEKqm0jeNTHrivT1eSa2u2TAFCRolgEAANPXRUnG+zIV9eTtUNHbgS46ZIgmmd8hCAEAwNQ1iWRqrjDelz7ZR/xkOK3niQGzEEEIAAAmrSUFtVdEKSErR7AGzEIEIQAAmLSKemm6/596RLVZ2M2W4j5CAACwfH6O1M/x7oVnKCEvhBhmLGfGI8IFCxZcuHBB7irMyYsvvvjTTz/JXYU5eeWVV44ePSp3FebkjTfeKCoqkrsKc/L+++/n5eXJXYU5WbFiRVZWlmHf04yD8PDhw3V1dXJXYU5KSkpqamrkrsKcHD9+/Pbt23JXYU5KS0t//fVXuaswJ6dPn75+/brcVZiTM2fOXL161bDvacZBCAAA0HEIQgAAsGpUktpY6lsWt27dSkpKunHjho7Pr6qqcnR05Dj0++gKn5i+qqur7e3tFQqF3IWYjZqaGltbW6VSKXchZuPOnTs2NjY2Nqa3Wa2pqq2tVSgUKpVKx+dnZGQEBwf/9XNMKAgJIb/++mttba3cVQAAgIXw8fF54Dcz0wpCAAAAI8McIQAAWDUEIQAAWDUEIQAAWDUEIQAAWDWzDMJ33nknPj6+d+/eGRkZctdiHp555plevXopFIrAwMBvvvlG7nLMwJNPPtmtWzeGYfz8/FauXCl3OWbj5s2bgwYNevzxx+UuxAy8+uqrAa3IXY55WLt2bZ8+fSilfn5+Bw8eNNTbmmUQEkLmzp3LMAyWWNORm5tbdnZ2Q0PDihUrnnzyyePHj8tdkambM2fOTz/9JIri119/vWTJksLCQrkrMg8LFixwcHCorKyUuxAzcPPmzZSUlJzfyV2OGdiyZctLL720Zs0anuf37Nnj4+NjqHc2yyB86aWXpk+fbmdnJ3chZuO1114LDAxkWTYxMbFv374IwgeKiopydXXVPggICLh06ZLcFZmBjIyMmpqa1NRUuQsxGy4uLv6/k7sWM/D222//z//8T1RUFMuyfn5+3t7ehnpnswxCaLeLFy+eP39+6NChchdiBkpKSr7//vvnn39epVIlJSXJXY6pq6mpWbRo0WeffSZ3Iebk008/7d69e2Rk5LZt2+SuxdSJonjy5Mk7d+6EhYX16dNnyZIlTU1NhnpzrLZlRerr66dPn75w4cKgoCC5azEDJSUlmZmZJ06ciI2Nxbp0D/T8888/88wzvr6+chdiNmbNmrVgwYIuXbpkZmZOnz599+7dI0aMkLso01VVVdXQ0LBx48atW7cKgjBx4kQHB4dXXnnFIG+OEaG10Gg0EydODAkJefPNN+WuxTzMmzdv48aNp06dOnz48OrVq+Uux6QdOHBg37594eHhxcXFV65cqampOXbsmNxFmbqIiIi+fft269Ztzpw506ZN27x5s9wVmbQuXbpwHLdgwQJvb+8ePXo899xzBmyWxPdcq9DU1DRt2jR3d/e0tDRK797oGf6CQqEYOHBgeXm53IWYtMbGRj8/v6VLlxJCLl++fOvWrVdffRWX+3SH/ysfiGXZwMDATnpzswzC8+fP37lzR61WX7x4sbi4ODAw0NHRUe6iTJckSVOmTLl27dpnn32m/Z7u7e3t6ekpd12mq7a2dvPmzTExMba2tvv379+0aRO+rf+1mJiYmJgY7eOPP/44JycHKfhAX3zxRWxsrL29fU5OzoYNG9A4+kBPP/30J598MmbMGEEQVq5c+fDDDxvqnc1y0e1nn322dTv7qlWrwsLCZKzHxDU3N0dERLQ+8uyzz86ZM0euekxffX39nDlzjhw50tDQEBAQsHDhwmnTpsldlNn49ttvi4qKPv30U7kLMXWPPPLIwYMHNRpN7969X3zxxUmTJsldkakTRfGVV1755ptvbGxsZsyYsWzZMkNtkWaWQQgAAGAoaJYBAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACrhiAEAACr9v8BBlEXCr2FiQ8AAAAASUVORK5CYII=", "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 }