{ "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.844838171836 NaN 1.97e-01 4.0 \n", " 2 -7.850309780061 -5.47e-03 2.93e-02 1.0 \n", " 3 -7.850617563449 -3.08e-04 2.78e-03 1.5 \n", " 4 -7.850646835438 -2.93e-05 1.25e-03 2.5 \n", " 5 -7.850647360445 -5.25e-07 5.67e-04 1.3 \n", " 6 -7.850647505437 -1.45e-07 3.68e-05 1.0 \n", " 7 -7.850647511425 -5.99e-09 1.49e-05 2.5 \n" ] } ], "cell_type": "code", "source": [ "scfres = self_consistent_field(basis, tol=1e-8, callback=callback);" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "... and show the plot" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=1}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nO3deUBU9f4//tf7fc4wMLgAAooCioqJCypouaCooIKCO65p2YKm1a1uH+12u2V1vWXdX/682ablrczqikuKouK+a+5mZi64obIpiDAsc5bvH6OEhgo4M2eW5+Of4Mz2YhrPc97nvF/vw1RVJQAAAFfFtS4AAABASwhCAABwaQhCAABwaQhCAABwaQhCAABwaQhCAABwaQhCAABwaQhCAABwaQhCAABwaQhCAABwafcMwjNnzkRHR0dFRU2ZMkWSpPvcrbi4uPqvpyhKzQp0XngrKiiKgqX+zFRVxVthhreigqqq2F1UsMZbcc8gfOGFF954442dO3eWlpYuWbLkXnd7+umnDxw4UP3XMxqN+HCblZaW4sNtVlZWJsuy1lXYhfLy8vt873QpkiSVl5drXYVdkGW5rKxM6yrsgqIopaWlFn/aqoNQUZT9+/fHxsYS0YgRI9LT0y3+wgAAAPag6iC8fv16vXr1GGNE5Ovrm52dbduqAAAAbKTqIKxbt25JSYn55+Li4vr169uwJAAAANv5IwgXLFjQp0+fgQMHrlu3Tq/X+/n5ZWRkENGePXsiIyO1qxAAAMCKRPN/fvzxx3feeWfRokW5ubmjR4/euXPnrFmzRo0aFR0dvXXrVpwjBAAAZ8XMczh79OjxxBNPJCcnE9G0adOI6JNPPsnKyrpw4UJ4eLiHh8e9Ht+1a9fIyMgWLVpU3hgREdGjR4+KX7dl0f/O0SfdiBEVFRV5enoyxs4U0qv7aUUMMWv9afbOaDS6u7tzjlZOKikp0el0oihqXYj2ysrKOOc6nU7rQrRnMpkURdHr9VoXoj1Jkkwm0332w67DPGvUYDBU/yGiKJrnu9zvPub/HDly5NNPPzX/3LVr188//5yIGjVq1KhRoweWdeXKlbteJigoqPJs+Cg/2nyFTd7FPnlMkWVZluWzN9n4HXxhD0Vx4Tnz5rcCzSREJMsy5/yBH1ZXYP5I4OsREcmyrCgK+mro9r4CbwURmT8SNXorBEGoVhDevHnTaDR6e3ubN/n4+GRlZVXzNTw8PF566aXo6Oj73+1fj9GbB+UXDwgfdXS/VKafsEtZ3Fdo5+3SOz5FUTAiNFNVFSNCM8YYRoRmgiBgRGgmSZIgCO7u7loXoj1z77XF3wpORJ6enqIoGo1G8yYrTRN9J1Jo5EFP7RZHbVa+6+3qKQgAAHaCExHnPDAw8OzZs+ZNZ8+ebdq0qTVe7IlQvjtPuFisGvDVHwAA7MOtRBo7duznn38eHx9vNBq//fbbt99+2+KvdPqGOmqz/FPv8neO69suk5bHiHFBGBQCgDPYsGFDxVjC4swnxlz8gLmvr+/IkSOt9OQiY+z06dPTp08fPHhwSEhIaWlpfHz8iBEjLPsy5hRc1FtoplNT+wtPbVdGbJKeb8Pf6yJwpCEAOLgvvviisLCwefPmWhfinK5du3bmzBkrBmHFrMXt27dfvnzZzc3Nz8/Psq9x+oY6erP8XW+hrTcrKiIi+m+04KWn784ov+ari3qL3jgdDgAOLjk52Xp7ahd3+PDhp556ynrPf8eUxSZNmlg8BYnIU0ff9RHa3jk7Zk5X4ateYjsf9tgq6ZfraCEAAABt2GLWSmMDa1xV+2NcIIsLFDr4KH3WSP9+THiyFRoJAADA1rSfvjm2BW/rzYZvlPfkqPO6CzqkIQAA2JBdxE64Dzs0TMwtpb5p0lWj1tUAAIArsYsgJKJ6OloWKyQE8YgVpi1XccoQAABsxF6CkIgY0YwOfGEvccxmafZRRetyAADAJdhREJrFB7H9Q8Rl55XxW+RiSetqAACs6WSBWlD+gPvcKKffCnCczIrsLgiJKLgO254geojUfZV0thD/+wHAaWXcpMR06T5ZeKOcEtKljJs2rMn12GMQEpG7QF/2FKa14V1XSSvO4zApADingUHs7x2Fgeuk62VV3HqjnBLTpVfb80G1XZDyP//5T25u7kOV+NBOnjz5v//9r8qbzp49+/XXX9u2nCrYaRCaJbfma+PEV/Ypr+2XFYwMAcAZxQWyNyOEhPV3Z6E5Bf/ang9pWvsd9YIFC/Ly8h6mvDFjxhw6dOhhnuHs2bNr1qyp8qZz5879+OOPD/PkFmHXQUhEnX3Zz0PEA7lqQnrV35gAABzdn7PQIil4F6PRWFhYeNfGwsLC+18e/MSJEzdv1uzIrKqq165dk6RbszwGDRr07bffVtx68+bNsjL72pvbexASkZ87rYsTw31YpxXS/lwMDAHACVXOQounoKqq06dPb968eYcOHXr27Gm+9Pq+ffvCwsI6d+4cEBAwZ84cIkpJSRk9enTFo3x9ff/+97+fOXMmOTm5c+fOP/zwQ3Ve6/jx423bto2KigoICPjss8+IaMWKFUOHDj127FiHDh2ef/75Tp06+fv7v/feexb50yxC+5VlqkPk9H4XoZu/kpAuvddFeAqLsQGA3dtyVS2p4dT3/oG8R6pERKOaMx1naZdq8NXfQ6Q+AVWfSvzhhx+2bNmSkZHh4eExefLk6dOnL1iwYPTo0W+//fYTTzxx/vz5yMjILl26lJeXF5kvjEBERPn5+a+88kpqaurHH38cHR1dnRpUVU1KSnrzzTfHjh175syZyMjIPn36mEwmo9EoSdKxY8dmzJgxb968M2fOhIeHJycnV/+vsyrHCEKzIU15ay82fIO8L0f9uLvghjQEADv27Wklu6RmB7Ekha4UqyrR7mzanyvX6LH+7qxPgFDlTT/99NPTTz9tMBiI6MUXX+zdu/fx48dLSkomTpxIRM2aNRs5cuSaNWvatWtXzdcyGo2nTp0iIp1O17Zt24rtJ0+evH79+pgxY4ioZcuWcXFx69ata9y4sflWHx+fcePGmW9q0qRJRkZGjf5A63GkICSiR+qzvUPEp7bLPVZJS2OFpnVwMUMAsFP/7VV1LN2L+Yjot70FvcDeOSSvHiD6WOgSddnZ2d7e3uafvb29CwoKrl275u3tzditXWiDBg2uXbtW+SGqqt7n3OGFCxdeeOEFIvL19V2xYkXF9pycHC8vr4qn9fHxyc/PrwjCevXqVdxTr9eXlz+og9JWHG9UVVdHS2KEkSH8sZXSpis4ZQgAzqDyecF7zSOtNVVVf/nlF/PPx44da9GiRWho6MWLFwsKCswbjxw5Ehoa6u3tXRGHp0+fNgehTqeT5bvHpmFhYTt27NixY0flFCQiRVEuXbp048YN869Hjx5t2bKlZf4Ga3KwEaGZeTG2x/zZ+C3yi2359A64yj0AOLA/z46JC2REQsJ6yVLjwsWLF4eGhvr7+7/88ssvvfRSSEjIoEGDJkyY8NJLL+3atevgwYOLFi0iot9+++3jjz8OCQmZP3++eWDXvn37zz77LCsrq0uXLqGhoQ98ITc3t2eeeWbatGmbNm26fPnyyJEjU1NTLfAHWJPjjQgr9A5g+4YIKy4owzbIhSatqwEAqJV7zRG17Ljw7bffPnHixOLFi//xj38899xzRPT999/HxMQsXLiwuLj4wIEDDRo0aNCgwbp16w4fPrxy5crZs2fPmDHDw8Nj7ty5ffr0+eWXX7Kzs6vzQs2aNXvmmWcWLFhw8+bNXbt2eXh4hIWFjRkzpmHDhlOmTKm429NPPx0YGGiBP8wSHHJEWCHQk20bJL6wR370J2l5P6GNF0aGAOBI7t8pYcFxYf369WfPnl15i06ne+mll+66W7du3bp162b+uaLDYerUqTV6rQEDBgwYMKDi1/bt27dv356IZsyYUbHx5ZdfJqLTp0/X6JmtxIFHhGZ6geZHCa935H3WSMvOYTE2AHAkS88p/xd+v37BuED2ekdhOVaatCbHHhFWmBjK23mzkZvk/XnqrM6CgJEhADiCpx958GgkIZgR1XKn9tFHHwUGBm7fvr12D6+pPn36HDlypPr3Dw8Pf+ONN6xXTzU5/IiwQoQv2z9UPJSnxqZJ2SVaVwMAYAf69etXt25drau4J39//6ioKK2rcKIgJKIGelobJ0YHsEdXSvty0FkBAAAP5lRBSEQCo5kRwrzufMgGacFJHFUHAIAHcLYgNEsM5tsTxP/8qkzcKtd0rT8AAHApzhmERNSqPtszWCxTKGq1dP4mDpMCAEDVnDYIiaiOjv7XV5jcmndbJW24jCwEAIAqOEn7xH0kt+atvdi4LfIzj7A3OwlYjQ0ArCEzM/PXX3/VugrnZO3rVDh/EBJRr0Zs32AhaZN8KE/+trfg5aZ1QQDgXCZMmPD6668vWLDASs+vqmrFJR1cU1hYmPWe3CWCkIiaeLLtCeIbB+XHVkrLYoV23i79kQIAyxoyZMiQIUOs9OSSJJlMJg8PDys9PzjzOcK7mC9z/49OPCZNWpKBzgoAACBynRFhhcdb8nbebMRGedMVdV53QedC3wQAAKAKrpgDHRuw/UPFi0VqbJqUhcXYAABcmysGIRH56CktThwYxB/9SdqDxdgAAFyYiwYh3b7M/WdRwpB0afZRnDIEAHBRrhuEZoOC2O7B4uKzyoStshGLsQEAuB5XD0IialmP/TxEFDn1SJUysBgbAICLQRASEbkL9N9ewnNhvMcqaV0mshAAwIUgCP+Q3JqnxIjP7JBf2y8rSEMAANeAILxDVCN2eJj4c446OF3KL9O6GgAAsD4E4d383Ck9Xmznwx5dKf1yHQNDAAAnhyCsgnkxtncieZ810ten0FkBAODMXG6Jteob24J3asCGbZT35Kgfdxfc8J0BAMAZYe9+P6292L4hYm4pxaRJV41aVwMAAFaAIHyAejpaFiskBPGIFaYtV3HKEADA2SAIH8y8GNvCXuKYzViMDQDA2SAIqys+iB0YKi47r4zbIhdjMTYAAGeBIKyBIE+2PUE0iNT5J+lkAQ6TAgA4AwRhzbgL9GVP4eV2vOdqacV5HCYFAHB4aJ+ojeTWPNKXjdwk78tV/9VZ4EzrggAAoLYwIqylSF+2f4h4ME9NSJeuYzE2AACHhSCsPV93WhcnhvuwTiuk/bk4ZQgA4JAQhA9FYPR+F+E/3XhiuvTV7zhlCADgeHCO0AKGNOWtvdjwDfLPuViMDQDAwWCfbRmP1Gd7h4jXy6j7KulCEQ6TAgA4DAShxdTV0ZIYYUJL3m2VtPEyshAAwDEgCC2JEf2lHf+hj/jENnn2UVzlHgDAASAILS86gO0bIqy4oAzbIN8o17oaAAC4LwShVQR6sm2DxIYe9NhK6QQWYwMAsGMIQmvRC/RFlPB6R95njbT0HDorAADsFNonrGtiKG/vw0ZslHdkqf/fY4KILx4AAHYGO2ar69SA7R8qnixQY9dKWSVUUtUlnIy4rhMAgEYQhLbQQE9pcWLvANZphRSVKt003XHrwTw1Nk2ScPQUAEALCEIbERjNjBDmR/GMm2rnn6TC21l4+Jr67A55cR8cNQUA0Ab2vjaVGMx/HiKWymqrJaacUnb4mvr0dnlZrBBSF1dyAgDQBibL2FpoffbrCN2AdVLbVUKzOkpanNi0DlIQAEAzGBFqoI6O5nYTPATKKKQi04PvDwAA1oMg1MDha2ryDnn7AGliK9b5J+nwNXTcAwBoBodGba3ivGBDgT7rzkWudl0lbUsQu/rhACkAgAbuOSLMzMycP3/+a6+9dvnyZVsW5Nx+K1Cf2SGv6PfH7JhPugvTwnjvVOlwHsaFAAAauGcQnjx5Mj8/f82aNdnZ2bYsyLk1r8tW9RPumh3zUVfh31153Hrp51xkIQCArd0zCGNjY2fMmOHv72/LapyeXqAmnlUcAn2+jbCwlzg4XdqTgywEALApTJaxF4OC2NfR4uB0afMVZCEAgO2IRPTpp59+9dVXlbc+++yzU6ZM0agk1xUXyJbGiKM3S4t6i/2aYO4MAIAtiEQ0derUqVOnal0JEBFFB7ClMeLwjdKXPYXBTTFeBwCwOv63v/2tyhuuXr2akpKSk5OzcePGbdu22bgsVxbViK2NEyfvlFecxzrcAABWx+fOnZuXl/fnG4xGY0ZGxuOPP66q6pUrV+71eEVRcnNzM+9048YNa9bs/CJ9WVqcOHWX/N0ZZCEAgHWx9u3bv/XWWyNGjKjd4zt27JiZmanX6ytvfOqpp6ZPn17l/YuLiw0GA2M4AUYlJSV6vZ7zex7//L2QDd2m/0d707hmsi0Ls72SkhKdTieKWN6BysrKOOc6nU7rQrRnMpkURblr3+KaJEkymUweHh5aF6I9RVHKyspq9FYYDAZBEO5/HzEgIOBhWubr168/d+7c6Ojoat6fMebp6YkgJCJBENzd3e8ThJ3r0qZBamwakwU+tY0zny8URRFBaObm5oYgNEMQVkAQVlAURafTGQwGyz6tqKoqYsluPVKf7UgQYtJkk0J/aefMWQgAoBUxOzs7ICBA6zLgnprVZVsGCTFpsqzSK+2RhQAAFsZPnz5d/QOboIngOmxHovjV78rbhzB3BgDAwvjzzz/v5+endRnwAI08aPMgcdl55bX9Tj5xBgDAxvgHH3ygdQ1QLQ09aOsgcfMV9f/2IQsBACwG55wciY+e0uPFHdnqc7tkBSuSAgBYAoLQwXi50fo48dh1dQqyEADAEhCEjqe+G6XHixmF6oStsoTZMwAADwdB6JA8RVrVX8wrVcdvlU3IQgCAh4AgdFQGkVIHiOUyDd8olWH2DABAbSEIHZgbp//FCHqBDdsolUhaVwMA4JgQhI7NjdP/+gp+7ix+vVRk0roaAAAHhCB0eAKjhb2E5nXZwPXSTWQhAEANIQidgcDoq15Cxwas7xrpepnW1QAAOBQEoZNgRHO7CVGNWL+1Ul6p1tUAADgOBKHzYERzugqDglj0aumqUetqAAAcBILQ2bwTKUwM5X3TpCtGLDwDAPBgCEInNKMDn9SK91kjXypGFgIAPACC0DlND+d/act7rZbPFiILAQDuR9S6ALCWqW24yCkmTd4QL4TWZ1qXAwBgpxCEziy5NfcUqU+avC5OaOeNLAQAqAKC0MmNb8kFRgPWymlxQgcfZCEAwN0QhM5vTAsucuq/VkrtLz7qhywEALgDJsu4hJEhfGEvcXC6tCcHc2cAAO6AIHQVg4LYN9Hi4HRp8xVkIQDAHxCELmRAIFsWK47bIm24jCwEALgFQehaejViy2LF8VukVRdwYXsAACJMlnFBPRqytXFiYrokqTS8Gb4JAYCrQxC6okhfljZAjF8nGSV6vCWyEABcGoLQRXVswDYPEvuvlSWFnmyFLAQA14UgdF1hXmzjQKFfmlws0bQ2yEIAcFEIQpf2SH22PUGIXStLCv2lHbIQAFwR9n2urlldtnmgMO+EMusI5pECgCtCEAIF12E7EsUfzipvH0IWAoDLQRACEVEjD9o8UFx2Xnltv6x1LQAANoUghFv8PWjbIHHzFfXVfTIWngEA14EghD946yk9XtyVrU7dJSsIQwBwDQhCuIOXG22IF38vUCfvRBYCgEtAEMLd6uho9QDxfJH6+FZZwuwZAHB2CEKogkGk1f3FIhON3yqbkIUA4NQQhFA1vUBLY4VymYZvlEoxkxQAnBeCEO7JjdOSGMFdYMM2SCWS1tUAAFgHghDuR8fpx76CvweLXy8VmbSuBgDAChCE8AACo4W9hBZ1Wfx6qRBZCABOB0EIDyYw+rKX0KkBi1kjXS/TuhoAAItCEEK1MKK53YSejVhsmpRXqnU1AACWgyCE6mJEH3UVEoNZ9GrpqlHragAALARBCDXzdqQwMZT3TZMuF2PhGQBwBghCqLEZHfhTrXjP1fK5m8hCAHB4uEI91Mb/hfM6OuqbJm+MF1rUY1qXAwBQewhCqKXnwrjAKHqNnB4vtPFCFgKAo0IQQu0lt+Z1dNR/rbx2gNDeB1kIAA4JQQgPZVwLLjCKWyenxQkdkIUA4IAwWQYe1ujm/PMoHrdW2peDuTMA4HgQhGABicH8q17ikA3S7mxkIQA4GAQhWMbAIPZjX3H4RmnzFWQhADgSBCFYTO8AtiRGHLdF2nAZWQgADgNBCJbUqxFbFiuO3yKtvIAL2wOAY8CsUbCwHg3Z2jgxMV2SVRreDN+0AMDeIQjB8iJ9WdoAMX6dVCzRhJbIQgCwawhCsIqODdjmQeKAtbKs0JOtkIUAYL8QhGAtYV5s40AhNk0ukuj5NshCALBTCEKwolb12fYEIXatLCn0UjtkIQDYI+ybwLqa1WWbBwqfnFD+eRjzSAHAHiEIweqC67AdieKPGcpr+2WtawEAuBuCEGyhkQdtHiiuvaQiCwHA3iAIwUb8PWjrIHHLFfX53TIWngEA+4EgBNvx1tP6ePFgnvrcTllBGAKAfUAQgk15udHGgeKpG2oyshAA7AOCEGzNU6TVA8TfCtTWS6XyO88YFplo3Ba5RNKoMgBwSQhC0IBBpM0DRc6odYpUdjsLjRIN2yglBDMPdLcCgA0hCEEbeoGODRf1IoUuMRklMko0ZIM0qRUf1wKfSQCwKex0QDNunI4NF+u6sfBUYcgmFSkIAJrAfge0pOO0K0E0ymxnNqVeVJefV3CCEABs7J5nYyRJ2rdvX25ubkRERHBwsC1rAtdhlChpszSns/xrofBzrrrwlDppm9wrgCWF8GHNeF2d1vUBgAu454jwxRdf/Prrr48dO9avX7+0tDRb1gQuouK84Ohm6nuRrE9j1tCdnR2tSwrhKeeUwO9NienSt6eVQpPWhQKAU7vniHDevHmccyIKDw9ftGjRwIEDbVgVOL/Ks2NKSoiIZkYIMw/JM36WF/QUJoby/DJKvaiknFNe2H1rjDikKa/vpnXdAOB07hmE5hQkon379rVv395W9YCrKDTRi215YvAdxyRmRgiLzigmhfQCeetpYiifGMoLymnVBWX1RfXlvabuDVliMB/WjPu5a1U4ADgbkYjS09O///77yltjYmImTJhARKtXr966deuWLVu0qQ6cVyMPuisFzSa0vHujl5s5EckoCZuuKCkZ6vSfTe28WVIIH9uC+3vYpFwAcF4iEUVGRgYEBFTe2qBBAyJKT09/991309LS3N3x9Ru0ZxApMZgnBlOJJGy8oqRkqG8dupWIY1rwhkhEAKgVMSUlJSkpyZx8le3YseNvf/vbmjVr/nwTgLY8/pSIMw+Z2nqzpBA+ugVvhEQEgJrgzzzzTFFR0Z9vWLp0KWMsISGhc+fOU6ZMsX1lAA9kTsRvewtXx+tmdOAH89SwFFNUqjT3uHLVqHVxAOAgWHh4+GuvvTZ27NjaPT4iIsJoNNapU6fyxjFjxtwrO4uLiw0GA2Osdi/nTIxGo7u7e8WkJFdWUlKi0+lE0QJrjJbKtCWbr7gorL3CWtdThwUpw4LlAMcZI5aVlXHOdTp0UJLJZFIURa/Xa12I9iRJMplMHh6O8zm2GkVRSktLDQZD9R9iMBgeuJsVQ0NDz58/X+uyDAZDcnJyp06dKm8MCAi4Kxor8/T0RBASEeccQWgmCIKlgrAOUVJ9SmpFZTKlX1ZWX1R7pist67HEYJ7UnDWva+8fPJ1OhyA0QxBWQBBWUBRFFMUaBWF1iAaDobCwsNaPFwQhLCzsscces2BNAA9PL9w6j/hpD2FPtppyTum2SvbVs6Tm7PGWvGU9e09EALAZnp+f37BhQ63LALAWgVFUIza3m3BlnO6LKCG/jKJSpbZLpZmH5NM3cGlgACC+f//+Dh06aF0GgNVVJOLl24nYa/WtRPwdiQjgwkQ/P7/evXtrXQaA7ZgTMaqRMKersDtbTTmn9FkjebuxpOZsTHPe2gtHTQFci7hmzRpMXQHXxP+UiDFpspcbJTVno5rzNkhEANcg4hJLAJUSkXZnq6svKYnrZXeBEpuyhCAe1QiJCODMLDBhHcBpVCTi+13o13w15ZwyabusqJQQzJJCkIgAzglBCFC1tt6srbcwM+JWIj61QzYpNDiYJYXwHo1wOgHAeSAIAR7grkR8ZodcKtOQpkhEACeBIASorrsS8dmdslGioUhEAAeHIASosbsScfJO+aaJhjVjSSG8e0PGEYkADgVBCFB7lRNx9UX1tf3yhSIaGMQSgll8IBexjiyAI0AQAlhAW2/W1pvN6MAzbqqpF9TZR5WntsvxgTypORIRwN7hHyiAJTWvy/7Sju9MFH8eIkb6stlHlYDvTRO3yqkXFZOidXEAUBUEIYBVhNxOxANDbyfi4qoT8aaJ0i9Xsdhp2iW1VLZRtQCuDEEIYF1N69xKxEPD7k7EcoWISGD04TF58Zk74nHhKWXeCcQggC3gHCGAjQTXYX9px/7Sjl8qVpefU2cfVZ7cJg8K4knNWUqMmLRJIqKRQUREC08pS88py2NFd0HjmgFcAYIQwNaCPG8lYmaxuuyc+p9flUnb5OgA/t5RuaiMOKcVF5GCALaDIATQTODtRLxiVJefV3NK2fP7eH03+qQ7V3CFRABbQRACaK+xgT3fhhlEYiRll7B/HpYn75T7NeHDm7FBwbyeTuv6AJwaghDALpjPC6b2VYmxkVvZC224u0grzqvP7zZF+LKEID6mBW/ooXWVAM4Is0YBtFd5doxBpJX9xJTzisBoSYxwaazuxbb8YJ7aOsUUlSrNPa5cMeKwKYAlYUQIoLGrRkq7pK6IFfUClUlERAaRVsSK47fKg4LJy40Sg3liMJXKwobLSkqG+vZhuY0XSwzmI0JYy3pY2BTgYSEIATQWYKClMXfPEK2jo5X97tjoLtxKRFkV9mSrKeeUqFS5gZ4lNWdjmvPWXkhEgFpCEAI4GIFRVCMW1Uj4qOutRIxJk73cKKk5Swzmkb5IRICaQRACOKqKRJzTlQ5fU1MvKmO3yJJCicG4RCJADSAIARweZxTpyyJ9/7hEYvJOuViioU1ZYjDvHcBw+QuA+0AQAjiVuy6ROPOQfOqGGhfIk5qzuECuQyIC/AmCEMA5VVwi8fxNdeUFdfZRZdI2eWAQT2rO+jfheqzfBvDCNFUAABgESURBVHAbvh8COLlmty8IdXi4GOnL/vOr4v+dKTFd+va0UmTSujgAO4ARIYCrqFjsO6+U0i4pKeeUF3bLvQJYUggf0pTXd9O6PgCNIAgBXI6vO00M5RNDeX4ZpV68lYiP+rOEID62BffHQm7gYhCEAK7LW38rEY0SbbqipGSoMw+Z2nqzpBCe1Jw1NqD/AlwCghAAyCDeWramRBI2XsFCbuBaEIQA8AcPEQu5gctBEAJAFbCQG7gOBCEA3A8WcgOnhyAEgGr580Juz+6UjRINbcqSQnj3howjEsExIQgBoMbuWsjttf1YyA0cGIIQAGqvYiG3czfVVXcu5DYgkLshEcER4HMKABYQcnsht0PDbi3k5rcIC7mBY8CIEAAsKbjOPRdyG9qM19NpXR/AnyAIAcAqKhZyu15Gq7GQG9gxHBoFAOvy0dPEUJ7aX7w6XvdiW34wT30kxRSVKs09rlw13nHPG+U0+6jy52f4/48r2SU2qhZcEIIQAGzEvJDbt72FK+N0Mzrwg3lq22WmqFRp9lHlbKFKRPXdqKBcfXmvrFZ61MxD8rHrqp+7VlWD88OhUQCwtYqF3MpkYUeWmnpR6XF7IbcnW/GvTymv7JVnRxIRzTwkXyyiL3sKaFIE60EQAoBm9ALFNmGxTf5YyK3vGtnLjXz09MR2alGPXSlBCoLVIQgBQHuVF3Lbla0uP6989ptqEFnOeKQgWB3OEQKAHeGMejZi9d1oVDOqp1Mjf5LUBz8I4KEgCAHAvpjPCy7oQfsT1CyjGr0aWQjWhSAEADtSeXZMPR0dHi4eu67Gr0MWghUhCAHAXhSUU5l8x+yYxga2e7C4M1tdcg5RCNaCIAQAe+HlRu91uXt2TBsvljZAfGGXdOQashCsAkEIAPauVyP2SQ8hMV2+UIQsBMtDEAKAA0gK4X/rwAeuk6+XaV0KOB0EIQA4hqlteHwQG7pBKpW1LgWcC4IQABzGh48Jzeqw0ZtlGYdIwXIQhADgMBjRl72EEkl9fjdGhWAxCEIAcCRunJbGintz1A+PVXHBJoBaQBACgIOpp6M1A4RPTijfnkYWggUgCAHA8TQ2sLQ4YfrP8obLOFsIDwtBCAAOqY0XWxIjjt+CRnt4WAhCAHBUvRqxL3sKg9FoDw8HQQgADmxwU/63jjwejfbwEBCEAODYngvjg4LYwPWSUdK6FHBMCEIAcHgfPCa0qsfGbkGjPdQGghAAHB4j+gqN9lBbCEIAcAa62432H6DRHmoIQQgATsLcaP/Zb8o3aLSHmkAQAoDzaGxgawYIM36W09FoD9WGIAQAp9LGi6XEiI9vkQ6j0R6qB0EIAM6mZyP2ZU9hSLp8/iayEB4MQQgATmhwU/56R95vrZxbqnUpYPcQhADgnKaE8WHNWGI6Gu3hARCEAOC0Zj8qPFKfjcEV7eG+7hmEr7/+emRkZERERJ8+fbKysmxZEwCARTCiL3sKpbI6bRca7eGe7hmEb7311sGDBw8dOtS7d+9PP/3UljUBAFiKudF+X646+yiaC6Fq9wxCvV5PRAUFBWfOnAkNDbVhSQAAlmRutP/8JBrtoWoiER05ciQzM7Py1g4dOgQFBc2aNWv58uWMsTlz5mhUHgCABTQ2sLQBQp81UiMPNiCQaV0O2BdOROfOnTtwp9zcXCL6+9//fvDgweTk5BkzZmhdJwDAQwnzYiv7iU9sQ6M93E28efPmsGHDhg0bdq97+Pv7FxUV2bImAABreMyfzY8SBqfLOxKEZnUxLoRbxKlTpy5atOjPNwwcODAoKEiW5Z07d3733Xe2rwwAwOIGN+VXjNRvrbx7sOjnrnU1YB+Ym5tbZmamn5/fXTcYjcbjx4+LohgWFubh4XGvx3ft2jU8PDwkJKTyxi5dukRHR1d5/6KiIk9PT8bwXYyMRqO7uzvnaOWkkpISnU4niqLWhWivrKyMc67T6bQuRHsmk0lRFPOsPYt7/SBty6L0AeTpCB86SZJMJtN99sOuQ1GU0tJSg8FQ/YeIovjAxBH9/PxOnDjx59wyGAyPPvpodcq6ceNGXl5e5Y1FRUWKUvXsLEVRFEVBENLtt0LrKuyCcpvWhWjP/CbgrSArfyre7UTP7mLjt1FKb1Ww+70R/oFUsNJbIXp7e+fk5NT68R4eHlOnTr3X+O/PTCaTXq9HEBKRLMt6vR4jQiJSFAUjwgoYEZpxzq03IiSihb0pYb30ygH2eZRgpZewFEEQOOfWeysciKIoqqpa/K3gRUVF9evXt+yTAgDYOXOj/c+56ntotHd5/OrVq82aNdO6DAAAW6urozUDxPknla9PIQtdGm/fvn2rVq20LgMAQAMBBkobILy2X16fieZC18W/+OILrWsAANBMRaP9oTxkoYviERERWtcAAKClx/zZgp5CQrp0phBZ6IowZREAgBKD+VsRQvw6OadE61LA5hCEAABERJNb8xEhLDFdKsYV7V0MghAA4Jb3ughhXmzMZknCNFJXgiAEALiFES3oKZgUmrYbV7R3IQhCAIA/6DilxIj7c9V/HcGo0FUgCAEA7mButP/ydzTauwoEIQDA3QIMlB4v/P2AgkZ7V4AgBACoQst6bEU/AY32rgBBCABQtUf92IKewqD10ukbyEJnhiAEALinxGA+M1IYuB6N9s4MQQgAcD+TW/ORaLR3aghCAIAH+FcXoY03Gu2dFoIQAOABGNH8KMGk0FQ02jsjBCEAwIOZG+0P5qmz0GjvdBCEAADVUldH6+LEb04r/0WjvXNBEAIAVJefO6UNEN44oKxDo70TQRACANSAudF+4lZpbw6y0EkgCAEAauZRP7awlzhsAxrtnQSCEACgxhKC2duRQjwa7Z0CghAAoDaSW/PRzVkCGu0dH4IQAKCW/tlZaOvNRm9Co71jQxACANSSudFeUtFo79gQhAAAtVfRaP/PwxgVOioEIQDAQzE32n97RlmIRnvHJGpdAACAwzM32kevlgM8WHwQ07ocqBmMCAEALKDiivZotHc4CEIAAMtAo72DQhACAFhMRaN9NhrtHQeCEADAkpJb8zHNcUV7R4IgBACwsHc7C+3QaO84EIQAABbGiL6IEmSVntuFRnsHgCAEALA8HadlseLxfPVdNNrbPQQhAIBVGERa1V/8Do32dg8N9QAA1mJutO+5WmqgpyFNMfCwU/gfAwBgRS3qsZX9xGd3yGi0t1sIQgAA6+rix/4bLQ7bIJ1Co71dQhACAFjdoCD2TqQwEI32dglBCABgC8+25mNbsIT1aLS3OwhCAAAbeSdSCPdBo73dQRACANgII/ocjfb2B0EIAGA7FY32bx/CqNBeIAgBAGzK3Gi/+Kzy6QlkoV1AEAIA2JqfO60dIPzziLzyArJQewhCAAANVDTa70GjvdYQhAAA2jA32g9Ho73WEIQAAJoZFMTe7YxGe40hCAEAtPTMI3xcC5awXioyaV2Kq0IQAgBo7G1zo/1mNNprA0EIAKAxRjS/p6DjbAoa7bWAIAQA0J7A6Ps+wq/56sxDyEJbQxACANgFg0ip/cUfzqpotLcxBCEAgL3wdae0AcKsI8pPaLS3IQQhAIAdaVGPrewvJKPR3oYQhAAA9qWzL5vfUxy4ropG+7WX1Bx0HFoaghAAwO4MDmZtvFjXVXc02i8/r8w+JruL2pXlpBCEAAB2hzPakSgGe1LET7ca7ZefV/7zq7Kqv1hPp3VxTgdBCABgjzijg8NET4G6rFaXX2RIQetBEAIA2CmB0fERYn4JJe8RngjlbthhWwfeVwAA+7X6khLmTeFe6qwjSpPvTVN2yruzMZvUwnDWFQDATpnPC66MZXpVmrZf8PNggQY2dZecX05jW7Dk1rx5XaZ1jc4AI0IAAHtUeXYMZ7Swl5BbQjdMdGS4uDxWKJGo2yopKlWaf1LBZSseEoIQAMDulEj041k1tdLsGHMWZharZwrVSF82t5twcYxuRge+8bLa5HvTqE3yxssqjpnWDg6NAgDYHQ+RlsQId23kjOZH/bFRL1BiME8MputlwtJzyoz9cs52Gt+SPfMIb1kPh0xrACNCAADH5qOn5Nb84FBxXZxARFGpUuefpLnHletlWlfmIBCEAABOoq03e7+LcHmc7v0uwq5sNeRH06hNcupFBdf7vT8cGgUAcCoCo9gmLLaJkF8mpJxTZh9VntuljGjGnn6Eh/vgkGkVMCIEAHBO3npKbs13Joob4gVvPQ1Ol82HTPNKta7MziAIAQCcXJgXmxkhZIwW3+8iHMxTW6WYEtOllHOKCYdMiQiHRgEAXAS/fcj0Rrmw8oIy/6Ty3E55RAifEsY7NXDpQ6YPGBF+/fXXH330kW1KAQAAG6jvRhND+YZ4cfdgMcBAIzbKbZdKs48qLnulw/sF4S+//DJv3rzVq1fbrBoAALCZVvXZzAjhzCjxiygh46baKsXUb62Uck4pd7FDpvcMQkmSpk+f/u6779qyGgAAsDHOKKoR+yJKuDpOl9yazz+pBCw2Td4p78xylZVq7nmO8N///vf48eP9/f1tWQ0AAGjFQ6SkEJ4Uwi8Vq9+fUSdtl3WcngjlT7TijTy0Ls6aRCIaMmRIZmZm5a0ffvjh3r17f/rpp4MHD2pUGAAAaCPIk83owGZ04Afz1PknlbAUU2c/NqElHxnCDc44w1IkopUrV961dc2aNZmZmZ07dzYajZmZmVOmTPn888+1KA8AADQT6cu+iBLmdhNSLyrfnlZe2SsPDOITQ3lME+ZM00z57t27/7x10KBBBw4cOHDgwKJFizp37owUBABwWe4CJYXw1P7i0eFipC/76z45LEWaeUg+d9NJTiIKWVlZ48ePv9fNnHMfH5/27dvf6w4LFy5s3bq1p6dnfiWqqnp4VH1Euby83M3NjTnVl4laMplMoijirSAiSZIEQeAcyzuQLMuMMUG4+7IDLkhRFFVVRdEZj8TVkKIoiqLodLoH39XK6rmxrv7suTDerSHbl6O+tFdeeUGVVHqkPnOzyWdWVVVJkiz+VjA3N7fS0tJa7447deqUlZWl1+srb5w0adJf//rXKu9fXFxsMBiw9ycio9Ho7u6OvT8RlZSU6HQ67PKIqKysjHNuD7s8zZlMJkVR7tq3uCZJkkwm071GFxoqk2lzNv/hnLA1m8U1Vsc2k3s3Uqy6c1cUpbS01GAwVP8hBoPhgbtZRkQ5OTl+fn61Kys6Ovqdd96Jjo6u5v2Lioo8PT0RhIQgrARBWAFBWAFBWMFug7DCVSMtyVC+Oa3kldK4luzZR3gL61wQsRZBWB2ciBBLAABQawEG+ks7fmiYuDZOIKLuqVLnn6T5J5WbJq0rqx6u1+sbNGigdRkAAODwzBdEvDJO934XYeNlNfD7WxdElO17Vg3v378/RoQAAGAp5gsiLokRzo/RxTZhs48qTX+UXtsvn7php3nI33jjDa1rAAAAJ1RxQcT0eIGIoldL9nlBRP7oo49qXQMAADizNl7s/S5C5jjdXRdElOxjdW/M0wMAAFsQbl8QsaBcWHVBmX9SeWmPMjKETWrFO2p6QUQEIQAA2JSXG00M5RND+e831B/OKsM2ygaBJobySa24vxZNImhiAwAAbTxSn82MEM7eviBi2NJbh0xNdx4y3Z9bxSyb32+ohRZqz0AQAgCAliouiHhujC4phM8/qTSqdEFESaF3DyvzT96Rjb/mq2M3y7kllpmGikOjAABgF+rpbh0yNV8Q8cntsp7TxFD+SQ/hhd2ypChTWhMRnSxQH98qf99HsNT6NQhCAACwL3ddEDF8mSnCj311Si2TWB9/NmmP/H0fIczLYvNrcGgUAADslPmCiJnjdE+05PV0NOOAMnSrsKSvJVOQEIQAAGDnPEWaGMrndRda1qWmnrTlqoVXqMGhUQAAsHfm84IpMbyJW9mTewRJUaa2sdhADiNCAACwaycL1LFbbp0XdOO0JEZIv6x+esJiy9IgCAEAwH5JCk3ZJVc+L+jG6X99hXWZ6uFraJ8AAABnJ3LaNFAU7pwcoxdoRT9BsNCMGVuPCC9dulRUVGTjF7VPly9fLiws1LoKu5CVlZWfn691FXYhNzc3Ly9P6yrsQl5eXm5urtZV2IX8/PysrCytq9BSReAVFhZevnz5ro0Pz9ZB+OKLL+7cudPGL2qfXn/99bVr12pdhV2YNWtWSkqK1lXYhblz5y5cuFDrKuzCwoUL586dq3UVdiElJWXWrFlaV2EX1q5d+/rrr1v8aXGOEAAAXBqCEAAAXNrDTpaRJOno0aPVv39+fv4vv/xiMBge8nWdQG5u7okTJ7Zt26Z1IdrLyso6ffo03goiunTpUkFBAd4KIjp//vzNmzfxVhDR6dOns7Ky8FYQ0YkTJ3Jzc2v0VnTu3NnT0/P+92Gq+lDTT3/44Yd58+aJYnUDtaioyN3dvfr3d2LFxcVubm46nU7rQrRXUlIiCIKbm5vWhWivtLSUMabX67UuRHtlZWWqqrq7u2tdiPbKy8tlWfbw0OJKfXbGZDKVl5c/MNgq++qrr1q2bHn/+zxsEAIAADg0nCMEAACXhiAEAACXhiAEAACXhiAEAACXZtPZm/n5+RkZGYGBgQ0bNrTl69obRVEOHDhw7ty5gICAqKgozl3360hhYeGBAwdycnIaNGgQHR2NiaNElJGRUVxc3L59e60L0UxeXt6FCxcqfg0LC3Pxhqvffvvt6NGjXl5ePXr0qFu3rtblaOP8+fPXrl2r+FUQhI4dO1rqyW0XhIMHD16/fr0oiv/85z9ffvllm72uHerevbvRaGzXrt2xY8fc3d23bNnish/uadOmXb16tXHjxr///vu1a9d27NgREBCgdVFaunTpUpcuXerUqVM5CVzN8uXL//GPf4SHh5t//eyzzx44/d2Jvfrqq4sXL+7Vq9eNGzcOHDjwxhtvaF2RNpYuXbp+/Xrzz2fPnvXw8Pj1118t9uyqrZw6daqsrCw+Pv6jjz6y2YvapxMnTph/MJlMbdu2/fjjj7Wtxx4oitKzZ8/3339f60I0Fh8f/8orrwQHB2tdiJa++OKLkSNHal2FXVixYkVgYGBOTo7WhdiXXr16zZ4924JPaLuDcqGhoTjwZRYWFmb+QRTFxo0bl5aWaluPnZAkydvbW+sqtPT111/7+/v3799f60K0V1hYuHnz5uPHjyuKxS6+6ogWL16cnJxcXl6+Y8eOgoICrcuxCxkZGXv27JkwYYIFn9N1z07Zgz179uzfvz8pKUnrQrS0evXq4cOHt2nTpmPHjpMmTdK6HM1kZWX961//+vDDD7UuxC5cunRp9uzZcXFxPXr0qHxmyNWcPXt2+/btCQkJH3zwQWho6IYNG7SuSHtffvnlwIEDLXsOBUGomTNnzowaNerTTz9t2rSp1rVoKTw8/Nlnn33yySeXL1++b98+rcvRzLRp02bOnOnn56d1IdqbNGnSiRMn1q9ff/bs2Tp16rz11ltaV6QZo9FoNBoPHDiQmpo6a9asadOmaV2RxiRJ+uabb5566inLPi3W/NTGhQsXYmNj33zzzbFjx2pdi8aCg4ODg4Pj4+OLiormzJkTFRWldUUaOHny5Lp163x9fbdt25aZmXn9+vXJkye/9957Pj4+WpemgYoFePV6/ejRo7/66itt69FQQEBAt27dBEEgopiYmMmTJ5eUlLjyoqPr1q2TZTk+Pt6yT4sg1EBmZmZMTMyrr7767LPPal2LHSkoKKjRWrrOxN/ff86cOeafPTw89u7dGxkZiXPqRHT06NEmTZpoXYVmevfuffLkSfPPGRkZPj4+rpyCRLRw4cInnnjC4tcqsN2i2ykpKQcPHkxJSWnWrFmXLl3GjBljwS4QxxIREVFYWDhy5Ejzr127dh06dKi2JWmlf//+3bt39/X1PXLkyJIlS7Zs2RIZGal1URpbv359cnKyK7dPPP300w0bNgwICDh48OCyZcu2bt3qsp+K7OzsiIiIxx9/PCQk5IMPPnj++edfeeUVrYvSTE5OTlBQ0JEjRyrmG1qK7YJw69atp06dqvi1b9++Ltsb9N133xmNxopf27Rp45rHA4low4YNe/bsKSgoCA4OHjVqVOPGjbWuSHuZmZnbt28fN26c1oVoZufOnZs3by4oKAgKCho1apQrjwiJ6PLly998801hYWG/fv1iYmK0LkdLGRkZe/futcY/DVyGCQAAXBpmjQIAgEtDEAIAgEtDEAIAgEv7fwmANF9mGs2XAAAAAElFTkSuQmCC", "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.3" }, "kernelspec": { "name": "julia-1.5", "display_name": "Julia 1.5.3", "language": "julia" } }, "nbformat": 4 }