{ "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 = [ElementPsp(el.symbol, psp=load_psp(\"hgh/lda/si-q4.hgh\"))\n", " for el in load_atoms(silicon)]\n", "positions = load_positions(silicon)\n", "lattice = load_lattice(silicon);\n", "\n", "model = model_LDA(lattice, atoms, positions)\n", "basis = PlaneWaveBasis(model; Ecut=5, kgrid=[3, 3, 3]);" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "DFTK already defines a few callback functions for standard\n", "tasks. One example is the usual convergence table,\n", "which is defined in the callback `ScfDefaultCallback`.\n", "Another example is `ScfPlotTrace`, which records the total\n", "energy at each iteration and uses it to plot the convergence\n", "of the SCF graphically once it is converged.\n", "For details and other callbacks\n", "see [`src/scf/scf_callbacks.jl`](https://dftk.org/blob/master/src/scf/scf_callbacks.jl).\n", "\n", "!!! note \"Callbacks are not exported\"\n", " Callbacks are not exported from the DFTK namespace as of now,\n", " so you will need to use them, e.g., as `DFTK.ScfDefaultCallback`\n", " and `DFTK.ScfPlotTrace`." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In this example we define a custom callback, which plots\n", "the change in density at each SCF iteration after the SCF\n", "has finished. For this we first define the empty plot canvas\n", "and an empty container for all the density differences:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots\n", "p = plot(yaxis=:log)\n", "density_differences = Float64[];" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "The callback function itself gets passed a named tuple\n", "similar to the one returned by `self_consistent_field`,\n", "which contains the input and output density of the SCF step\n", "as `ρin` and `ρout`. Since the callback gets called\n", "both during the SCF iterations as well as after convergence\n", "just before `self_consistent_field` finishes we can both\n", "collect the data and initiate the plotting in one function." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using LinearAlgebra\n", "\n", "function plot_callback(info)\n", " if info.stage == :finalize\n", " plot!(p, density_differences, label=\"|ρout - ρin|\", markershape=:x)\n", " else\n", " push!(density_differences, norm(info.ρout - info.ρin))\n", " end\n", " info\n", "end\n", "callback = DFTK.ScfDefaultCallback() ∘ plot_callback;" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Notice that for constructing the `callback` function we chained the `plot_callback`\n", "(which does the plotting) with the `ScfDefaultCallback`, such that when using\n", "the `plot_callback` function with `self_consistent_field` we still get the usual\n", "convergence table printed. We run the SCF with this callback ..." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) α Diag\n", "--- --------------- --------- --------- ---- ----\n", " 1 -7.846724494063 -0.70 0.80 4.0\n", " 2 -7.852285512770 -2.25 -1.54 0.80 1.0\n", " 3 -7.852621281480 -3.47 -2.53 0.80 3.2\n", " 4 -7.852621932194 -6.19 -3.35 0.80 2.2\n", " 5 -7.852621941581 -8.03 -3.99 0.80 2.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+gvaeTAAAgAElEQVR4nO3deWAUVbr38XOqqpMQkkBYE0KAIPsWQ0iEEBZlJ5CwmKgoIiogeu/V0XFw5l6X6zbg3HH0HUcFxxmdUWckEAIBWUJYI1sIBFB2CSIgqyyBrFXV7x+NjDIsWTpd3anv56/uprr6Sdv2r586p+pIp9MpAACwK8XqAgAAsBJBCACwNYIQAGBrBCEAwNYIQgCArRGEAABbIwgBALZGEAIAbI0gBADYmieCsKSkZPLkyYmJiQ888MClS5c88IoAAFSSJ4Lw97//fWRkZG5ubo8ePV5//fUbbXb06NGTJ09WfreGYbijOlSHaZpcnM9CfPitZZqm1SXYV218+XgiCJcsWXLfffcJIe67777FixffaLO33nrrk08+qfxui4uL3VAcqqW8vJzvYgvx4beQ0+nk/bdQaWmp23+IeCIIT5w40axZMyFE8+bNT5w44YFXBACgkjwRhEFBQSUlJUKI4uLikJAQD7wiAACVVKMg3Lhx41tvvTV9+vStW7f+9PFNmzY98MADY8eO/fzzz4UQ0dHR27ZtE0Js27ate/fuNXlFAADcS6vJk1988cUWLVosW7Zs4MCBvXr1cj34zTffDB06dObMmREREY899pjD4fjlL3/50EMPff311+np6e+++647ygYAwD1qFIQrVqwQQtx+++0/ffC9994bN27c448/LoQ4ceLEH/7wh/Xr1y9dujQ/P3/ixIktW7a80d5M09y4ceP777/vuhsVFTV48OBrtjldKpoGXLltGIZrvsZlXUghAmv0p6BqDMOQUkoprS7Epq5++OF5TqeT999CVX3zFUW55TeV+9Nj69atDzzwgOt2YmLiU0895XQ6w8PDR40adfMnmqb53Xff5eXlue5WVFT069fvmm0eWaeltjHTWptCiPLy8rKyssu6SF2rPtvVeWcYE5o9p6ysjO8CC7k+/FZXYVNOp7OsrMzhcFhdiE25PvmqqlZy+4CAAAuC8OTJk40aNXLdbty4cWlp6fnz50NDQ2/5RE3T0tLSnnnmmZtskz5EjM3WVYcysZ1iGIbwC5ywVp/SWUlqyyVyPEpRFE3TNI023BqGYQQGBlpdhU05nU6n08n7byF/f//KB2FluP+LLDAwsLS01HW7pKRESunGT0ygJhYM0cZm60KIYU3k/dn6Q+2V+9uRggCAanJ/hERGRhYWFrpuFxYWNm/e3N/f3437d2XhX/eb/VY4HulACgIAasT9KZKWlvbZZ5+5Thz8y1/+kpaW5vaXcAphOsXJUuXFbcaBC1zoCwBQfTU6NDp+/PjVq1dfvHjxkUcemT59+qJFixITE1NTU9PT07t06RIaGlpRUZGdne2uWl2KdTEmW5/SURncpPiulf6xC/X/jlaf7aEoTGAE4EGLFi2q0lUhUUMjRoyYPHlybey5RkH48ccfV1RUXL0bHBwshHA4HAsWLDh48GBJSUnXrl0VxZ1NZ7EuUn4cFywqcuaNcQxfqn90wMg6Yn40QG0XQhgC8JAvv/zS398/OTnZ6kJsITc3d+XKld4YhEFBQTf6p3bt2tVkzzdyzyp9Skcl7cc5ooGaWDZCS8nWu4UqfRbp/3O7+l/d6AwBeEj37t1TU1OtrsIWKioqlixZUks797Hp7x8N0Br/fOZNoCYWDtGkEFM7yUlrjSXfmR/2VyPrk4YAgErxsSmXja83/zRQE/U00bmh3DBaG9RCiV2gz9nLyfUAgErxsSC8OU0RM6KVVUna7L3myOX6sctMKAUA3EKdCkKXbqFyc7I2IEzpmUlrCMDrmLf6iX7LDeBedTAIxY+tYfYI7b09ZmqOcbrU6oIAQAghxKvbzYlrDOPGUWc6xZT1xn9v5UK+nlM3g9ClRyO5JUXrGiqiMyoyDtMaArDeyEi57Kg5YbWhX+87ySnEExuMzw+ZQyLq8pezt6nj77VDES/1VDMGa7/JM9NyjLNcrx+ApXo2kdkjtJXHzPvXXJuFTiEe/9L4+wFz8TBtYHh1pr6Hh4cXFRW5pU5d11euXOl01vQo7aVLl5o0aXLz/XTp0uXIkSM1fKGaqONB6NK7mdw+VmsbInrM1xd+S2sIwErXzcKap6AQoqysrObR5VJcXDxkyJCaL7Xm5+f3wgsv3HwhJDeWXT0+dh5htdXTxMw4NbmV86F1xqcHnbMT1VB3XgkcAKrAlYVDlur3rxGfDlRVxQ0peA2n0/nZZ5/t2rWrffv2EydO9PPzcz2+fPnyDRs2hIWFTZgwoUGDBkKInJycFi1adO7cWQhx5MiRgoKC5OTkzMxMIcSf//xnRVHGjBnTrFmzyrzopUuX5syZ88MPPyQkJIwcOVIIIaUMCLiynHpGRkZ8fPzKlSsPHTo0ZMiQf19x1iq26AivSmguC8ZqbUNE9wx9yXdMzAJgmZ/2hY/lujkFhRDTpk376KOP2rVrN2/evPHjx7sefOGFF55++umwsLD8/PyYmJgLFy4IIT744IO1a9e6Nti5c+cbb7xRvVcsLi5OTEw8ePBgZGTkU0899eabbwohSkpKpk2bZpqmEOK111679957Dx06FBQUlJyc/OWXX7rh73QHu3SEVwVqYmacOqiF89H1xvCW8vd3qEEsNA3AHXRTjFiuV+kpt4XIeYWmECKmsXytwHitoArPzRyi1b/BV/jBgwczMjKOHDkSGBg4ceLE1q1bb9mypX379r/73e927drlugTmsGHDZs+e/atf/eq6exgzZowQ4tFHH638+tuffvppixYt3n33XSFEXFzcXXfd9V//9V/XbDN69OgZM2YIIU6dOrVw4cK+fftWcue1ynZB6DIkQu4arz272eiRoX/YX73Tfb/CANiWIsWMHlVYOd3pFG9+ZToUpyJFQz/5bPeqraLjf+Mjert27erevbtrUXR/f/+YmJivvvpKCBEaGnr1QtD9+vVzPVglx48fdyVZQEDABx988NN/2rlzZ+/evV23Y2JiysrKvvvuu8aNG/90m9tvv911IyIiYufOnVV99Vpi0yAUQoQ4xOxEddlR56Q1xohI+WZv9Ua/rQCgMhQpBkdUNspcs2PWnzCXDddCHGLIUv3D/eLTgarmjgGrS5culZb+6wTqkpKSwMDAwMBA10qxLsXFxa6kVFX16qSY4uLim+85JCTEdaD13zvFy5cvX12JQdf1iooK1/5/SlWr8EPBY+w1RvjvhreUO8drQoge8/V1Jxg1BOAJ18wRvck5FdWWn5+/f/9+IURhYWF+fn5CQkL79u0DAwOzsrKEEJcvX54/f/7AgQOFEG3atNm2bZsQwul0ZmRkuJ5ev359h8Nx9uzZa3YbFBQ0ZsyYMWPGjBo16t9fdOHCha6snTt3bvv27Zs3b+6eP6aW0QSJhn5idqK65DvnhNXG+DbyjXjV3xt/sgCoI657psQ180hr3hdGR0fff//94eHhmzdvfvnll1u1aiWE+Mtf/jJ58uR333137969AwcOTEtLE0JMmTKlb9++d95554ULF7p06eJ6uqqqjz/+eHR0dMuWLT/++OOuXbtW5kXbtm3bp0+fiIiIrVu3pqen1/Rv8BSC8IqkSLl9rDb9SyM2U/9ogNqrCaOGANzvJucLujcLmzdvPn/+/IKCgtatW4eFhbkeHDZs2MGDB/fs2RMWFhYREeF6sE2bNvv27du9e3fr1q0bNWp09ZjqW2+99cYbb1y+fNm16HplDBw4cPLkyQcPHuzSpYvrWcHBwSdPnnSt0L569ep69eq5tpw+fXrNT1J0F4LwX5oGiHmD1PRCc9Ry/aEOysuxqp/djxwDcLMZW4xPDprLhmuJYdf5td2ziVw6XBu6VP/FJvHHhJoem/L397/jjjuueTAwMDA2NvaaB4OCguLj4123HY5/zaT38/O7egJiJTVu3PinE2SklFfPQQwJCbn6+NWTC70BQXit1Cilf5jyWK7RK1P/eIAa05jWEIDbPNpRSWmt9G1+wy+W+KZyxYgbnhdxc5GRka7eq1u3bp7vt0aOHFnJ8+6vERERUfmTNGoDQXgdzeuJBUPU9EJz+DJ9cgfllVjVQWsIwB06NJAdGtxim/im1fz9vWPHDteNmJiYmJiY6u2k2u6+++7qPXHdunXuraSq+IK/odQopWCs4+tzzsQsfc95JpQCQN1EEN5MeKBYNFR7pKMyYLE+a4d5kyXEAAA+iiC8BSnE1E7KlhRt+VGzX5a+/wJhCAB1CkFYKW2CZU6S9lAHJTFLn7XDNElDAKgrmCxTWa7WcHCEnLzWWPyd+df+arsQJpQC9tW2bdtXXnll7ty5VhdiCz/88ENiYmIt7ZwgrJq2wXJ1kvbHr82ERfqrvdQpnap0jVwAdceUKVN69epldRU2EhkZWUt7JgirTJHiyW7KkJbyobXG/ELzw/5qy/qkIWA7iqL8+5np8EWMEVZTl4Zyw2jtrhZKzwX6nL1uukouAMDjCMLq0xQxI1pZlaTN3msmLdePFzOFBgB8D0FYU91C5aZkrX+YEkNrCAA+iCB0A4ciZkQrK0Zo7+4203KMM6W3fgoAwEsQhG4T3UhuStHahogeGRULDtMaAoBvIAjdKUAVM+PUjMHar/PMtBzjbJnVBQEAboUgdL/ezeT2sVrbENFjvr7oW1pDAPBqBGGtqKeJmXHq3EHqM5vNtBzjHK0hAHgrgrAW9W0ud4zTwgNFz0w95zgnVwCANyIIa1egJt7uo85OVB9eZ0zLNS5VWF0QAODnCEJPGBohd47ThBA9MvQ139MaAoAXIQg9pIGfmJ2o/ilBnbjGmJZrXNatLggAIIQgCD1sRKTcNV4TQkRn6OtP0BoCgPUIQk9r6CdmJ6oz45S0HP25PKPMsLogALA3gtAad0cpO8Y5DlwQsZl6/hlaQwCwDEFomWb1xPzB6os9laTl+nN5Rjln3gOAFQhCi6VGKQXjHHvPi7hMveAsrSEAeBpBaL2weiJziPo/McqwZfpL2wyDNAQADyIIvUVqlJKXon150tk3S997njAEAA8hCL1IqyC5YoT2cAel/2J91g6T1hAAPIAg9C5SiKmdlC0p2vKjZv/F+oELhCEA1C6C0Bu1CZY5Sdqk9kpClj5rh2mShgBQawhCL+VqDTcla0u+M4ct049cIgwBoFYQhF7tthC5JkkbFan0ytTn7KUzBAD3Iwi9nSLFk92U1UnaB3vNEcv0o5dJQwBwJ4LQN3QNlRuTtTvDlZ4L9Dl7uQgNALgNQegzNEXMiFZyRmrv7zGTluvHi2kNAcANCEIf072R3Jyi9Q9TYhbonxykNQSAmiIIfY9DETOilayh2msFZlqOcabU6oIAwJcRhL4qvqncPlZrGyJ6ZFRkfktrCADVRBD6sABVzIxT5w/WZmwx03KMH8qsLggAfBBB6PP6NJMFY7W2IaL7fD3rCK0hAFQNQVgX1NPEzDj187vUpzeZD64xiiqsLggAfAdBWHckhsn8sVo9TfTI0Fcd5+QKAKgUgrBOCXGI2Ynq+33Vh9Ya03KNS7SGAHArBGEdNKyl3DVeE0JEZ+hrv6c1BICbIQjrpgZ+Ynai+scE9YE1xrRco1i3uiAA8FYEYV02MlLuGq+V6CIuU887TWsIANdBENZxDf3E3waqL8UqySv05/KMMsPqggDAyxCEtpAapewY59h/QfTK1LedoTUEgH8hCO2iWT2RMVh9oacycrn+XJ5Rzpn3ACCEIAjtJjVKKRjn2HNexGXqBWdpDQGAILSfsHpi4RD1me7K0KX6rB2mQRoCsDeC0KYebK/kjdGyj5mJWfq+C4QhAPsiCO2rdZDMHqlN7qD0y9Jn7TBN0hCALRGEtiaFmNpJ2ZyiLTtq9l+sH6A1BGA/BCFEVLBclaQ92F5JoDUEYD8EIYT4sTXclKwt/s4cvkz/7jJhCMAuCEL8y20hcvVIbVALJXaBPmcvZxoCsAWCED+jKWJGtLI6SZuz1xyxTD9GawigriMIcR1dQ+WmZG1guNIzk9YQQB1HEOL6XK1h9gjtvT3mqOX698VWFwQAtYMgxM30aCS3pGi9msrYzIp3dxtpOddepPSHMjFyuc56hwB8F0GIW3Ao4qWe6oLB2h93Ow8VOUct16+u5XS+XCSv0B/vrAZqlpYIADVAEKJS7mgmt4/VBkfIvNMiMauizBDny8Wo5fpz0eqoVtLq6gCg+vglj8oKUMXMOHVES2fqSqPLQiUi0PlCT40UBODr6AhRNQPCZf5YrdiQX50TbYKtrgYAaowgRNWcLxf3rdY/7GPeGyViF+jphZxcAcC3cWgUVXB1XHBwMz25jYwIkpPWGFtOOWfFqwqHSAH4JoIQlfXT2TGlpUII8UovNdRfvl5g7DkvPrtLDXFYXSIAVB2HRlEFv427do7o092Vvw7QmtcT8Zks8AvAJxGEqKyGfqJf2HUOgI5uJT/srz7dXRm4WM8+RhYC8DEEIdxjaidl3mBt8jpj1g6mzwDwJQQh3KZvc7k5Wc04bN632uCiawB8BUEId4qoL9eO0gJU0TdLP1zEYVIAPoAghJsFqOKv/dXpnZU+i/RVx8lCAN6OIEStmNpJ+XyQ9sAanSFDAF6OIERt6R8mN6do6YXmxDVGCUOGALwVQYhaFFlfrknSSg3Rb7F+5BKHSQF4I4IQtSvIIeYOUlOjlL5ZxqZTZCEAr0MQotZJIWZEK3MS1THZ+of7GDIE4F0IQnjIiEi5fpT25i5zWq5RQRoC8BoEITynfQO5KUU7VSLu+kI/WWJ1NQAghCAI4WHBDpExRB0VqcQv1PNOM2QIwHoEITzNNWT4ToIyaoX+8QEOkgKwGEEIa4xupawbpf22wJyWa+ikIQDrEISwTMcGckOyVljkTFqunyuzuhoAdkUQwkqN/MXS4VpMExm/UP/qHEOGACxAEMJiqhQz49RXYpVBX+jzCzlICsDTNKsLAIQQ4t7blE4N5biVRt4Z5+u9VEVaXRAA26AjhLe4vbHckqLlnXaOXqGfL7e6GgC2QRDCizQJEMuHa90byTsW6nvOM2QIwBMIQngXTREz49TnY5SBS/SF3zJkCKDWMUYIb/RAO6VzQzl+pbHxlPO3cSojhgBqDx0hvFRsE7khWV37vfOeHOMy6/oCqDUEIbxXi0C5JkkL8RMJi/RDRQwZAqgVBCG8mr8q/txPfaKL0neRvvIYWQjA/QhC+ICpnZT0QdqktcasHUyfAeBmBCF8Q2KY3Jyizj9sTlhtFDNkCMB9CEL4jJb15bpRmp8iErP0by9xmBSAexCE8CUBqvhogPpYZ6XPIn3192QhADcgCOF7pnZSPhqg3ZOjM2QIoOYIQvikoRFyY7L2yUFzWq5RThoCqAGCEL7qthC5MVk7WybuXKJ/X2x1NQB8FkEIHxbkEOmD1ORWSvxCffMphgwBVAdBCN8mhZgRrcxOVFOy9b/s5yApgCojCFEXjIyU60Zp/7fTnJZrVJCGAKqCIEQd0aGB3JyinSwRg77QT5ZYXQ0A30EQou4IdogFQ9SkSCV+ob71DEOGACqFIESd4hoy/L87lOFL9b8d4CApgFtjYV7UQalRSscGcmy2kX/G+WZvlYV9AdwEHSHqph6NZN4Ybfd5Z9Jy/VyZ1dUA8GIEIeqsRv5i2XDt9sYyfqH+9TmGDAFcH0GIukyVYmac+nKsctcXesZhhgwBXAdjhKj77rtN6dRQjltpbDntfL2XqjBkCOAn6AhhCzGNZV6KtuWUM3mFfqHc6moAeBOCEHbRJECsGKF1ayTjF+p7zjNkCOAKghA2oiliZpz6THel/2J90bcMGQIQgjFC2NDUTkq3UJmaY2w763yxJycZAnZHRwg7Smgu88aoy44671tlXNatrgaApQhC2FSLQLk2SauviYRFemERQ4aAfRGEsC9/VXzYX32ii5KwSM85ThYCNkUQwu6mdlLmDtIeXGPM2sH0GcCOCEJA9AuTm1PUeYXm/auNEoYMAZshCAEhhGhZX64frWmK6JulH7nEYVLARghC4IoAVXw8QJ3UXolfqK/5niwE7IIgBH7myW7KxwO0tBz9/33NkCFgCwQhcK1hLeWGZG3OXnNarlFOGgJ1HUEIXEe7ELkpWTtTKu5aop8osboaALWJIASuL8gh5g1WR7dS4jP1LacZMgTqLIIQuCEpxIxo5b1EdfQK/a/7OUgK1E0EIXALSZFy/SjtdzvNablGBWkI1DkEIXBrHRrITSnaiRIx+Av9FEOGQN1CEAKVEuIQmUPUO1vIuIV6/hmGDIG6gyAEKksK8VJP9Y14ZdhS/e8HOUgK1BEszAtUzT1tlc4N5dhs4+tzztd6qSzsC/g6OkKgyno0klvGaFtPO0ct18+XW10NgJohCIHqaOwvlg3XohvL+IX61+cYMgR8GEEIVJOmiJlx6ks9lbu+0BccZsgQ8FWMEQI1MuE2pXNDOW6lsfm08/VeqsKQIeBr6AiBmoppLLekaJtPOVOy9QsMGQK+hiAE3KBpgMgeobUNlncs1PeeZ8gQ8CUEIeAemiLe7qM+3V3pt1jPOsKQIeAzGCME3GlqJ6VrqLxnlbH7nPhVNCOGgA+gIwTcrG9zuTlZzThsTlhtFOtWVwPgVghCwP0i6st1o7R6qkhYpB8uYsgQ8GoEIVAr/FXxl/7q412UPov0VcfJQsB7EYRALZraSfl8kDZxjTFrB9NnAC9FEAK1q3+Y3JSipheaD6wxShgyBLwPQQjUusj6Mne0JoVIXKwfucRhUsC7EISAJwSo4m8D1bQopW+WsfEUWQh4EYIQ8BApxIxo5YN+6ths/c/7GDIEvAVBCHjU8JYyd7T2h13mtFyjnDQEvABBCHhauxC5KUU7XSoGfaGfKLG6GsD2CELAAsEOMX+wOipSic/U804zZAhYiSAErOEaMny3rzpqhf7Rfg6SApYhCAErjWol143SZu00p+UaFaQhYAWCELBYxwbyy9Ha4SLnkKX66VKrqwHshyAErNfIX3wxXOvdTPZZpO/6gSFDwKMIQsArqFLMjFNfjVUGL9XnFXKQFPAcFuYFvMi9tymdQ+XYbGPrGefrvVQW9gU8gI4Q8C7RjWTeGC3vtHPUCv18udXVADZAEAJep7G/WD5c69FIxi/Ud59nyBCoXQQh4I00RcyMU1+MUe5comd+y5AhUIsYIwS81/3tlM4N5biVxqZTDBkCtYWOEPBqPZvIDcnquu+d3TP0P+25tjXMPuZ8cqNhSWFAnUEQAt6uRaBcnaT1bip/k2e+sv1fsZd9zPlCvvF8jGphbUAdwKFRwAf4q+LD/mpsU/PpTcbpUvFatyspmDVUaxJgdXGAjyMIAZ/xeGelY0M5apmef9LfVEhBwD04NAr4kkHh8r2+av45ea5cXNY5swJwA4IQ8CXZx5yz95rbR5apUnSdp8/aYZqkIVAzBCHgM66OC7YMdG4fq/Vurvx5nzlgiX7gAmEIVB9BCPiGa2bH+Cnii2Fq51DZqr5MyKI1BKqPIAR8QwM/sXjYz2bH+Cli7l3qo52UzSnasqNmv8X6PlpDoOoIQsA3xDeVjf2vfTBAFXeGy7bBclWSNqm90i9Ln7XDNEhDoCo8EYRLlizp2bNnTExMTExMbm6uB14RsBspxNROyuYUbflRs1+WvpdLdQOV5okg7NSp0/r167dv3/7OO+/853/+pwdeEbCnqGCZk6Q91EHpv5jWEKgsTwThbbfdVr9+fSFEcHCwqnI5KKAWuVrDLSla9jEzMUvfQ2sI3Io7ryxz9uzZjRs3/vSRRo0aJSQkuG6XlpY+8cQTL7/8shtfEcB1tQmW2SO1D/aaiVn6lE7KK7Gqg/kAwA24MwgvXry4devWnz7SqlUrVxCWl5enpqZOmDBh5MiRbnxFADfiag2HtZRT1htxmfpfB6gxjVnGCbiOSgXhuXPn5s6dm5+fX1RU9I9//OPq4+Xl5a+88sqKFSuaNm36/PPP33HHHS+99NK/P90wjAcffLBfv37Tp093V90AKqN1kFwxQksvNIcv0yd3UF6OVf1oDYGfq9T/E4cOHVq1alW9evXmzZv308dffvnllStXvv/++6NHjx4xYsTZs2ev+/Q333wzOzt769ataWlpkyZNckPVAKoiNUopGOvYc17EZerbzjBqCPyMdDor+3/Frl27evbsWVFR4bqr63p4ePiCBQsSExOFEEOHDh05cuRTTz317088d+7cuXPnXLcVRWnTps119/+LX/xi3759rr0JITp27JicnHyTeoqKioKDgytZPNyrtLRU0zRNY/USa1T7wz//W/HUZvFgO/Hi7YLWsHqcTufly5eDgoKsLsSmiouL/f39Kz/vUtM0KW8xKFD9L7Ljx4+fOXMmPj7edTc+Pn7nzp3X3TI0NDQ0NLQy+ywuLr4amUVFRaZ57XrcP2Wa5s03QO0xf2R1ITZV7Td/bKRIaCr/c5O8I0vMSTBjG7u9tLrP6XTy4beQ682/ZbZVSfWD8NSpU4GBgX5+fq67oaGh27dvr0kpqqqOHj36mWeeqeT25eXl/v7/dqUNeITT6aQjtFBNPvyt/MXCYSK90By3WjzUQfnfnqo/5zRVhdPp1HWdLx+rGIZRpY6wMqp/cCQkJKS0tNQwDNfdS5cuNWzY0E1VAahdqVHKjnGOAxdEbKaed5pRQ9ha9YMwIiJCVdVDhw657h48eLB169ZuqgpArWteT8wfrL7YU0leoT+XZ5QaVhcEWKT6QVi/fv2UlJR33nlHCHH48OHFixdPmDDBfYUB8ARXa/jNRdF9vr7uBK0h7KhSYzynTp1q3ry567aUMiIi4ujRo0KI3/3ud8nJyVFRUefPn//Vr37VrVu3WqwUQO1oVk+kD1LTC820HD2ltfJmb7U+g7+wk0p93ps1a3bdsyzatGmzc+fOY8eONWzY0HU1UQA+KjVKGdxCeS7PiM7Q/9xPHRjOZWhgF244kygiIoIUBOqAUH8xO1F9q7c6cY0xLde4VGF1QWf688sAABHvSURBVIBHcEotgJ8Z1UruGq8JIaIz9NXfM2qIuo8gBHCthn5idqL6xwR1Eq0hbIAgBHB9IyOvtIY9MvSc47SGqLMIQgA31MBPzE5U3+2rPrzOmJZrFNEaoi4iCAHcwvCWcue4K61h9jFaQ9Q1BCGAW3O1hu/3VR9db6TlGOfKrC4IcB+CEEBlDWspd43XwgNF9wx90bcsv4A6giAEUAUhDvF2H/Ufd6rPbDbTcowfaA3h+whCAFXWL0zuGKe1DRHd5+uZtIbwcQQhgOoI1MTMOPXzu9QZW8y0HOMsrSF8FkEIoPoSw2TBWK1tiOgxX884TGsIn0QQAqiRepqYGaemD1J/k2em5RhnSq0uCKgighCAGyQ0l9tdrWFGxbxCWkP4EoIQgHu4WsP5g7Xn8820HOM0rSF8BEEIwJ36NLvSGnaZVzFnL60hfABBCMDNAlQxM05dPFR7+ytz9Ar9eDFXZYNXIwgB1Io7msmCcVpicyVmgU5rCG9GEAKoLQ5FzIhWVozQ3t9jJi3Xj16mNYQ3IggB1K7oRnJzitY/TInNpDWENyIIAdQ6V2u4coQ2Z685Ypn+Ha0hvAlBCMBDujeSm5K1geFK7AJ9zl6TMISXIAgBeI6miBnRyqok7YO95ohl+pFLpCGsRxAC8LRuoXJjsnZnuNIrU3/7K5PeENYiCAFYwNUark7SPvvGHLhEP3iRMIRlCEIAlukaKr8crY1vo/RZpM/aQWsIaxCEAKykKeLJbsqmZO2L78wBS/T9FwhDeBpBCMB6t4XI1UnaxHZK3yxaQ3gaQQjAKyhSTO2kbE7Rlh01+y3W99EawlMIQgBepG2wXJWkTWqv9MvSZ+0wDdIQtY8gBOBdpLjSGq44ZvbL0veeJwxRuwhCAN4oKliuHKk91EHpv5jWELWLIATgpVyt4ZYUbeVxs2+WvpvWELWDIATg1doEyxUjtIc7KP2z9OfyjAqWr4C7EYQAvJ2rNcwfq20744zL1LefpTWEOxGEAHxD6yC5YoT23zHK8GX6c3lGOa0h3IQgBOBLUqOUgrGOvedFXKaef4bWEG5AEALwMeGBInOI+j8xStJy/bk8o8ywuiD4OIIQgE9KjVIKxjn2XxC9MvWttIaoAYIQgK8KqycyBqsv9FRG0xqiBghCAL7N1RoevChiM/Utp2kNUWUEIQCf17yemDdIfbGnkrJCf3KjUaxbXRB8CkEIoI5IjVJ2jHMcLxbRGfq6E7SGqCyCEEDd0ayeSB+kvtlbmbDamJZrXKY1RCUQhADqmtGtlF3jNCFEdIa+5ntaQ9wCQQigDgr1F7MT1bf7qBPXGNNyjUsVVhcEL0YQAqizkiLlrvFXWsPVtIa4AYIQQF3W0E/MTlT/mKBOojXEDRCEAOq+kT+2hj0y9JzjtIb4GYIQgC008BOzE9V3+6oPrzOm5RpFtIb4EUEIwEaGt5Q7x11pDbOP0RpCCIIQgN24WsPZieqj6420HOOHMqsLgtUIQgB2NDRC7hqvhQeKHhn6om9Z5NfWCEIANhXiEG/3Uf9xp/rLLSatoZ0RhABsrV+YLBirtQ0R3efrmbSGtkQQArC7QE3MjFPnDlJnbDHTcoyztIY2QxACgBBC9G1+pTXsMV/POExraCMEIQBcUU8TM+PU9EHqb/LMtBzjTKnVBcEjCEIA+JmE5nK7qzXMqJhXSGtY9xGEAHAtV2uYMVh7Id8cvUI/Xuy87ikWBWed317irHyfRxACwPX1bia3jdW6hsroDP3dPeYvNhk//deCs86H1xkGOej7CEIAuKEAVcyMU7NHaCdLRNYRc8r6K1noSsF5g9W2wdLaClFzBCEA3MLtjeWWFG1KR/WzQ+aQpcbOc5IUrEsIQgC4NYciZkQra5K0/DPOISv9Xo1VSME6gyAEgMpySNEmSPZoZKatMtJyjN3nGSGsCwhCAKiUH8cFlexBFdM6K6dKnYO/0NNyjP0XiEPfRhACwK1dMzvm93eovZvJcW2U2CYyMUt/cI3xzUXi0FcRhABwC04h/nebmfHz2TG/jVM1RdzRTH5zj6NrqOy9SJ+Waxy7TBz6HoIQAG5BCrFgiNrm57NjpBBv9VYHhstgh5gRrexNdYT6i27z9Wm5xokSqypFdRCEAOAGjf3FzDh1f5oj1F90m1fxXJ5xjlUsfARBCABu0zRAzIxTt4/TSnTRfm7Fc3nGhXKra8KtEIQA4GaR9eXbfdT8sdq5MtEhveKlbUZRhdU14cYIQgCoFa2D5OxEdU2SduiiaD+3YtYOs0S3uiZcD0EIALWoc0P5t4Fqzkgt/4yzY7r+9ldmmXHrZ8GTCEIAqHVdQ+XcQeq8werK42bHdH3OXlNnoUOvQRACgIfEN5VZQ7XP7lQ/P2R2SNfn7DVZxckbEIQA4FEJzWXOSO1vA9RPD5o95uvphSZpaC2CEAAskBgm147S3u6jztphRmfo6YUcKrUMQQgAlhkcIfPGaK/1Ul4vMPss0lceozm0gGZ1AQBga1KI0a2UpEhl/mHziQ1G0wDxai91YDiLHXoOHSEAWE+RIjVK2XO39mQ35dH1xpCl+tYzdIceQhACgLe4GoepUcqYbGPIUn37WeKw1hGEAOBdHIqY2kk5dI+WGqWMWm6k5Rj7WPu3NhGEAOCN/BQxtZOyP03r21wOXKyn5RgHWfu3dhCEAOC96mviyW7KgTRHbBOZsEh/cI1RWEQcuhlBCADeLujHtX/bhoj4hfq0XOP7YqtrqkMIQgDwDY38xUs91T13O0L9Rff5FU9uNE6WWF1TnUAQAoAvaRIgZsapBeM0IUTXeRXP5RnnWfu3ZghCAPA9LevLt/uo28Zq58pEx/SKl7YZF1n7t7oIQgDwVa2C5OxEdVOy9n2x6MDav9VFEAKAb4sKlrMT1VVJWv4ZZ9TnFbN2mKWs/VsVBCEA1AVdGsq5g9TsEVr+GSdr/1YJQQgAdUf3RnLuIPXzu9SsI2Z71v6tHIIQAOqa3s1k1lDt7wPUf3xjdp+v/+2AyeK/N0EQAkDdlBgmVydp/6+P+s7uK2v/kobXRRACQF02OEJuSdH+0Fv9bYHZZ5GedYSRw2sRhABQ9w2OkPljtWe6K89uNvtm6auO0xz+C0EIALYghUiNUnbfrT3VTZmWayRm6etOEIdCEIQAYCuutX93361N7aRMXmsMWapvO2P3OCQIAcB2HIp4sL2yJ1VLjVKSs43RK/QdP9g3DglCALCpK2v/pmqDWygjlulpOcaBC3aMQ4IQAGwtUBNPdlP2pzpim8i+WfqDa4xDNlv7lyAEAFxZ+/ebexxdQ+UdC/VpucbxYrvEIUEIALgi2CFmRCt7Ux2h/qLrPH1arnHCBmv/EoQAgJ9p7C9mxqn70xyh/qLbvIrn8oxzZVbXVJsIQgDAdTQNEDPj1O3jtBJdtJ9b8VyecaHc6ppqB0EIALihyPry7T7q1jHauTLRIb3ipW1GUYXVNbkbQQgAuIU2wXJ2oromSTt0UbSfWzFrh1miW12T+xCEAIBK6dxQ/m2gunLklbV/3/7KLDOsrskdCEIAQBV0C5VzB6nzBqsrj5sd0/U5e03dxxe0IAgBAFUW31RmDdU+u1P9/JDZPUP/2wHT8NnTDglCAEA1JTSXOSO1DxLVD/f58Nq/BCEAoEYSw+TaUdpbvdVZO8zbM/T0Qh87VEoQAgDcYHCEzBujvdpLeb3ATFik5/jO2r+a1QUAAOoIKcToVkpSpDL/sPn4l0bTAPFaL3VAuLS6rlugIwQAuJNr7d89d2tPdlMeWW8MWapv9e61fwlCAID7XY3D1ChlTLYxZKm+/ayXxiFBCACoLQ5FTO2kHLpHS41SRi030nKMfd639i9BCACoXX6KmNpJ2Z+mxTaRfRfpaTnGwYteFIcEIQDAE+prYka0cvheR2wTmbBIf3CNUVjkFXFIEAIAPCfox7V/24aI+IX6tFzj+2KLSyIIAQCe1shfvNRT3XO3I9RfdJ9f8eRG42SJZcUQhAAAazQJEDPj1IJxmhCi67yK5/KM81as/UsQAgCs1LK+fLuPum2sdq5MdEyveGmbcdGza/8ShAAA67UKkrMT1U3J2vfFooNn1/4lCAEA3iIqWM5OVFclaflnnFGfV8zaYZYaYs9552sF117Iu0QXT2wwyt1xfW+CEADgXbo0lHMHqdkjtPwzzo7p+roTzu8uOZ/LM65uUGaIe1YZ8U2lnztCjCAEAHij7o3k3EHq53ep8wrN5cecm045Z+QZQogyQ6TmGOOj5KT27okwVp8AAHiv3s1k9ggt94Tz+Xzjgz1mwWnVoTpT2yruSkFBRwgA8H6JYXJ1kvZhP23zGcWhON2YgoIgBAD4hDJD/PWA+UasEVZP/nS8sOYIQgCAt7s6LvhAlPnHPvJiuXBjFhKEAACvds3sGCnEn/qqbsxCghAA4NWOFTvvu+1nc0SlEO8kqE0DpK3PIywtLf3++++trsK+zp49e/HiRaursK8jR46Ypju+AFB1xcXFJ0+etLoKe2kbLO+77UpanTlzpqioSAihSPFMd8XW5xHm5+dPmjTJ6irs69VXX/3000+trsK++vXrd/78eaursKl169ZNnz7d6irs6ze/+c2CBQvcu09fDUIAANyCIAQA2BpBCACwNel0Oq2u4YpPP/302WefrVevXmU2rqiouHz5csOGDWu7KlxXUVGRpmmV/I8Ftzt79myjRo2klFYXYkfl5eUlJSUNGjSwuhCbKioqcjgcAQEBldx+8eLFnTt3vvk2XhSEpmkePnzY6ioAAHVHy5Yt/fz8br6NFwUhAACexxghAMDWCEIAgK0RhAAAWyMIAQC25nsr1BuGsWfPnoKCAsMwuMqa5x07dmzJkiX79u1r2rTp/fffHxkZaXVFNuJ0OufPn79jx47Lly936dJlwoQJgYGBVhdlR9u2bdu6deu9994bEhJidS02snHjxl27dl29+8gjj6iq6pY9+15HuGTJkuHDh//pT396+umnra7Fjh566KHc3NyWLVsWFhZ26dLl66+/troiGzFN85NPPgkMDGzTps3f//73gQMHGoY7lydFZVy8eHHixInTpk07ffq01bXYS3p6+scff3zoR2485cH3Tp8wTVNRlA0bNowePfrs2bNWl2M7paWlV09lHT9+fLt27WbNmmVtSfZUXFwcEhLy1VdfderUyepa7GXq1Km33377E088cfDgwdtuu83qcmzk6aefDgoKevnll92+Z9/rCBXF92quS356QYfS0tKgoCALi7GzDRs2NGjQICIiwupC7GX16tUHDhx4+OGHrS7EprZv3z5r1qx//vOfpaWlbtwtoYJq+uKLLzZv3jxlyhSrC7Gd8ePHh4WFjRs3bu7cucHBwVaXYyPFxcX/8R//8d5773FxO0uEh4dHRERcvHjx7bffjo6OPnfunLv27HuTZeANNm/ePGnSpE8++SQsLMzqWmzno48+unjx4rx58+69996dO3eGh4dbXZFd/PrXv544cWKnTp3KysqsrsWOnn32WdcN0zT79ev3zjvvPP/8827ZMx0hqmz79u0pKSkfffTR8OHDra7FjoKDgyMiIp588snIyMjs7Gyry7GRTz/99J///GevXr0SEhKEEGPGjFmyZInVRdmRoih9+vQ5dOiQu3ZIR4iq2blz58iRI2fPnp2UlGR1LbZTWlrq7+/vOi53+vTpw4cPt27d2uqibCQnJ0fXdSFEeXl5QkLCq6++2rt3b6uLspGSkhLXijelpaXZ2dkPPvigu/bse7NGjx49OmbMmEuXLn3zzTfR0dFRUVHp6elWF2UjvXr1OnjwYLt27Vx3hwwZ8tvf/tbakuwjOzv7sccei42NlVKuXr06OTn5gw8+YLzK88rKygICApg16mFRUVFdu3Zt2LDh+vXro6Kili5d6q6V4HwvCMvKyr766qurdwMCArp27WphPXaze/fukpKSq3dDQ0Pbtm1rYT12s2fPnt27d0spu3Xr1qFDB6vLsSmn07lt27Zu3br5+/tbXYuNHD16dNu2bcXFxe3atevVq5cb9+x7QQgAgBsxWQYAYGsEIQDA1ghCAICtEYQAAFsjCAEAtkYQAgBsjSAEANgaQQgAsDWCEABgawQhAMDWCEIAgK39f3j8K0UXvmeuAAAAAElFTkSuQmCC", "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" ], "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" ] }, "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.7.3" }, "kernelspec": { "name": "julia-1.7", "display_name": "Julia 1.7.3", "language": "julia" } }, "nbformat": 4 }