{ "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 Δtime\n", "--- --------------- --------- --------- ---- ---- ------\n", " 1 -7.774283520076 -0.70 0.80 4.5 \n", " 2 -7.779029236222 -2.32 -1.52 0.80 1.0 37.2ms\n", " 3 -7.779318521323 -3.54 -2.59 0.80 1.5 37.6ms\n", " 4 -7.779350260492 -4.50 -2.92 0.80 2.8 46.2ms\n", " 5 -7.779350714885 -6.34 -3.24 0.80 1.0 36.8ms\n", " 6 -7.779350850454 -6.87 -4.85 0.80 1.0 42.5ms\n", " 7 -7.779350856105 -8.25 -4.87 0.80 3.0 49.1ms\n", " 8 -7.779350856131 -10.59 -5.27 0.80 1.0 53.6ms\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+gvaeTAAAgAElEQVR4nO3deUBU5f4/8Oc558wAAyggIrKoLCqKZIIoIuUCaCooamiaWlbqLTMt++at+8tsvXZvWd7UXEpLy+uegeACKpqi4nVLXBAVFTcEBdkGmLP8/hhFRJMZZpgzy/v113CcOedjx+bNc57Pcw6VJIkAAADYKkbuAgAAAOSEIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJuGIAQAAJtmRkFYXl4+e/ZsuaswCxqNRu4SQG84axZHkiRBEOSuAvQjCILR7wxqRkFYVFS0atUquaswC1VVVXKXAHrDWbM4giDU1NTIXQXop6amRhRF4+7TjIIQAADA9BCEAABg0xCEAABg00wUhIIgHDp06Pjx43j8IQAAmBXONIcZNWqUm5tbaWmpp6fn/PnzG72fMg1RMMSerb+9sIq0tDeoQgAAsE2mCMLDhw8XFxdv3LhRkqQuXbpcvXrVx8encbvadV1celbcGMPVzcJlZ8VdN6T/9nskHgEAABpiikujx48f79GjByGEUhoWFvbnn382elfD2jLP+zHD03g1f2/Ljzni5sviimeRggAA0BimCMLS0lKVSqV9rVKpSktLDdnbxA7MKH9mRDqv5smPOeKmS/UHiAAAALprqkujgiAUFxe7uLhwHNeqVatDhw5pt9+6datVq1YG7nxiB0YipMsmjZ8TTRnI2SEFAQCgsQwaEaalpcXGxnp6eg4ePLju9n379rVr1y40NNTHx2f79u3R0dE7duyoqKgoLCw8ceJERESEYTUTQoggEl4gmQXSiWK0oQIAQOMZNCJ0dHScPHny6dOn09PTazeKovjSSy99+umnL7/8clJS0oQJE65cufKPf/xjwIABlNIFCxY4ODgYWLR2XjBnlGLmIeGZJH5lH3Z0ABZEAoDMdu3atXjxYrmrsFoDBw589dVXm2LPBgVhZGQkIWTZsmV1N2ZmZpaWlo4fP54QMnTo0LfffnvHjh0TJkyYMGHCk/dWWVlZUFDQv3//2i2jR49+8cUX671t5UU26SrzS28NryZfPkVcWfblveRIQfXsp6zn5rkVFRWUUrmrAP3grFkcnuc1Go0R77u9b98+jUYzduxYY+0QamVmZm7btm306NFqtVqpVLKsrlNi9vb2HNdA0hl/jvDixYuBgYG1VXbs2PHixYu6fNDBwaF58+YffPCB9kdKaUhIiJOTU933bL4sbr0p/j6Qs2PttFs+60laOIqfHaP51crlz7AOJloY2bQkSar3Fwfzh7NmcbRBaPg1qlpKpTIoKCgxMdFYO4Raoij+9ttvTk5OLMvqFYS6MH5u1O0RJYQ4OTmVlJTo8kFKqb29fUxMzBPeE+PFDPJh6nXHvN2FGe3PzMoSolP532K5Vkb7Vw0AANbP+FNrLVu2rJt8xcXFhreJ1nJSkMf2iHqpyMq+7LC2TPhm/mgR2mcAAEBXxg/CLl26nDlzpqKighDC8/yxY8dCQkKMfpRHUUJmdWXmRTCDtvObLxv5aVUAAGCtDArCoqKi9PT0M2fOlJSUpKena28ZExwc3KNHj5kzZ+bm5n7wwQc+Pj69e/c2UrUNe96P2f4cN+OAOOeo9fTOAABA0zEoCC9fvrx06dKrV6926tRp6dKlaWlp2u1r1qypqKh4/vnnL1++/PvvvxujTj083YJmDmVT86Wxu4XaO7EBAJitCp7U/PVlLImQkhoTVmN7DGqWCQsLW7du3aPbPT09V61aZcieDeSlonvjuFf3on0GACzAhAxBkMi6aFb5yNhEImT6AWHfTel/CRyD5TlNw2rXoduz5Jd+7LC2TI/f+WO30T4DAObro1Bmf4E4PI2vfnhKRyJkxgFh6Vnx0+5s41LQ39//1q1bRilSFMX09HTDl11qNBp3d/fKysonvKdbt265ubkGHkh3VhuE5H77zFc9mee28b+jfQYAzNVTbnTnYC6rUBqR/iALtSm45Ky4MYYb4tvIwWB1dbWxHofO83xsbKxarTZwPyzLzp49W6lUPuE9NTU1omi6L22rWH/+RIl+TGAzmpAmHLstzQnF/bkBwBxpszA6lR+Rzm+K4ZSsEVLwUevWrTt69Ki/v/9LL71kZ3fvtiQ7d+7ct2+fu7v7mDFj3NzcCCEZGRktWrTQNvxfu3bt0KFDI0aM2Lx5MyFkxYoVdnZ28fHxrVu31uWIlZWVS5cuLSwsjIiIiI+P1260t7/3IPWkpKSQkJA//vgjNze3f//+/fr1M9bfVC/WPCKs1a0FPTCUTbkivbhbqEIzKQCYpbrjwmmZxk/B6dOnL1q0KDAwMDk5OS4uTrvxiy++eOONN1q2bJmdnd2tW7fbt28TQn7++eedO3dq33D27NkvvviicUesqqrq06fP6dOn27Rp8/e///3zzz8nhPA8P2XKlOrqakLIv/71r+HDh+fk5Li4uCQmJu7atcsIf0/9Wf+IUMtLRf+Iv98+E8N5oH0GAJpe+jXpyz/1++07oBndfk2SJOmpFvTbbOHbbD0+O6E9Mz7w8cOb/Pz8lStXXr58uVmzZi+99FJAQMDevXtDQ0M/++yzw4cPBwcHE0ISEhIWLlw4e/bsx+4hISGBEDJx4kTdbyW4YcMGZ2fnpUuXEkKeffbZ8PDwmTNnMsxDFSYkJMyZM4cQUlxcvHnz5rq3mzYZWwlCcr995uOjQq8kPmkAG+yKBiwAaFph7nTWU3rMyEiELDotMkRSsMSRpTNDGE6fL6pOLn/5RydPngwKCmrWrBkhRKFQdO/e/eTJky4uLnZ2dtoUJIRERUVlZWXpcTxCCCGFhYXvvPMOIYRhmJ9//rnuH/3555+1z93r1KkTx3F5eXkBAQF13/P0009rX/j4+OzZs0ffoxuFDQUhIYQSMieU7ewi9k/ll0WxQ9vaxJVhAJCLqx2J8dY1yrTdMVuvir/Fcr6OJDqV/+6UtCnGOM8er6ioqKqqqv1RrVarVCqVSqXtptE+OEW7kRDCsmxtd+iT2zsJIY6OjiNHjiSEPPr0lYqKitqmGFEUq6ur696JWsu4t89uHFtMglH+zO+x3BuZ4pcn0EoKAGahXo/oY/tIDXTy5Mns7GxCyNWrV/fv3x8VFdWuXTsPD48NGzYQQtRq9bp16/r27UsIadu27dGjRwkhkiRt3LhR+3GlUqlSqe7cuVNvtyqVKiEhISEhYdiwYY8eNCkpSXvHzd9++83Ly8vX19c4fxmjssUgJIREeNBDQ9n1eeJrfwhPuKEDAIAJPHalhNGzMCQk5LXXXhs6dGh4ePj777/fvn17juOWL18+Y8aMQYMGhYSEhISEaB8l+8orr/zxxx99+vTp3r173WUM06ZNCw8P7969+7Fjx3Q8aMeOHXv37h0XF/f6668vWbKk3gShmbCtS6N1eTvSPXHc+Ayhfwq/Ce0zACCfdw4Ky3LE5AFc7MPXUWvXVIzeJWyKaeSa+lotWrRITU09fvy4j4+Pl5eXdmP//v3Pnz9/5syZli1b1g7XvL29z549e+rUKR8fn5YtW9ZeHZ07d+6nn35aXl7u7Oys40EjIyOnTp167ty5Tp06aWcolUplQUGBo6MjIWTr1q21SykmTpyojWHTs90gJIQ4cmRjDPvxUSEymf89Fu0zACCPXh50sG/9FNTSZuGeG5JR7q+mVCp79OhRb6ODg0NoaGi9jSqVKjw8XPtaG2BaCoXC1dVVr4O6urr27Nmz7hYPDw/ti7qBamdnV7u00cRsOgjJ/faZTi5i/1T+h2fY+DbmOGwHAOs2yv9J3zxPudGn3BoZgz4+PtpulKCgoBdeeKFxO2m02NjYuiGqOy8vL4VCYfR6/oqtB6HWaH+mrRN9fqdwupjM6oosBAArcejQIe0L7RSgiY+uXXrYCLXPMjINfOnfU9s+MwntMwAAtgRB+IC3I80YwhVVkehUvrCq4fcDAIAVQBA+xElBNsWy0V60VxJ/ugQPbwIAsH4Iwvq07TOfhTH9UvgtV5CFAABWDs0yj/dCANPWmT6fLrwVzKB9BgB0wXHcTz/9ZOJGDxtRXFxce9tSo0MQ/qVeHvTQMHZYmnChTFoYySqQhgDwRFOmTJHriXq2wMfHp4n2jCB8Eh9HumcINy5DiE7lN8ZwLe3lLggAzJizs3NYWJjcVYDeMMxpgJOC/BbL9veikUn8GbTPAABYHQRhw7TtM5+EMX1T+JR8ZCEAgFVBEOpqTADzWww3+Q8BD28CALAmCEI9RLaih4axay+KU/YJGqQhAIBVQBDqx8eR7o3jCtRk8Ha+uFruagAAwGAIQr1p22d6t6I9fufPon0GAMDCIQgbQ9s+8/7TTN8Ufud1ZCEAgAVDEDbeKx2YjTHchAzhP6cwYQgAYKkQhAbp3Yr+Ec8uPYv2GQAAS4UgNJS/Mz0wlLtRSYZs50tq5K4GAAD0hCA0AmcF+S2WDXVH+wwAgOVBEBoHS8nccHbWU0zfFH4X2mcAACwHgtCYXu3IbIzhxmXwC05jwhAAwDIgCI2sdyu6L577/jTaZwAALAOC0Pj8nenBYdyNShKH9hkAALOHIGwS2vaZbu60x+98zl1MGQIAmC8EYVPRts+89xTTdwvaZwAAzBeCsGm91pFZH82Ny+AXon0GAMAsIQibXJQn/SOOW3hanLJP4JGGAABmBkFoCgHN6KFh3PVKKW4H2mcAAMwLgtBEnBVkcyz3dAvaE+0zAADmBEFoOtr2mXdDmL5b+N03kIUAAGYBQWhqk4KYddHci7v5RWifAQAwAwhCGTzjSffGcQvQPgMAYAYQhPIIbEYPDuOuVUhxO/i7aJ8BAJAPglA2zRTk9wFcx+Y0KpnPK8OUIQCAPBCEcmIpmd+LnRbMRCbxK3PFCr7+GwSJHL+NjAQAaEIIQvlNDmJW9eXeOiCEb9bUzUJBIuMzhJ24PRsAQFNCEJqFGG/6vwSuVEOC1mu0K+4Fiby8R+jiSmeG4BwBADQhfMmai8Bm9OQIhauSBKzV5FcyL+8ROrnQD57GCQIAaFr4njUjrnbkyHBFqDsNTVX6OCIFAQBMAV+15oWhpKU9DXQWF54WCtRyVwMAYAMQhGakdl7wwMCaZzxp0HpNYZXcNQEAWDsEobnQpmDtvGDyAK6zK+28QYOnVQAANCkEobn4+qQY4vZgXpChZPcQztuRRCU/sroQAACMB0FoLmaGMO899dDpUDIkM17hrCAfHRHkqgoAwOohCM0FSx+zUcWR5AHcujzp22zcnBsAoEkgCM2duz3ZMYj9NltcmYssBAAwPgShBfB1pCkD2feyhG1Xcbs1AAAjM1EQfvzxx2FhYYMHDzbN4axPsCv9PZabkMFnFiALAQCMyURBGBcXt2LFiqtXr5rmcFappwdd0Ycbmc6fKUEWAgAYjYmCMCwszMXFxTTHsmJDfOlXPdmBW4Ur5chCAADj4Iy7u4yMjOrq6rpb+vTpY29vb9yj2LIXA5k71WTwdmFvHOdmJ3c1AACWz8hBmJmZWV5eXndLREQEgtC4pgUz+RXS4O38zsGco5FPIACAzdH1ezQ1NfXgwYPnzp2bMWNGRERE7fZt27bNnz9frVaPGTNmypQpH3zwQdPUCQ/5sgc76Q8hIY1PGcgp0fkLAGAAXb9EV6xYUVFRsXfv3vz8/NqNJ0+eHDVq1MSJEz/55JMvvvhi9erVf/XxefPmvfHGG1euXBk1alRGRoaBRQMlZEkU20xBJ+4RREwXAgAYgEqSHt+jISEhs2fPTkxM1P745ptvCoLw/fffE0J++OGH5cuXZ2ZmPvaDBQUFFRUV2tctW7Z0dnZ+9D3nz5/v2bPnzJkza7fExMR069ZN9/KsRllZ2WP/E9Wj5smQdBLiSub3NEFR0AAdzxqYD57nNRqNg4OD3IWAHtRqtVKpZFlWx/ezLMswDQz5DJpiOn78+KRJk7SvIyIi3nrrrb96Z6tWrRrcmyAIgiDcuXOndktRUZEg2OJtNrX/KRp8m5KSDX3IgDTmX3+SmcG474zMdDxrYD6E++QuBPSg7ylrMAWJgUF469at2kURbm5uarX67t27zZs3b9ze7OzsXF1dv/rqK0NKsg4ajUbHDiN7e7J1kBSVLHg6ca92xGyhnHQ/a2AmeJ5nWRZnzbJIkqTXiFAXBgWhs7NzZWWl9nVFRQXLsk5OTsaoCvTgpaJpg9hntwguSjLSD1kIAKAfg74327Zte+HCBe3r8+fP+/j4GDelQUcBzejW59g3M4W9N9E5AwCgH4OCcOzYsatWrSorK5MkafHixWPHjjVWWaCvp9zo6v5c4k7+2G1kIQCAHnQNwtjYWEppdnb2qFGjKKV79+4lhAwfPjwyMjIgIKBdu3YlJSXvvfdeU5YKDejXmi6JYods53PvIgsBAHSl6xxhWlraoxtZll2xYkVRUVFNTY2Xl5dRC4PGSGjLFKjJoO3CvnjOEz3hAAA6MMIdutzd3Q3fCRjLlCDmZiUZuJXfE8e5KOWuBgDA7KHJ0Ap9FMr086LD0/gqrI8CAGgIgtA6zevJtnKgL+wSeKyzBwB4IgShdWIoWdWXrRKkqZkYFQIAPAmC0GopGLIphjt5R5p9BFkIAPCXEITWTMWRpAHc+jzpm2xcIQUAeDwEoZVztydpg9j52eLKXGQhAMBjIAitn48jTRnIvpclbLuKhfYAAPUhCG1CsCv9PZabkMFnFiALAQAegiC0FT096Io+3Mh0/kwJshAA4AEEoQ0Z4ku/6skO3CpcLkcWAgDcY4RbrIEFeTGQKa4hg7cJf8RzbnZyVwMAYAYwIrQ5b3Zm4trQQdv4co3cpQAAmAEEoS2a24MNcaPD0/kaLKkAAJuHILRFlJAlUWwzBZ24RxAxXQgAtg1BaKNYSn7py+ZXSG8dwA3YAMCmIQhtlwNHkgdw+wukuSdwhRQAbBeC0KY1V5KUgeyys+IPOchCALBRWD5h67xUdMcg9tktgquSjPTDL0YAYHPwxQckoBnd+hz7Zqaw9yY6ZwDA5iAIgRBCnnKj/+3PJe7kj91GFgKAbUEQwj19W9OlUeyQ7fy5u8hCALAhmCOEB4a1ZW6qyeDtwr54ztNB7moAAEwCI0J4yJQgZkJ7ZsBWvqRG7lIAAEwCQQj1ze7GRHvR4Wl8FZbaA4ANQBDCY8yLYNs40dG7BB7LCwHA2iEI4TEoIT88w9YI0tRMjAoBwMohCOHxFAzZGMNlF0sfHkEWAoA1QxDCX1Jx5PdYbmOe9E02rpACgNVCEMKTuNuTHYPY+dniz7nIQgCwTghCaICPI00ZyM7KErbmY6E9AFghBCE0LNiV/h7LvbSHzyxAFgKAtUEQgk56etDV/bjnd/JnSpCFAGBVEISgqxhv+lVPduBW4XI5shAArAfuNQp6GBvAFFeT2K3CvjjOAzcjBQCrgBEh6GdqZ2ZEOxq/gy/XyF0KAIAxIAhBb/8MZ0Pc6PB0vhpL7QHA8iEIQW+UkCVRbHMlnbhXEDFdCAAWDkEIjcFSsqoPe61CmnYAo0IAsGwIQmgkB44kD+QOFEj/PIGbzgCABUMQQuM1U5Btz3E/nRN/yEEWAoClwvIJMIiHA0kdyPZJEVyU5Hk//F4FAJYH31xgqIBmNHUgOzVTSLuGzhkAsDwIQjCCp9zo2v7cuAz+2G1kIQBYGAQhGEff1nRpFDtkO3/uLrIQACwJ5gjBaIa1ZQrUZPB2YV8854kbsAGAhcCIEIxpchDzUntmwFa+uFruUgAAdIMgBCP7sBsT7UWHp/NVWGoPAJYAQQjGNy+CbetER+8SeCwvBACzhyAE46OE/PAMWyNIr/4hoHMGAMwcghCahIIhG2O486XSh//DFVIAMGsIQmgqKo4kDeA2XZLmncQVUgAwXwhCaEIt7MiOQezcE8Ki0/Wz8FqFlJqP66YAID8EITQtH0f6RTg746Cw9uKD2LtaIQ1LE9zsZKwLAOAeLKiHJvdaR0YtkPEZvLOSG+xDr1ZICWnCgkg2woPKXRoAAIIQTGJaZ6ZSQxJ28Oui2c+OiUhBADAfCEIwkVldmcvl0vPpQqIfY8fKXQ0AwH0IQjCRqxVSVqH07lPsL+eFQ+kSS8kofzrKn+nqhqEhAMgJzTJgCrXzgnPDmeXPcq0cyE99WAVDRu0U/Nby0w8I+26igxQA5GGKEaEoisuWLdu/fz/HcRMnTnzmmWdMcFAwH/W6YwZ4U0LYdw8JWwZyc0LJqWJpfZ74yh8CL5L4NjTRj4nyxBgRAEzHFEGo0Whu3bo1c+bM8vLyF198cefOnQEBASY4LpiJpWfFhZFszzrdMQO8qSixK3PFt7swwa402JWtTcSJewUlQxL96Sh/prMLEhEAmhyVJJNekhoxYsSUKVMGDhz46B9dunSpX79+eXl5pqzHPJWVlTk7O8tdhWy0ifjLecmOIYn+dLQ/08kSEtHGz5ol4nleo9E4OODhmZZErVYrlUqWNWbHnUmbZc6cOXPmzJnevXub8qBgceqNEeO2C/YsSfSnL/gzQZaQiABgWYwchImJiRcvXqy7ZfXq1R07diSE5OfnjxkzZuXKlU5OTsY9KFgrbSLO7kYyC6T1eWJ0quCiJIn+dEwA07E5EhEAjEPXINRoNKdOnTp//ny3bt3qzvCVlpauXbu2tLR08ODBnTp1Wr9+/WM/fv369WHDhi1cuDA8PNwIVYMtYSiJ8qRRnuw3EfcSsV8K76qkif50bADTAYkIAIbRdflEjx49EhISJk2alJaWVruxsrIyIiIiNTX11q1bkZGRGRkZj/1sVVVVdHR0ZGTk9evX169ff+nSJYPLBlukTcT5vdirYxRLotjiatI3hQ/ewM85KuTexeoLAGgkXUeEe/fudXZ27tevX92Nq1evdnZ23rRpE6W0devWn332Wd++fR/9rCiKL7/8MiFEe9U0KCjor47C83zdK6uenp4qlUrHCsF21Bkjstox4rNbeDc7muhPxwUygc0wRgQAPegahI9th0tLSxsyZAillBASFxc3c+ZMjUajUCjqvU2lUs2aNavBQ1RWVhYUFERHR9dueeONN15//XUdK7QmFRUV2v+q0KCnncjTIeSTLiSriPktn+mdxLgpyfA2wph2op+TSYeJOGsWR9s1Kgh4drQl0bdr1N7enuMaSDqDmmWuX79em1utW7cWRfHmzZu+vr6N25tKpfL29sbyCUKIJEloKdJXrDOJ9SPfSeRAgbQ+j43ZKbrb0QntmXGB1NvRFPmEs2ZxsHzCErEsa17LJyh9sAxR+wK/EYO82PtXTedFsAcKpPV5YthmIbAZTfRjEv2plwr/PgGgPoOC0MvL6+bNm9rXN27cYBimVatWxqgKwFCPJmLob/cScZQ/0xpTzwBwn0E33R4wYEBKSoooioSQ5OTk6OjoRycIAeTF3u81vTZWMSeUPVIkBW/URCXz87PFm2q5iwMAM6DrLda+/fbbzMzMjIwMX1/fgICAd999t0ePHpWVlT179mzbtm2HDh1WrFiRlJRkyA21cYu1WrhZV5OqFsiOa+L6i9KWfLGzC030Y0YHMJ4GzxPhrFkczBFaIjlvsRYREeHt7Z2YmKj90dvbmxCiUqkOHDiwadOmkpKSrKys9u3bG7EygCZix5L4Nkx8G1IlsGnXxPUXpTlHNcGuNNGPeSGAaYVvRQAbY+qbbj8BRoS1MLYwsSqBaBMx+YqoTcQxAYyHnomIs2ZxMCK0RBZ/020A82T/yBjxo6OaLq400Y8ZG8i0tJe7PgBoSghCgAdqE1HNs+nXH0rEFwMZ90cSUSMSxSMNZxIhgkg4gxrRAMB08D8rwGM4cCS+DbOyL3tjrGJWV+ZIkdRxvSZ+B78yV7xb8+Bt4zOElPyHJhckQt7MFP57UTR1xQDQWBgRAjyJNhHj25BKnt15XVx/UXr7oCayFU30YxLaMUui2LgdvCAyQ9syhBCJkLcyBQeWjA/Er5gAFgNBCKAT1cOJuCpXmpapebY1HRfAzP1TJIT0dSNvZQp2LPmqpzGn8QGgqSEIAfRTm4glNezvl8V1F8VTxdIre4U2jso+XkhBAMuDIARoJBcleak981J75k41GZ7OHylkckrFC6VkRDsa35ZpYSd3fQCgG8xkABhEIuSjI0K4O80ZWv10C9rVjey6LgWuvXcXt2sV5rJOFwD+CoIQoPG03THaecFmCil1IJdxU3rej954UTGrK3O6RArdzAdv4OccFXLuIhEBzBQujQI03tsHBRVHvuxxb16wuZIkD+DitvMOHNXOIy7qfe/ZF/1TBBcliW9L43yZKE88DQrAjCAIARrvtY5MF9eHUs1FSbYM5Mo198Z/tU+D+iaCHLstJV8RX/1D0Igkvg1N9GN6e+IBngDyQxACNF69FNRyURIXZf3tDCVh7jTMnZ0TSk4VS+vzxBkHhSvl0nM+TKI/fc6HefQONQBgGghCAFMLdqXBruycUHKxTEq+LH15Qpy4RxjsyyT60wHejB3WXwCYFn4LBZCNvzOd3oXZF88dG8FFedKlZ0WPX+7dyK1UI3dxADYDI0IA+fk60slBdHIQc7uapFwR1+eJ0zKFZ1vT+DZMQlu9HwgFAHpBEAKYkRZ2ZEJ7ZkJ7poInu66L6y9Ks7LuPTT4eT/q7YjeGgDjQxACmCNHrv4DoT45Jvg507g2dGwA06E5EhHAaBCEAGat9vEXgnRvSWLfFN5VSRP9aXwbJswdiQhgKAQhgGWosySRzSyQtuSLY3YLPJYkAhgMQQhgYZj7iTg3/N6SxCn7hDINGeRL49pgSSKA3hCEABbs0SWJr+wVBvlgSSKAHvCrI4A1qF2SeCSBC3On/zn1YEliGZYkAjwRRoQAVqWNE53ehU7vwhRVkdR8LEkEaBiCEMA6udvfW5JYUkPSronJl6X3sjRdXGmiH5PoT71U6K0BuAdBCGDlXJQk0Y9J9HuwJPHjY4I/liQC3Ic5QgBboV2SuKEJr6UAABgBSURBVLIvWzhO8W0EW1xN+qbce27wkaKHnhssEfJellD+yOTiylzx4C08YRisDYIQwOZolyTO78VeHaNYEsVWCeSFXULAWn76AWHfTUkihBLSsyWN28HXbbRZc0Fcc1Hs6oYRJFgbBCGA7dIuSZwbzuaO4pIGsK52ZPI+od0afso+QcmSNzox8fezcM0FceV5cWM054DpFLA6+EcNAIQ8bklizl0p2IVGJvPTOzObLiMFwWphRAgAD6ldkpg1jItvyxSppWkHhBXPIAXBaiEIAeDx/Jypt4p0c6eBziQqhcfCfLBWCEIAeLzaecFNsdytSmngNmQhWCcEIQA8Rt3umPbN6dTOjIKS+B3IQrBCCEIAqE8iJOcuqdsd82E3Nr+C9G1Nz5ZgHSFYGwQhANRHCfkolKnbHePAkW8imPUXpadbYB0hWBsEIQDoZFhbJqAZ/SZblLsQACNDEAKArhZEMl/9KVwqw9VRsCoIQgDQVRsnOi2YfTcLg0KwKghCANDDrK5M9h0pJR+DQrAeCEIA0IOSId9FstMPCFWC3KUAGAmCEAD0E+tNQ1vQL0/gAilYCQQhAOhtfi924Wnh3F1cIAVrgCAEAL21VpFZXdlpmbg8CtYAQQgAjTE9mClQk02XcIEULB6CEAAag2PIgkj27YNiOe4+ChYOQQgAjRTlSfu2pp8dxwVSsGwIQgBovH/3YJfniCfvoGsGLBiCEAAaz8OBzAllp2YKSEKwXAhCADDI3zoxNSJZfR5dM2CpEIQAYBCGkoWR7HtZYkmN3KUANAqCEAAMFeZOh7alHx1B1wxYJAQhABjBF93Z9XnisduYKwTLgyAEACNwtSOfd2en7BNERCFYGgQhABjHyx0YFUeWn0PXDFgYBCEAGAclZEEk+//+JxRWyV0KgD4QhABgNF1c6dgA5oPD6JoBS4IgBABj+jiM3X5VOnALU4VgMRCEAGBMzgry757MlH0Cj7lCsBCmCEJJkg4fPvzrr79u2bJFrVab4IgAIKPR/oyPI1l0BkkIlsEUQVhTU7N48eLCwsL09PQePXpUV1eb4KAAIKP5Eexnx4TrlbhAChaAM8Ex7OzsfvzxR+3rnj17XrhwoXPnziY4LgDIpX1zOjmIeS9L/KUvK3ctAA0wRRASQnieX758eX5+vr+/f8eOHU1zUACQ0T+eZrts5Hddl/p7UblrAXgSYwahRqMZN25cvY1r166tfV1VVXXr1q2qqipHR0cjHhcAzJADRxZEsq/vF/4cwdlhWAhmzAhBePnyZYZhfH19OY77/PPPH38Yjps8eTIhZPz48du2bRs5cqThxwUAMzfIl35/hn6bLc7qigZ1MF86/evkeX7QoEEeHh6U0kOHDtVur6ioiImJ6dOnT+/evePi4qqrqwMfQQgpKioqKCgghBQWFh47dky7EQBswYJI5uuTwqUydM2A+dIpCCmlEyZM2L9/v7Ozc93tixYt0mg058+fP3/+fFFRUW1HTD0lJSXjxo0LCwsbPnz4rFmzunbtaoTCAcAStHGibwWz7xzCUgowXzpdGmVZdsyYMYQQSh+a9F6zZs2MGTM4jiOEvPLKK6tWrZo6deqjHw8MDExLS2vwKJWVlTdu3OjWrVvtlpdeeunVV1/VpUIrU15eLncJoDectb/yuj+JzLVbn1P9nJd53XqN53mNRsPzvNyFgB7UarVSqWRZXaed7e3tFQrFk99j0BzhpUuXaq9zBgYGXr582ZC9qVQqd3f3H374oXZLhw4d6o1BbYfN/sUtGs7aX1nYW/rbfiY+kLM3p64ZbRA6ODjIXQjogeM4vYJQp30a8uGKigp7e3vta5VKVVpaamA1CoUiLCzMwJ0AgLmJ8aZh7nTuCWFOqDklIQAhxMA7y7Rq1aq4uFj7+s6dO56ensYoCQCs0Pxe7KLTYs5ddM2A2TEoCLt163bw4EHt6wMHDoSGhhqjJACwQp4O5O9d2bcyzWuaEIDofml07dq1d+/eramp2bx584kTJ8aMGePs7Pzmm2++8MILXbt25Xl+4cKFycnJTVorAFi0t4KZlbnixjxxpB+WFYIZ0TUIs7Ozb926NWHChDt37ty5c0e7Ij4mJub7779fsGABwzA//fRT7969m7JUALBsHEMW9mbH7hYG+jBODfTxAZgOlSRzuWR/6dKlfv365eXlyV2I/MrKytB/aHFw1nQ0ca/gYU++7CF/1wy6Ri2RvssndIELFABgUl/1ZH/OFf+8Yy6/ggMgCAHApFrYkY9C2TczBSQhmAkEIQCY2pQgpkYkv57HfdfALCAIAcDUGEoWRrKzssSSGrlLAUAQAoAswtzpsLb0w/9hWSHID0EIAPL4Zzj722UpqxBzhSAzBCEAyKO5knzenZm6XxARhSArBCEAyGZCe8ZJQX7IQdcMyAlBCACyoYQsiGRnHxEKq+QuBWwYghAA5BTsSscFMu8fRtcMyAZBCAAy+yiU3XFVyizAVCHIA0EIADJzVpCvejJ/2y9oMFcIckAQAoD8RvkzPo5k0WkkIcgAQQgAZuE/vdjPjwvXKnCBFEwNQQgAZiGwGZ3Sifm/LAwKwdQQhABgLj7oyh66Je28jkEhmBSCEADMhQNHFvZmX98vVGMxBZgQghAAzMhzPrSTC52XjQukYDoIQgAwL9/1Yr45KeSV4QIpmAiCEADMSxsnOqML+85BDArBRBCEAGB23n2KybkrbbmCQSGYAoIQAMyOkiHf92bfzBQqeblLARuAIAQAc9SnNe3Vin55Av2j0OQQhABgpr6JYL8/I+bcxQVSaFoIQgAwU54O5P2n2WmZGBRC00IQAoD5eiuYKawiG/LQQQpNCEEIAOaLpWRBJDv9gFiqkbsUsF4IQgAwa71b0YE+9NOjuEAKTQVBCADm7t892V/OiyfuoGsGmgSCEADMXQs7MieMfTNTQBJCU0AQAoAFmNSR4UWyKhddM2B8CEIAsAAMJQsj2b8fFoqr5S4FrA6CEAAsQ6g7Hd6O+fAIumbAyBCEAGAx/hnObr4sZRVirhCMCUEIABajmYJ80Z2Zul9A2wwYEYIQACzJ+PaMs4L8kIOuGTAaBCEAWBJKyHeR7Owjwi213KWAtUAQAoCFCXalEwKZvx9G1wwYB4IQACzPx2Hs7hvSnhuYKgQjQBACgOVRceRfPZg3MwUN5grBYAhCALBIiX5MGyey8DSSEAyFIAQASzW/F/vFceFaBS6QgkEQhABgqQKb0b91Yt7NwqAQDIIgBAAL9n5XNuuWtO0qBoXQeAhCALBgDhxZ1Jt964BQjcUU0FgIQgCwbAN9aLAL/fokLpBCIyEIAcDi/SeS+TZbyCvDBVJoDAQhAFg8X0f6dgj79kEMCqExEIQAYA1mhjDn7krJV5CFoDcEIQBYAyVDFkex0zLFCl7uUsDSIAgBwEo860kjW9G5J9A/CvpBEAKA9ZgXwS49K54tQdcM6AFBCADWw9OBfNCVnXYAg0LQA4IQAKzKm8FMURVZdxFdM6ArBCEAWBWWkiVR7DuHxLs1cpcCFgJBCADWpkdLOtCbfnoMF0hBJwhCALBC/+rJ/npePHEHXTPQMAQhAFihFnbk4zD2zUwBSQgNQhACgHV6rSMjiGRlLrpmoAGmC8KcnJzWrVuvWbPGZEcEAFvGULIgkn0vS7hdLXcpYN5MFISiKL777ruRkZHV1fgnCQAmEupOR/kxH/4PXTPwJCYKwu+++2748OHe3t6mORwAgNbn4WzSFenQLcwVwl/ijLivgoKCnTt31t3SsmXL2NjYS5cupaenJyUlTZ8+3YiHAwBoUDMF+Wc4MzVTODSMY6nc1YBZMmYQVldX37hxo+4WhmEIIbNmzerUqdOyZcuys7PLysp69erVoUMHIx4XAOAJxgcyK3LEZWfFv3VCeyA8hk5BeOvWrZUrVx45cqS4uHjbtm2129Vq9fvvv799+3YPD4+PP/64b9++M2fOfPTj48ePv379utFKBgDQ04LebP8UfkQ7xsNB7lLA/OgUhDdu3Dh79qyfn9+GDRvqbv/oo49OnDiRnJyclZU1bNiw3NxcDw+PRz8eFxenfZGdnR0WFobhIACYWGcXOqE9M+uwsOJZVu5awOxQSdJ1DvnkyZOhoaEajUb7o0ajadWqVUpKSq9evQghgwYNio6Ofvfdd5+wh8rKSo7jlErlY/80Nzc3PDx80qRJtVsGDx4cGRmpY3nWpKyszNnZWe4qQD84a2aukifdkpllkdKzre596fE8r9FoHBwwSLQkarVaqVSyrK6/0CgUCu0k3RM0fo7w+vXrxcXF3bt31/7YvXv37OzsJ39EpVI94U8ppQzDuLi41G7R/a8KAPBkKo7MDRVnZDGHhkgKzBVCHY0PwsLCQpVKpVAotD+6uLgcPXrUoFI4rnnz5v/4xz8M2Yl1qKmpsbOzk7sK0A/Omvl7oQP5JY9fel7xTghDCGFZlmEYnDXLIoqiXiNCXTT+9yIXF5eqqipBuLdStayszNXV1UhVAQA0ifm92LknhGsVWFYIDzQ+CL28vBQKRW5urvbHc+fO+fn5GakqAIAmEdCMvtGZmXkINyCFB3QKQkmSiouLS0tLCSHFxcV3794lhKhUqpEjR3799deSJOXk5KSkpIwbN65piwUAMEyVQCp5crRI2pr/0KBwfrZ4vhTDRBulUxDevn07ICAgPj7e2dk5ICCgtpPz3//+9+nTpz08PCIiIj7//POOHTs2ZakAAIayZ4mXinqq6FsHhar7tyD9+qSYeUtq52TZN54priY/nXvMSHfFObGkxvTlWBKdmmXc3d3v3Lnz6HYvL6/9+/eXlZU5ODhwnDFvUgMA0ERmdGEIEb86KX2VTWZ1Jl+fFLMKpV/7spyFt5I2U5Kd16XrleIHTz/4m2j/duMDZazLAhghvbB2CgAsy4wuTEmNNPeEWF7N5FdZQwoSQlhKfurDvrxH+OL4vSy0moxvanosqG9qly5d6tevX15entyFyA9Lsy0RzprFCdnAZ5fc+wJ0URJ6/8qoiqN295ODY4iz4sFHmisJc/9tDiy1v9/DzzKkWZ23NVMQ9v4e7FnicP9u3zq+jaGkeZ37jjgrSG2S2TFExd17G6XE5ZHbkwgSeXmP0MmF2rHEKlNQ3wX1usD1TACwRV+fFDu5knH+4v+KudX92Aqe1A4KKnmp+v5cGy+SMs2DT92tIaKeb1MLUnHNvR8EkeTWeVuphgj391AlEPX9H0SJ3K0zq1emIfz9t1WLpJK/tzdJInUn/5wU5N6NAiSy7iJRceSHZ3BPEp0gCAHA5mivGa58hkiCuOQCHbtbWN2Prb3djKudjl0z5tVcU64hGpEQQhaeFrOKpDtV0ruHxOkHxUkdmdc6Um9H86rWrCAIAcC21M6cEZHXCPd6Z+ploSVyUhBCyNcnxZPF0qZollLy8h6hpT1RC1K334SuLejkIGZEOwYPZXyUJZ92AAA9VQuktEZa3e+hmbMZXZj+XvRKubk0TDRa3e4Ybe9MYRVppqBXxigmBzHfZotB6/kvT4iFVXIXamYQhABgQ+xY8nEY++io6PVOTEAzyx4rFVaRC6UPZTxLyYpn2asVUrmGJPox++O5Nf3Zi2VSx/WaUTuF9GsWH/zGgq5Rc4T+Q0uEs2ZxbPYxTHdryNqL4nenRF4iL7dnJgcxrpZz4/Gm6BrFiBAAwLY0V5LJQczJkdwvfdlTxVLgOs2UfcLx2+YyKDI9BCEAgI0Kc6cr+7Knnlf4O9Ph6UL3zfzSs2IlL3dZJocgBACwaZ4OZFZX5sIobm44m35NavNfzZR9wukSGxogIggBAIAwlMR403XR7IkRXGsV6Z/Cx27l1+eJGht4YhWCEAAAHvB2pHNCWe2Ki6VnxXZr+L8fFqxgbckTIAgBAKA+JUMS/Zi0QdyuISwhJGwzH7uVT74iWmUeIggBAOAvdWxO54az+WMU4wOZOUfFDuv4L0+IRda1JB9BCAAADbBnyYT2zJGEe0vyO1jXknwEIQAA6CrMnS6JYvNGK2K86dsHhc4b+C9PiMXVcpdlGAQhAADop3ZJ/iqrWJKPIAQAgEayjiX5CEIAADCIpS/JRxACAIARWO6SfAQhAAAYk8UtyUcQAgCA8VnQknwEIQAANKG6S/I/OmKOS/IRhAAA0OS0S/KPDjfHJfkIQgAAMJ1Hl+TPzxbLNXKWhCAEAABTq12SvzSK3V8gtV2jmbJPOHFHngEiJ8tRAQAACCFRnjTKk72pZn8+JyakCS3syOQgZnwg42DCdMKIEAAAZFZvSb7vw0vyc+5Knx6rvxqxSiBTM4VqwQhHRxACAIBZqF2Sf/zhJfn+zvSmWnov60Ho1Yhk1E4hxJXascY4rhH2AQAAYDw+jyzJd1aQAjXRZmGNSJ5PFwb70r91Mk6EYY4QAADMkXZJfqIfk10sfX9GTLkiNrej5+8yGlEa0oYxVgoSjAgBAMDMdXGlCyPZSy8o3u7CbL3GFFdLRkxBgiAEAACLoGTJjqvSl6HiUy1o3flCwyEIAQDA3NXOC05qL/4nglbwxIhZiCAEAACzVq87hhKyIJI1YhYiCAEAwKxdq5BG+z/UI6rNwlYOFOsIrda8efN4npe7CtBDaWnp4sWL5a4C9JOdnb1lyxa5q4CG+TnTFwPvpVVycvKpU6cIIZSQmSEM1hFarW+//basrEzuKkAP165d+/HHH+WuAvSTlZWVmpoqdxWgn9TU1KysLOPuE0EIAAA2DUEIAAA2DUEIAAA2jUqSWTwgmBCSm5vbpUsXHx8fuQuR3+XLl319fRkGv6ZYDI1GU1BQgH+9lqW8vLyqqsrd3V3uQkAPRUVF9vb2Tk5OOr5/7Nixn3766ZPfY0ZBSAi5ePGi3CWYherqajs7O7mrAP3grFkcURQFQVAoFHIXAnrQaDQsy+o+TmjdurWDg8OT32NeQQgAAGBiuPgGAAA2DUEIAAA2DUEIAAA2DUEIAAA2DU+oNyOSJB0+fDg9Pf3OnTtdunQZO3asUqmUuyjQVW5u7u7du+Pi4ry8vOSuBRrG8/yaNWuOHDni5uY2cuTIzp07y10RNECSpOTk5P3799vb2w8dOjQsLMxYe8aI0Izk5+ePGjWqpKTE19d38eLFMTExgmDMh09C0+F5fvz48dOnT8/JyZG7FmgYz/ODBw9etGiRj48Pz/MHDhyQuyJo2CeffPLOO+8EBASoVKp+/fqlp6cba88YEZoRLy+v8+fPcxxHCJk4cWLLli2zs7O7du0qd13QsLlz5w4cOPDSpUtyFwI6Wb58+c2bN48ePar93w0swoYNGz777LMXXniBEJKfn79x48aYmBij7BkjQjPCcVzt/5Y8zwuCoPvdE0BGZ8+eXb9+/fvvvy93IaCr1NTU8ePHb9269ZtvvsFw0FJ07tz52LFjhBCNRpOdnR0cHGysPSMIzdT06dNHjBgREBAgdyHQAFEUJ02atHDhQnt7e7lrAV3l5eUtW7Zs48aNZWVlI0eO/O677+SuCBq2ZMmSnTt3+vv7e3l5+fn5TZ061Vh7xmUBc/Thhx+eOHEiIyND7kKgYfPmzQsJCYmKipK7ENADwzAdO3b86aefCCHh4eHjxo2bNm2a3EVBA2bMmOHr6/vTTz+VlpZOnjz5xx9/fO2114yyZ9xizex8/vnnq1ev3r17t4eHh9y1QMNCQ0PVarWjoyMh5M8///Tz8/u///s/Y/3/CU0kPj6+c+fOX375JSHkypUrbdu2LS0tdXZ2lrsu+EtVVVUqlers2bMdOnQghCxfvnzx4sXGekIvRoTmZd68eStXrszIyEAKWopff/21srJS+3rAgAHvvPPOkCFD5C0JGpSQkLBq1SpJkiilBw8e9PX1RQqaOTs7OycnpwsXLmiD8MKFCy1atDDWzjEiNCM5OTlBQUF+fn5ubm7aLV9//XWfPn3krQp05+np+d///rdfv35yFwINUKvV0dHRdnZ2AQEBSUlJy5YtGzZsmNxFQQO+++67OXPmjBgxori4ePfu3UlJSb179zbKnhGEZkStVp8+fbruloCAABcXF7nqAX1pL41ibGERNBrN7t27y8vLIyIicA8ES5GXl3f06FF7e/vIyEhXV1dj7RZBCAAANg3LJwAAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKYhCAEAwKb9f0067pzZO3H5AAAAAElFTkSuQmCC", "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.8.5" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.5", "language": "julia" } }, "nbformat": 4 }