{ "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 `AtomsBuilder`\n", "to build a bulk silicon lattice,\n", "(see AtomsBase integration for more details." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using AtomsBuilder\n", "\n", "system = attach_psp(bulk(:Si); Si=\"hgh/pbe/si-q4\")\n", "model = model_DFT(system; functionals=LDA())\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 `ScfSaveCheckpoints`, which stores the state\n", "of an SCF at each iterations to allow resuming from a failed\n", "calculation at a later point.\n", "See Saving SCF results on disk and SCF checkpoints for details\n", "how to use checkpointing with DFTK." ], "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. This example is a bit artificial, since the norms\n", "of all density differences is available as `scfres.history_Δρ`\n", "after the SCF has finished and could be directly plotted, but\n", "the following nicely illustrates the use of callbacks in DFTK." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "To enable plotting we first define the empty 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 = 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`. The latter is the function\n", "responsible for printing the usual convergence table. Therefore if we simply did\n", "`callback=plot_callback` the SCF would go silent. The chaining of both callbacks\n", "(`plot_callback` for plotting and `ScfDefaultCallback()` for the convergence table)\n", "makes sure both features are enabled. We run the SCF with the chained callback …" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) α Diag Δtime\n", "--- --------------- --------- --------- ---- ---- ------\n", " 1 -7.775304062666 -0.70 0.80 4.8 93.4ms\n", " 2 -7.779131302266 -2.42 -1.51 0.80 1.0 136ms\n", " 3 -7.779315988994 -3.73 -2.61 0.80 1.0 18.7ms\n", " 4 -7.779349970093 -4.47 -2.80 0.80 2.5 23.8ms\n", " 5 -7.779350231296 -6.58 -2.84 0.80 1.0 18.3ms\n", " 6 -7.779350841374 -6.21 -4.07 0.80 1.0 18.5ms\n", " 7 -7.779350855836 -7.84 -4.68 0.80 2.0 22.6ms\n", " 8 -7.779350856114 -9.56 -5.11 0.80 1.5 20.7ms\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+gvaeTAAAgAElEQVR4nO3deVxUdds/8Ot7zplhR9lEEFQgFRcEUcwFFxB3cRfLMjNzy1+3Zd2Ptpl322Npi/lYaptZeVeuYeLGokmuiSmuqaCilriwr3OW3x9jhLgwA8OcWT7vV38Mh+Gci0b48D3nOtcwRVEIAADAXnFqFwAAAKAmBCEAANg1BCEAANg1BCEAANg1BCEAANg1BCEAANg1BCEAANg1BCEAANg1BCEAANg1BCEAANg1CwrC4uLi+fPnq12FRdDpdGqXAEbDq2Z1FEWRJEntKsA4kiSZfDKoBQXhjRs3vvnmG7WrsAjl5eVqlwBGw6tmdSRJqqysVLsKME5lZaUsy6bdpwUFIQAAgPkhCAEAwK4hCAEAwK6ZLwgvXLhw5coVsx0OAADAEIJ5DjN16tQbN24UFRX16tXr9ddfr/N+CnVUKZG3Y83t2UVKkBurV4kAAGCXzLEiPHbs2IkTJzZu3Lht27Y1a9Zcu3atzrv673m5X5J4487uvLkHpUm70QMNAAB1YY4gPHToUHR0NBEJghAVFXXkyJE672p6KDeiBev9s/hX2e0tr/wm7biibOxvpqUtAADYGHPkR35+vqurq/6xm5tbXl5effb2RmeeSIrdIqYOFZaekJJylOQhgpeDKQoFAAD7Y44g9PHxycnJ0T++efNmkyZN6rnDNzrzkiKFrRe9HSk9HikIAAB1V69To+np6Y899lj79u3Hjx9fffvRo0fDw8OdnJzatm27b9++Pn36JCcnV1ZWFhYWHj58uGvXrvWrmYhIlElUlKxC5VSeiWftAACAXanXirCioqJnz55+fn4HDhyo2qgoyoQJE6ZMmfL888+vWrVq/PjxWVlZU6ZMiY6OlmX5zTffdHNzq2fR+uuC5xI009PFmC3iV32Fx0PQMgoAKktNTV2+fLnaVdisgQMHTpkypSH2zOo/vfSzzz5bvXr1nj179B8eOHBg0KBBubm5Go1GUZQWLVqsXLly0KBBsiwzxhi7b2KdOnUqKipq4MCBVVtGjx49atSoGk9bcJTbfpXb0k/y1CpENHM/vyabey1cfrGd7TSOFhcXV11VBWuBV83qiKKo0+mcnJxMtcP3338/IyNjwoQJptohVNm7d+/ly5e//vrrsrIyrVbL87yBX6jVagWhliWf6a8Rnj17tk2bNhqNhogYY+3atTt79uygQYM4rpbTsA4ODi4uLtXPskZERDg63nHP4IrTys6/lOQhnJfD7cq/6ksuDvK7x7kblfwHD3OcTawMdTpdjW8cLB9eNasjiiLP8yZ81QRBCA0NHTdunKl2CFVkWd6wYYOjo6OiKEYFYa3RQw0RhPn5+S4uLlUfuru737p1y5Av5DjO2dk5ISHhAc+Z2IrGh5Dnnd0xn/TkZrRVntsnJaQp3/blna3/TgqO4wx58cCi4FWzOtzfTLXDB5zxgvpjjHHVmHDPpv+59fb2LiwsrPowPz/fx8fHVDt31dRMQb2OnmzHYMHXiXokipdL0D4DAACGMn0QtmnT5vTp0/r3ZpNlOTMzs23btiY/yt0Ejj7tyU9uzfVIlI7cRBYCAIBB6hWEBQUFhw8fvnTpUnFx8eHDh7OysoioU6dOoaGhb7zxRkFBwaJFi9zc3Pr06WOiams3uwP3YTdu0Dbx50vIQgAAqF29gvDYsWPTp0/funUrz/PTp0+v6hv+4YcfMjIy2rRps3379k2bNpn5wsmYIC6xvzA9XVp6wsTvYgwA0BDKRHrwX+6lopkqsU/1aizp1avXb7/9dvf24ODgbdu21WfP9fRwE5Yezw/bLp0rVD7oxvO4gA0AFmxauuQk0Iroe/+u+t+j8uaL8t7h1t8HaKlstsktyI3tHyGcK1SGbRcLdWpXAwBwf2934VKuKE//Isl3LQwXZ8rvHZU+6m7o3QI1BAcH5+bm1rc+IiKSZTk5OVmS6nvHtk6n8/b2Li0tfcBzOnXqdPbs2XoeyHA2G4RE5Kahn/oLwe6s12YxB62kAGCpmruytKH8rj+VqXvuyMLFmfLbR6Ttg4WuPnU8r1VRUVH/qSl6oij279+/rKys9qc+EM/z8+fP12q1D3hOZWWlLJvv2paNr7UFjpb14Jccl3skSpv68529cZIUACyRPgtjtkjT0qWV0TzHTJCCd/vxxx8zMjKCg4MnTZrk4HD7XrSUlJT09HRvb+9HH33U09OTiHbt2uXl5RUWFkZEV65cOXDgwOjRozdt2kREX331lYODQ3x8vJ+fnyFHLC0tXbly5fXr17t16xYfH6/fWDXEIDExMSwsbM+ePWfPno2NjY2JiTHVd2oUW14RVpndgfu4Ozdku5h4Ee0zAGCh9FmYdlWZli4tOmb6FJw9e/Ynn3zy0EMPbd68ediwYfqN77zzzjPPPOPj43P8+PFOnTrdvHmTiL7++uuUlBT9E06fPv3OO+/U7Yjl5eV9+vQ5efJk8+bN582b9/bbbxORKIrTp0+vqKggovfee2/UqFFnzpxp3LjxuHHjUlNTTfB9Gs/GV4RVRrXkgtzYyJ3SqXyaG24X8Q8Aqku+orx7zLiLan7Oyupziiyzzj7sld+M+9onWnETH7r377ecnJzVq1dfvHjR3d190qRJISEhv/zyS2Rk5FtvvXXo0KH27dsT0ciRI5ctWzZ//vx77mHkyJFENHnyZMNn6q5bt87NzW3lypVE1Lt376ioqBdeeKHGfQQjR45csGABEeXl5W3atCk2NtbQ79Z07CUIiSjCi/0az8fvkLKKlGU9eAFpCAANrLM3m9vRuD6X/2bJv9+UHbTk48BmdzBufHLbxvf9VGZmZmhoqLu7OxFpNJouXbpkZmY2btzYwcFBn4JEFB0dffDgQaOqJaLr16/PmTOHiDiO+/rrr6t/6tixY926dbtdW9u2giBkZ2eHhIRUf05ERIT+QUBAwO7du409uknYURASUTMX9ssw4ZFUceh2cW2c4K5RuyAAsGkeDhTXzIgsW5wpb8iWU4cKTZ0oZov0Q5asv15YfyUlJfqBX3plZWXOzs7Ozs76bhr9lFT9RiLieb6qO/TB7Z1E5OLiMmbMGLrXqNWSkpKqphhZlisqKvT7r87w8dkNx+6WRa4a+mmA0LoRi94sXipGKykAWIrq3THVrxfefU9F3WRmZh4/fpyILl++/Ouvv0ZHR7ds2bJJkybr1q0jorKysh9//LFv375E1KJFi4yMDCJSFGX9+vX6L9dqtc7Ozne/iYKzs/PIkSNHjhw5YsSIuw+amJhYUlJCRBs3bvT39w8MDDTNN2NSdheERMQzWtqDn9Ka67FZ+u0GshAA1Hd3j2hVFta4p6LOwsLCnn766eHDh0dFRb300kutWrUSBOHLL7987rnnBg8eHBYWFhYWNnHiRCJ66qmn9uzZ06dPny5dulS/jeHZZ5+Niorq0qXLkSNHDDxomzZtevbsOWzYsJkzZ65YscIy36HFvk6NVje7A9e6EYvfLn4azY9sYYmvDQDYiXePyu8ek5KHCDVu8WruylKG8DFJ0oxfpZXR9T2F6OXllZSU9PvvvwcEBPj7++s3xsbGnjt37tSpUz4+PlXLtWbNmp0+ffrEiRMBAQE+Pj5VZ0cXLlz45ptvFhcXu7m5GXjQHj16zJo1648//mjbtq3+CqVWq7127Zr+3fq2bt1adSvF5MmT9TFsfvYbhEQ0OJBtHyzE75B+v6ksiFT/PDUA2KcQd9o5uGYK6rV0Y2lD+G2XTXPuSqvVdu3atcZGJyenyMjIGhudnZ2joqL0j/UBpqfRaDw8PIw6qIeHx8MPP1x9S5MmTfQPqgeqg4ND1a2NZmbXQUhEHT3ZvuF8/A7pSon0SU9eg5UhAJjd2KAH/epp6cZmtK1jw0xAQIC+GyU0NPSRRx6p207qrH///tVD1HD+/v4ajfm6Ge09CInI35ntHipMSJOGbhfX9hMaPWjuDwCANTlw4ID+gf4SoJmPrr/1sA527txp2koeDCsgIiJXDW3sz7dtzKI3ixfRSgoAYE8QhLfxjJZ0559tz0Vvlg5dRxYCANgLBOEdpoVyn/fih+0Qvz+PqaQAAHYBQVjTwACWPFiYd0hekFHft90CAADLh2aZewjzZPuGC/E7xJxiaXk0WkkBwCCCIKxatcrMjR52Ii8vr2psqckhCO/Nz5l2DxMeS5MGbxPXxQmN0UoKALWZPn26Wu+oZw8CAgIaaM8IwvtyEWhjf/4/GVLXn8QtA/hWjfCmvgDwIG5ubp07d1a7CjAazvo9CCNaEMn/O4zr/bO45y+0kgIA2CAEYe2mhnKr+wpjU8Q1aCUFALA5CEKD9G/GUoYIr/yGVlIAAFuDIDRUBw+2N15IylEm/yJVYmUIAGArEIRG8HOmX4YJ5RLFbhFvlNf+fAAAsHwIQuM48rQmho9rxronimcK0D4DAGD1EIRG07eSzgvn+v4s/oJWUgAAK4cgrKMpbbhv+goJKeK353DBEADAiiEI6y6uGUsZIsw/LC/IkLAwBACwUgjCemnvwQ6OEFKuKo+mSuW4sQIAwAohCOvL25GSBws8o35J4nW0kgIAWBsEoQk48PRtDN+/GeueKJ7Ox1lSAABrgiA0DX0r6SsRXN8t4q4/kYUAAFYDQWhKk1tz6/oJj6VJq8+ilRQAwDogCE0suinbPYx/+3d59j5JxsoQAMDiIQhN7yF3tne48PtN5ZFUqUxUuxoAAHggBGGD8HKgHYMFB55ik8TcMrWrAQCA+0MQNhQHnlb35QcGsO6J4im0kgIAWCoEYQPSt5K+E8X1SxLT0EoKAGCREIQNbnywvpVUXH4KraQAABYHQWgOPXzZnmHCR8fRSgoAYHEQhGYS4s72DheO3lISUqVStJICAFgMBKH5eDrQjsGCM0+xSeI1tJICAFgGBKFZaTn6ui8/qgXXI1E8iVZSAAALgCA0N0Y0N5xb2JWL2SIm5SALAQBUhiBUx7ggblN/4ek94icn0UoKAKAmBKFqujdh6fHC0pNoJQUAUBOCUE3BbmxvvJB5SxmbIpXcp5UUEQkA0KAQhCrzcKBtgwV3DYX8qJt3SKrx2TXn5Um7am4EAAATQhCqT8vRqj78/2vLf3BMnvXrP7H3Y5b84gHpf8LxGgEANCBB7QLgtlc7cX7ONCNdKtTRJ53pxyz5uf3SjsFCBw+mdmkAALYMQWhBprThmjjR6J3SiZvavyqknYOF9khBAIAGhtNuliW+OfdqJHcsj2vXmEMKAgCYAYLQsvyYJa84JX/WXbf7T3ki2mQAABoeTo1akKrrgi2Ecr9GTvHbxUZa+r8evNp1AQDYMqwILcV/z8vP7Zd2/t0dMyiALevJrTglv3TXPRUAAGBCCEJL4ePIanTHTAvlXwznfsiSi3Uq1gUAYOMQhJYirhm7uzvmf7vwgwK4hFRRxERSAICGgSC0dEu685JMs/biBCkAQINAEFo6DUfr4oT9ucpHx7EqBAAwPQShFXDTUNJA/sPj8sYLyEIAABMzXxCeP3/+5MmTZjucjWnmwjbG8dPTpf25eDsKAABTMtN9hMOHD79161ZhYeGxY8fMc0TbE+nNVvURRu0U0+OFEHcMnQEAMA0zrQgTExPXrFljnmPZsCGBbEFnPn6HlFehdikAALYC1witzPRQrn8zNipZrMTlQgAAUzDxqdE5c+YUFhZW3/Lee+95enqa9ih27sNu/JhkaUa69GVvTF8DAKgvEwRhfn5+RUWFr68vEU2ZMkWnu2MOipubW/0PAdVxjNbE8DFJ4tu/y69EYE0PAFAvhgbhzJkz09LSzp49+/33348bN06/UZblmTNnrl271tHRsV27dhs2bGjfvn2DlQr/cBJoU3+he6LY3JUmPoQsBACoO0N/h3br1m3VqlWtWrWqvvGnn35KSUnJzs6+fPmyVqtdtGjR/b583rx5o0eP1mq1Xbp02bx5c71KBiIiaupEWwbyLx6QUq/ihgoAgLozdEU4adIkItJoNNU3fvfdd0888USjRo2IaNasWc8+++ybb755zy9fuHBhrYcoLS39888/O3XqVLVl6tSpEydONLBCW1JcXGzI0wJ5+qobNyFVSYqtbOWGOFSZga8aWA5RFHU6nSiKahcCRigrK9NqtTxvaIeEo6NjjeS6W72uEWZnZ48dO1b/uHXr1pcuXZIkyfD6anB2dvby8vr888+rtrRo0cJuLzEa+I0PcaP/leRxe7h9w4UmTg1dFNTCbv+5Wil9EDo54SfHmgiCYFQQGrTP+nxxUVFR1b8hFxcXSZKKi4v1C8S60Wq1nTt3rk9Jdmhya+5coTJsh7hrqOCMN1oGADBSvfosfHx88vPz9Y/z8vIcHR3d3d1NURUY560ufJtG7MndkozzowAARqpXEHbs2PHgwYP6xwcPHuzYsSNjGP2lAkb0RW/+Rrnyym94tyYAAOMYeirtl19+uXbtWkFBwf79+4mob9++Pj4+M2bM6N2794ABA3x9fd96663XX3+9IUuFB9FytDZO6JEoNneVZ7bFDRUAAIYyNAh3796dmZnZrVu3nJycnJycDh06+Pj4hIeHf/vtt0uWLCkvL3/++eefeOKJBq0VHszLgbYO4nttlkLc2YBmWJoDABiEKYqlXFa6cOFCTExMdna22oWor6ioqM79h+l/KaOTxZ1DhHBPZKFZ1edVA1Wga9QaGXv7hCFwDs3WRDdly3ryw7ZLl0ss5U8cAABLhiC0QeOCuJltuRE7pWJd7U8GALBzCELb9HIE19WHjU8VRbxbEwDAAyEIbdbH3XmdTLP24oYKAIAHQRDaLA1H6+KE/bnKkuNYFQIA3BeC0Ja5ayhpIP/BcXnjBWQhAMC9IQhtXDMXtjGOn54u7c9FEykAwD0gCG1fpDf7qo8wLkW6WIwsBACoCUFoF4YGshfDuCHbpPxKtUsBALAwCEJ7MbsDF9eMjdopVuJyIQBANQhCO/JhN76xls1Ixw0VAAD/QBDaEY7Rmhj+ZL7y9u9YFQIA3IYgtC9OAm3qL3x+Rv7mHLIQAIAIQWiHmjrRloH8iwektD/RRAoAgCC0S+0asx9ihQmp4pkCZCEA2DsEoZ3q68fejuKHbJNyy9QuBQBAVQhC+/VUa258MIvfIZaKapcCAKAeBKFdezuKb92IPblbknGKFADsFYLQrjGiL3rzN8qVV3/DzYUAYKcQhPZOy9HaOGH9BeXTU7ihAgDsEYIQyMuBtg7i3zoi77iCM6QAYHcQhEBEFOzGfojlH08Tj91CFgKAfUEQwm3RTdn/9eCHbpculyALAcCOIAjhHwnB3Iy23IidUrFO7VIAAMwFQQh3eCWC6+rDxqeKEpaFAGAfEIRQ08fdeZ1Ms37FDRUAYBcQhFCThqN1ccK+XOXjE7ihAgBsH4IQ7sFdQ0kD+fcz5U0XkYUAYOMQhHBvzVzYhjh+2h7pQC6uFgKALUMQwn119mZf9RHGpkgXi5GFAGCzEITwIEMD2Qth3JBtUn6l2qUAADQMBCHU4rkOXFwzNmqnWInLhQBgixCEULsPu/GNtWxGOm6oAAAbhCCE2nGM1sTwJ/OVd37HqhAAbA2CEAziJNDGOGHlafnbc8hCALApCEIwlJ8zJQ3iXzggpf2JJlIAsB0IQjBCu8bsh1jhsTTxjwJkIQDYCAQhGKevH3urCz94m5RbpnYpAACmgCAEoz3VmksIZmNSxHK0kQKA9UMQQl28E8W3cGVP7JJknCIFACuHIIS6YESf9+KvliqvHcaqEACsG4IQ6siRp58GCOuyleWncEMFAFgxBCHUnZcDJQ3k3zwi77yCM6QAYK0QhFAvIe7sh1j+sTTx2C1kIQBYJQQh1Fd0U/Z/Pfhh26XLJchCALA+CEIwgYRgblooN2KnVCKqXQoAgJEQhGAar3biorzZ+BRRwrIQAKwKghBMZmkPvlKmFw7ghgoAsCYIQjAZDUfr4oTUq8rHJ3BDBQBYDQQhmJK7hhL78+8elTddRBYCgHVAEIKJtXRjiQP4aXukA7m4WggAVgBBCKbX2Zt91UcYmyJdLEYWAoClQxBCgxgayF4I44Zsk/Ir1S4FAOCBEITQUJ7rwPVrxkbtFCtxuRAALBiCEBrQR934Rlo2Mx03VACA5UIQQgPiGK2J4ff8Jb/0W80sPFugpFzFFUQAUB+CEBqWs0AJwfz7x+SlJ//JwrMFSmySdAmtNABgAQS1CwDb91YXrkJWnt8nN3Vi44I4fQq+EsFNbo2/wwBAfQhCMIdFXfliHU1Ik0RJ+Z9DyisR3Iy2SEEAsAj4ZQRm8mlPPsaPTdwtt3QjJ4FO5isyzowCgAXAihDM5GyBciqf+vix47fop4vy27/TtVIl0pt18WZRPizKhwW5MbVrBAB7ZKYgPHPmTHp6uiAI8fHxnp6e5jkoWI6q64Iz2nKLM+XPTstpQ3g3LTt2Uzl8Q/n5kvKfDPlKqdLBg3X2ZtFNWbQv5+esdtEAYB/MEYRlZWWzZs0aMWJEUVFRly5dDh065OXlZYbjgoWonoJE9GIYR0SxSVLqED66KYtuenshmF9Jv11X0q/Jq8/K/2+vxDPq7M30/3Vvwnk7qvktAIANM0cQOjk5JScn6x8fPnz44MGDgwcPNsNxwUK8/bv8n87cU9V6RF8M43QyLTkhv9eVr9rYWEtxzVhcs9tbrpYqh28oh28oK0/LT+6WHPjbuRjty3X3ZS44qQ8AJmLWXyf5+fmZmZnh4eHmPCio7qs+/N1X/14K5x7cK+PvzPybs/jmtz/MKlLS/1IO31AWZEhHbynNXVjVerGLD3PkH7gvAID7M3EQTp8+/eLFi9W3LFu2LCQkhIgqKioeffTRl19+2d/f37QHBQt3vx4Yo3pjgt1YsBt7ohURkSjTmYLb68W12fLvN5UWrrdDMbopi/Bi9wheAID7MC4I8/LynJ2dHRwcqm88d+5cfn5+eHi4RqNZuHChKIrVP+vh4UFEOp1u/PjxQ4YMefLJJ+tdM9g7gaP2Hqy9x+1c1Ml07Nbt9eLK0/LlEiXM85/1YjsPhlgEgAcwNAgnT568adOm/Pz8Tz/9dMaMGfqNsixPnDgxPT3dz8/v+vXrycnJQUFBd3+tJEmPPvpoQEDA0KFDs7KyfHx83NzcTPYdgN3TcLcvH+o/LNLR0ZvK4RtK8hXl3aNoRgWAWhgahBMnTnz99dcnT55cfeO2bdv27t2bmZnp7u7+7LPPvv7666tXr777a0tKSjiOy83NnTdvHhE988wzffv2rXflAPfmpiE0owKA4ZiiGDHeIyYmZvz48VUrwsmTJ3t7ey9atIiIjhw50rNnz+LiYo6r47SakydPRkRENGvWrGrLzJkzn3nmmbrtzaoVFxe7urqqXYVt+rOMfs/jjtxiR26xAzc4B54iPOROnkp3b+VhH8WZr/u0G7xqVkcURZ1O5+TkpHYhYISysjKtVsvzhjbIOTo6CkItS756NctcunQpKipK/zgoKKisrOz69eu+vr5125uzs7Ovr29KSor+Q8ZYs2bNtFptfSq0Uoqi4FdqA2nlSq18aNzfH1Y1o757Sjn6q9HNqNsuK339bj+t6lU7ma8IjFo3wqVJS4cgtEY8zxsVhIaoVxCWlZVVNc44OjoSUUlJSb2qEYTg4OD67AHAKPVsRv32nPxhpvLTAKEqMk/mK3FJ4pLuPIIQwFrUKwh9fX1v3bqlf3zz5k0i8vPzM0FRAGqoQzPq1334J3dLw3eIP/UXiOh0vjJwq/RuV35ckNWPs1fudX/LPTcCWLt6/bhGRUWlp6frH6enp3fo0AEnGcBm6JtRZ3fgVvflT4wVLk/QLIzig91Y8hUlIUXyWK2L3izO2S/182cuAhuxUzyWx/pvld6J4iY+ZPUpqJOpZ6J4+IZSY2NCivTTRVmtqgAaiKErwpSUlPPnz1+9ejU9PZ3juAEDBrRs2fKpp55avHjxokWL2rRpM3fu3Pnz5zdorQAqqtGM+lcZHbou/3Zd+TFbPnRdKRZZ7F8O4R6087KS/pfkJJAjT84Cc+DIVUMajtw0JHDkrmE8o0Za4hg11hJj5KFljFFjC7sUruHotU780O1i4gChqw8jIkmhSbulYp0yMABTfMDWGBqE165dy8rKGjVqFBFlZWXprwU2bdp0165dH3/88f79+996663HH3+8ASsFsCRNnSi+OaefAHcyX+m7RfTUKpUy9fHjJIVKRaqQqERUSkXKKSGdTEU6EmUq1MmSQgWVJCuUX0mKQnmViqJQfiVxjBppiWfkrmECR24a0nDkqiEtRy4Cc+DJWSBHnpwEcuKZI0/OAjnw5CKQlidXgVXLWuI5aqRlHFFjh7qfyRwcyL7qLQzfISYOEDp7s4m7pLwKZWN/AdPswPYYd/tEg7pw4UJMTEx2drbahaivqKgIMwesxel8RX9GNL5JybMZztfKlJ/6C07GX3yXFSqoJEmhQp0iylSkI51MxTqqlKlEVCokKhWpXKIykcokpVyqylqqlKhYVKplLUkyFVQqMlF+BSlEHg7EiBpr2T9ZqyWBkZuG/ZO1GnLgyFlgt7P27xVtxnX5w+Nye0+ukUZZHyc42FYKomvUGhl7+4QhMMMfoO6qUnDiQ1xREa3qwz+5WxqxU6xDFnKMPByIiLwd717F1atDJa+CFKL8SuWfrK0kUaEinfJP1uqoQqa8SuV21opULlGpKJdLJCps71/yTwNsLQUBqiAIAepuxq/Swijusb+7Y3hGq/rwE3dJy0/Lz3ewlJYZfb56Ohidr5JCE3dJPX0VL0d+dLK4a6jQwxdNo2CDEIQAdbdzsKC5M+94Rt/0tYV3v9CnoP66oANPWYVK/61i2tDbvTMAtsRS/mgFsEaae/0A8Yw4Kw8LUabxqVJhpbKpv+DIEyPaMpB319CArTXvqQCwAQhCAKiJ52hgM1a9O8bDgbYNFhjRrQoEIdgaBCEA1MSIpoZyNbpjwj3Zwq78CwfkUi4je4EAABisSURBVPE+XwZgnRCEAGCo6aFcZ282LV1SuxAAU0IQAoARPunBn8xTPj+DQWtgOxCEAGAEJ4F+iOVf+U1C1wzYDAQhABinVSO2MpofmyLdrFC7FABTQBACgNFGtOBGtmBP7haxKgQbgCAEgLpY1JUvqKT3juJiIVg9BCEA1IXA0Y/9hKUn5Z1XsCwE64YgBIA6aupE3/Tln9wtXSlBFoIVQxACQN3F+LFn2nHjUiQdTpGC1UIQAkC9vBzB+TqxV37DXfZgrRCEAFAvjOjL3vz6bGV9NlaFYJUQhABQXx4OtD6Of2avdDofFwvB+iAIAcAEIrzYfyL5hFQJI7nB6iAIAcA0ZrTlIr0wkhusD4IQAEzm0578iTzlC4zkBquCIAQAk3ES6MdY/mWM5AargiAEAFPCSG6wOghCADAxjOQG64IgBADTe68rn19Ji47hYiFYAQQhAJiehqO1/YSPT8jJGMkNFg9BCAANQj+Se9Ju6WopshAsGoIQABpKjB+b2ZYbm4yR3GDREIQA0IBe6cQ1wUhusGwIQgBoQIzoK4zkBsuGIASAhoWR3GDhEIQA0OAivNgCjOQGS4UgBABzmNmWi/Ri0zGSGywPghAAzOTTnvzxPOXLP3CxECwLghAAzEQ/kvulQ1IGRnKDJUEQAoD56Edyj0mRbmEkN1gMBCEAmNWIFtyIFmwSRnKDxUAQAoC5LcJIbrAkCEIAMDeM5AaLgiAEABVgJDdYDgQhAKgjxo/NaMuNS8FIblAZghAAVPNqJ87Hkb2KkdygKgQhAKhGP5J7Xbay4QJWhaAaBCEAqMnDgb6P5aenS2cKcLEQ1IEgBACVRfmwNzrzCSkYyQ3qQBACgPpmtuU6YSQ3qARBCAAWYVlP/shNjOQGFSAIAcAiuAi0MQ4juUEFCEIAsBStGrGl3TGSG8wNQQgAFiQhmBvenD25W8KqEMwGQQgAlmXxw3xepbIYI7nBXBCEAGBZNBz9N4b/8Li0+08sC8EcEIQAYHECXNi3fYUJaRjJDeaAIAQASxTrz2a05SakSSJOkUIDQxACgIV6JYJzFujVw7jLHhoWghAALBTH6Lu+wtosjOSGhoUgBADLpR/J/cyv0vlCXCyEhoIgBACLFuXD5nfiRydjJDc0FAQhAFi6Z9pxEV7sX/twsRAaBIIQAKzAJz35/bkYyQ0NAkEIAFbARaANGMkNDUMwwzFkWf7mm2/Onj3r4uIyYcKEFi1amOGgAGBjWjdiS7vz41OlQyOFxlq1qwEbYo4VoSRJN2/ejI6ObtKkSZ8+fYqKisxwUACwPQnB3JBA9sQujOQGUzJHEGo0mjlz5gwaNGjKlCne3t5//vmnGQ4KADZp8cP8rQrl/UxcLASTMcepUSLS6XSvvvrqH3/8ERsb27p1a/McFABsj4aj72P5rj+JXX1Y76ZM7XLAFpgyCCsrK7t3715j4+HDh4mI5/kxY8acOnVq8eLFL7zwgq+vrwmPCwB2RT+S+9FU6dBI3t8ZWQj1ZUQQ5ubm5uTkhIaGuri4VG0sLS3dtWsXYywmJsbR0VEfe3fjOK5r165du3bdtWvXnj17xo4dW9/CAcCOxfqz6W25CWlS8mBBQPM71I9B/4J0Ol3z5s1btmzZtWvX48ePV22/du1aWFjYRx99tHjx4k6dOt26deueX37x4sXU1NSrV6/u2bPnl19+6dSpk2lqBwA79ipGcoOJGBSEPM///PPPhYWFrq6u1bd//PHHkZGRO3bsSElJeeihhz755JN7frlGo9m4ceOMGTO+/vrrb7/9NiQkxASFA4B9w0huMBWDTo1yHNexY8e7t2/atOmNN97QPx4/fvySJUteffXVu5/m7++/dOnSWo9SWlp68eJFxm6f8WeMvfnmm7NnzzakQhtTUlJS9f8BrAVeNfPTEH3ZjSXs0bRy1AW5Gn1LhSiKOp1OkrCmtCZlZWVarZbneQOf7+joKAi1JF29mmUuX74cGBiofxwYGHj58uX67M3Z2blFixbZ2dn12YltUBSlxuIbLB9eNVX0caX5JfLje7n9wwUnI3+f6YPQycmpYUqDBsHzvFFBaIh6XWUWRbGqGo1GU1lZaYqSAACMMAsjuaF+6hWETZs2vXHjhv5xbm6un5+fKUoCADDOJz35fbnKVxjJDXVSryCMjo5OSUnRP05JSenVq5cpSgIAMI5+JPc8jOSGOjH0nPrixYtv3LhRUVGxbNmyjRs3zp0718PD4/nnn+/du7enp6ckSatXr96/f3+D1goAcD9VI7l/Gyk0wkhuMIahK0J3d3cPD4///Oc/7du39/Dw4DiOiCIiIvbs2XPz5s3CwsK9e/e2bdu2IUsFAHiQhGBucCCbiJHcYCSmKJbyb+bChQsxMTHoGiWioqIiNzc3tasA4+BVswQ6mWK2iCNbci+G1f5XPrpGrZGxt08YArOJAMB26Edyf5Ap/fKXpfyJD5YPQQgANqVqJPfVUmQhGARBCAC2JtafTQvlHkuTRNxPAQZAEAKADXqtE+ck0GsYyQ0GQBACgA3iGH3bV/ghS9mIkdxQGwQhANgmTwf6Ppaf+at0vhAXC+FBEIQAYLO6+rDXOvGjk6UyUe1SwIIhCAHAlmEkN9QKQQgANg4jueHBEIQAYOOqRnIfuYmLhXAPCEIAsH2tG7GPu/MJKVIB3jUV7oIgBAC7MD6YGxTAntiNkdxQE4IQAOzFB934m+XKB5m4WAh3QBACgL3QcPTfWH7RsXuM5L5RrkpFYBEQhABgR9w0TJTZ6J13jOQ+X6h03iQeyMVJUzuFIAQAO9JYS5sH8mWSMmjr7ZHcF4qU/luluR25h5swtasDdSAIAcC+dG/CdgwWzhQoj+9WLhZTbJL0Yhj3TDv8MrRfeO0BwO709GWb+gubLtHDScK8cKSgvcPLDwD2qHUjauzAikVan63kVahdDagKQQgAdkd/XXBBBG3pJ+25Jj/0o7guG/dU2C8EIQDYlwtFiv664IxQ6uOrpAwRZEWZs1+O3yFeKUHjqD1CEAKAHSmopD5bpOrXBbs3YVsGCuWSEuTKIjeJK0/LCEN7I6hdAACA+TTSUuIAPtzzjjsleviyXcOE1u5saig3ZY+05rz8WTTfqhHuprAXWBECgH2pkYJ67RozgaMwT7Y3XhjVguu5WXz3qIyxpHYCQQgA8A+Bo9kduP0jhB1X5F6bxZP5CEPbhyAEAKgp2I0lDxGebM313izOOyRVoqXUpiEIAQDugRFNC+Uyx2j+KKAum8RD17E0tFkIQgCA+/Jzpg1x/GuduBE7xdn7pBJR7YKgASAIAQBqMS6IOzFGUy5R+AYx9SqWhrYGQQgAUDsPB1oRzX/cnZ/8i/TELglT2WwJghAAwFBDAlnmGMHDgcI2iBsvoIXGRiAIAQCM4K6hJd3572P4lw7JCSnSdby1vfVDEAIAGC26KTsySgh2p/brdCtPY2lo3RCEAAB14STQwih+x2BhxWl56HYxBwO7rRaCEACg7iK82IHhQu+mXJdN4pLjMiZ2WyMEIQBAvQgczQ3n0ocJGy/KfbaIZwoQhlYGQQgAYAKtGrG0ocLEh7heGNhtbRCEAACmoZ/KdnCEkHJV7rJJzLiBMLQOCEIAAFNq6cZ2DBZejuCGbBfnHZIqJLULgtogCAEATG9cEHd0tOZ8IXVYL+7+E0tDi4YgBABoEL5OtLYf/0E3buIuaXq6VKxTuyC4DwQhAEADim/OHRsjEFHHDeLOK1gaWiIEIQBAw2qspRXR/Kc9+al7pIQU6SYGdlsYBCEAgDkMDGAnxwrB7tRxvbguG1PZLAiCEADATJwFWhjFr4vjXz8sx+8Qr2Aqm2VAEAIAmFX3Juz30UK0Lxe5SVx5GkPZ1IcgBAAwNw1Hc8O55MHC52fkmC3iWUxlUxWCEABAHWGebG+8MKoF1xNT2VSFIAQAUI3A0ewO3L7hwo4rcq/N4sl8hKEKEIQAACoLcWfJQ4QnW3O9N4vzDkmVaCk1LwQhAID69AO7M8dozuRTl03ioetYGpoPghAAwFL4OdPG/vxrnbjhO8TZ+6QSUe2C7AOCEADAsowL4k6O1ZRLFL5BTL2KpWGDQxACAFgcDwdaEc1/3J2f/Is0PV0qxMDuhoQgBACwUEMC2bHRgiNP7daJGy+ghaahIAgBACxXIy0t6c5/H8O/dEhOSJGul6tdkC1CEAIAWLropuzIKCHYndqv0608jaWhiSEIAQCsgJNAC6P4HYOFFaflodvFHAzsNh0EIQCA1YjwYgeGC72bcl02iUuOy5jYbRIIQgAAayJwNDecSx8mbLwo99kinsHA7npDEAIAWJ9WjVjqEGHiQ1wvDOyuN/MFYWVl5fz58zMyMsx2RAAAG8YxmhbKHRwhpFyVu2wSf74kD9gq5lfe8Zz8Suq/VbxYjJx8EPMF4TvvvLNu3brMzEyzHREAwOa1dGM7BgsvR3BP75HyKil2i3ir4vanCipp4Faxgwdr4cpUrdHSmSkIjx07duHChbi4OPMcDgDArowL4n4frWnpyrKLle6J4q0KKqikAVvFHr7sw2682tVZOnMEoSiK//73v999910zHAsAwD41daK1/fgV0fyVUiV0ra7vz0hBQwkm3NehQ4dmzJhRfUtERMQXX3yxePHiRx991NfX14THAgCAuyUEcV28uY4bdJl5SlwAu1ZGvk5q12TxDA3C4uLijIyM3NzcsWPHVt9+8uTJnTt3+vn5jRgxIioq6vDhw3d/7cWLF5OTk7/77rszZ86kpqb6+fkNGDDABLUDAMCdCirp0VRxahsuv5JWn5U/PyPHB3Jzw7n2HrhMeF9MUWrvJvr1119jY2ObNWuWk5Oj0/0zBf3nn3+eNGnSk08+eeTIkcrKyt27d/P8g5bh//rXvzp37jxp0qR7fvbChQsxMTHZ2dnGfg+2p6ioyM3NTe0qwDh41ayOKIo6nc7JyXZWTNWvCypEz++X0q4qI1twX/wht21M/2rPxTe3+lvmysrKtFrtg7PGWAb9T4mMjMzLy/vpp59qbF+wYMGiRYvef//97du35+bmJiUlPXg/w4cP79y58wOeIMtyXjXVQxcAAB6gRncMI/qwGx/jz37OkQ+NECY+xM09KEduFFeflUUMK72TQStCvczMzMjIyKpwun79epMmTW7cuOHl5UVEs2fPFkVx2bJldS7l5MmTHTt2rP439dy5c2fNmlXnHVqv4uJiV1dXtasA4+BVszo2tiIs0LHvsvlnWt/xrvYK0ad/CI8FSY00iqzQ9j+55WeF80VsSog05SHJXWN99xcauyJ0dHTUaDQPfk7dm2WuXr2q1Wr1KUhEfn5+Bw8erPPeiMjZ2TkwMBCnRvVwks0a4VWzLjYWhG5Ecz3vsX1utdNwCe6U0IYybigfHec7bpEff4j7d0cuwMWaLh8KgqDOqVEDGb64BAAAtUR6s9V9+d9GCkTUcYP4xC7pRJ5d//auexD6+flVVlbm5eXpP7x27Zqfn5+JqgIAgIYV5MaWdOezx2s6e7OB26TozeLmS3Z68bDuQdikSZPw8PAtW7YQkSRJ27dv79+/v+kKAwCABtdIS7M7cFnjhWmh3LyDcie77KYx6BphQUHB1KlTCwoKJElKSEjw9PRcvnw5Ec2fP3/atGl//PHH4cOHXVxc4uPjG7haAAAwPS1HT7TiJrbiUq4oS05Ir2fIM0K56W25xlq1KzMLg7pGKyoqEhMTqz50cnIaNmyY/nFGRkZycrK3t/cjjzzi7Oxcn1JwH2EV3JFmjfCqWR0ba5YxoSM3lQ8z5S05lthN0xD3ERpx+0RDQxBWwa9Ua4RXzeogCB/sQpGy/LT85Rl5UAD3P+FcB8uYTaPaDfUAAGBvWrqxhVH8ufGazt5s8N/dNJaycjIpBCEAANyXu4Zmd+DOjxemhXIvHZI7bRBXn5V1ttVNgyAEAIBa6LtpMscIix/m12bLLb7XLciQ8ivVLstEEIQAAGAQRhTXjG0eICQNFLIKKeQH3ex9Uk6J1Z8uRRACAIBxIrzY6r784ZGCk0CdN4oJKdKh61YchwhCAACoi6pump6+bHSyFXfTIAgBAKDuqnfTvHxIjtggrjwtl0tql2UMBCEAANRXVTfNsh785kty8A+6BRlSXoXaZRkGQQgAACYT3ZRtHiBsHSRkFVLIj7rZ+6RLxZZ+uhRBCAAAJhbuyVb35U+M0Xg4UORGMX6HeNCCu2kQhAAA0CD8nGlBJJ/9iCbOnxtrwd00CEIAAGhAbn9308zuwL15RA63vG4aBCEAADQ4DUfjgriDI4RPevCbL8lB3+sWZEi3LKObBkEIAADmo++m2T5Y+LOU2qzVzd4nXVS7mwZBCAAA5tbRk62I5jPHaDwcqPNGMX6HeCBXtThEEAIAgDqaOv3TTZOQqlo3DYIQAADUpO+mOZcgzO7AvXVEbrNWXHJcLhPNVwCCEAAA1KfvpjkwQviyF5989fZsGn03zZGbytQ9knTnUvFGOY1NkQp1Jjg0ghAAACyIvptmx9/dNNPTJQ2jv8qUR1Ml8e83BL5eTv2SxGA3cteY4IgIQgAAsDhhnmxFNH98jMbPmWKSRFmhq6XKhDRJlOl6OcUliQMD2HtdeZMcC0EIAAAWyteJFkTyWeM1/fy5S8WU+qfcb6cwYJtswhQkBCEAAFg4Nw3NCeOyxgsLIvnfb1GZxP43ymQpSAhCAACwCnmV9NlpeWYbuU2jO64X1h+CEAAALF3VdcGFkfKPMVyZdPt6oUkgCAEAwKLV6I5x4GldP8GEWYggBAAAi5ZfoUwI4ap3xzjwtLafEOxOFQhCW/XBBx+IohnHKkC9FRYWLl++XO0qwDjHjx//+eef1a4CateqEZsbfjutNm/efOLECSJy5GlhFO8imGD/CEJL9NFHHxUVFaldBRjhypUrX3zxhdpVgHEOHjyYlJSkdhVgnKSkpIMHD5p2nwhCAACwawhCAACwawhCAACwa0xRVH5r4Cpnz57t0KFDQECA2oWo7+LFi4GBgRyHP1Oshk6nu3btGv71Wpfi4uLy8nJvb2+1CwEj3Lhxw9HR0dXV1cDnT5gw4c0333zwcywoCIkoKytL7RIsQkVFhYODg9pVgHHwqlkdWZYlSdJoTPH+BWAuOp2O53nD1wl+fn5OTk4Pfo5lBSEAAICZ4eQbAADYNQQhAADYNQQhAADYNQQhAADYNVOMaQMTURTl0KFDycnJt27d6tChw4QJE7RardpFgaHOnj2blpY2bNgwf39/tWuB2omi+P333x8+fNjT03PMmDHt2rVTuyKohaIomzdv/vXXXx0dHYcPH965c2dT7RkrQguSk5OTkJCQn58fGBi4fPnyuLg4SZLULgoMIorixIkTZ8+efebMGbVrgdqJojhkyJBPPvkkICBAFMV9+/apXRHU7o033pgzZ05ISIizs3NMTExycrKp9owVoQXx9/c/d+6cIAhENHnyZB8fn+PHj4eHh6tdF9Ru4cKFAwcOvHDhgtqFgEG+/PLLv/76KyMjQ//jBlZh3bp1b7311iOPPEJEOTk569evj4uLM8mesSK0IIIgVP1YiqIoSZLh0xNARadPn167du1LL72kdiFgqKSkpIkTJ27duvXDDz/EctBatGvX7siRI0Sk0+mOHz/evn17U+0ZQWihZs+ePXr06JCQELULgVrIsjx16tRly5Y5OjqqXQsYKjs7+7PPPlu/fn1RUdGYMWOWLl2qdkVQuxUrVqSkpAQHB/v7+wcFBc2aNctUe8ZpAUv02muvHT16dNeuXWoXArX74IMPwsLCoqOj1S4EjMBxXJs2bVatWkVEUVFRjz/++LPPPqt2UVCL5557LjAwcNWqVYWFhdOmTfviiy+efvppk+wZI9Yszttvv71mzZq0tLQmTZqoXQvULjIysqyszMXFhYiOHTsWFBT073//21Q/n9BA4uPj27Vr9+677xLRpUuXWrRoUVhY6ObmpnZdcF/l5eXOzs6nT59u3bo1EX355ZfLly831Tv0YkVoWT744IPVq1fv2rULKWgtvvvuu9LSUv3jAQMGzJkzZ+jQoeqWBLUaOXLkN998oygKY2z//v2BgYFIQQvn4ODg6up6/vx5fRCeP3/ey8vLVDvHitCCnDlzJjQ0NCgoyNPTU7/l/fff79Onj7pVgeGaNm363//+NyYmRu1CoBZlZWX9+vVzcHAICQlJTEz87LPPRowYoXZRUIulS5cuWLBg9OjReXl5aWlpiYmJPXv2NMmeEYQWpKys7OTJk9W3hISENG7cWK16wFj6U6NYW1gFnU6XlpZWXFzcrVs3zECwFtnZ2RkZGY6Ojj169PDw8DDVbhGEAABg13D7BAAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2DUEIQAA2LX/D/qNQLN75LprAAAAAElFTkSuQmCC", "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" ], "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" ] }, "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.10.5" }, "kernelspec": { "name": "julia-1.10", "display_name": "Julia 1.10.5", "language": "julia" } }, "nbformat": 4 }