{ "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 Input and output formats 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.846541914835 -0.70 0.80 4.6\n", " 2 -7.852232192886 -2.24 -1.54 0.80 1.0\n", " 3 -7.852621531955 -3.41 -2.53 0.80 3.2\n", " 4 -7.852621949638 -6.38 -3.38 0.80 2.2\n", " 5 -7.852621958664 -8.04 -4.43 0.80 1.6\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+gvaeTAAAgAElEQVR4nO3deWAU5f348WdmnpCT+8oB4YwgBEIOjnBoQBAQcgGilZJ6Qq3V+KO1ag+NVlvxatN6fIFaFUWtAkmEhCMgCIFwhjuIRFBOAQGBkHNm9vfHKioGyLHJ7Gber7+SZbP7oCHvPJ+dmVUcDocAAMCuVKsXAACAlQghAMDWCCEAwNYIIQDA1gghAMDWCCEAwNYIIQDA1gghAMDWCCEAwNYaKIR79+598803d+3a1TBPBwBANTVECNetW5eSkuJwOKZNm7Zy5cor3S0zM3PNmjXVf9jKykpXrA61YRgGF+ezEN/81tJ13eol2Jeu6y7/4dMQIfzHP/7x7LPP3n333S+99NJLL710pbvl5eVt3ry5+g9bVlbmitWhNiorKw3DsHoV9sU3v4UcDgf//S1UUVFhmqZrH7MhQlhYWNinTx8hRJ8+fQoLCxvgGQEAqKaGCGFZWZmXl5cQwtvbu7S0tAGeEQCAaqprCB0Ox9mzZysqKi67fd++fdu2bXNOz4KDg48fPy6EOHr0aEhISB2fEQAAF6pTCCdNmtSyZctWrVplZGRcurGysjIxMXHMmDH33HNP3759v/766wkTJsyZM8c0zdmzZycnJ9d5zQAAuEydQnjvvffu2LEjIiLixzd+9NFHRUVFhYWFBQUFkZGRzz333EMPPeTn5zd8+PDy8vJHH320Ls94qqqXqC/qooRjuAAAtaLU/TjUfv36Pf7447fddpvz06SkpAEDBvzxj38UQqxateqOO+5wzkWv6eabb966dWuzZs2cn0opc3Nz27Rp8+P7TF7jNamTMbmTKYQoLi4OCAgoMZTJa+TvrjeGB7r4OCJcRVlZmZRSSmn1QmzK+c1v9SpsyuFwlJSU+Pv7W70QmyopKfH29tY0rZr39/HxueZPKtf/IDt06NClKHbp0uXEiRNlZWU+Pj7X/MJevXrFxMTce++9361MytDQ0MvuM/9mkZyre3mrU7urDodD9QmYkqvfd70a393btX8LXJ38ntULsSmHw0EIreJwOBRF4b+/VVRVrVEIq8P1P8hKS0u9vb/Lko+Pj8PhKC0trU4IpZStW7fu2rXrVe7jJ0XGKJmcqwshRrdRpuTqd4apU7pzoTgAQC25PiHt27c/c+aM8+PTp0/7+Pi0aNHChY/vbOGcz8yYnCZ3X0cFAQB14vqKxMTE5OXlOT9et25dTEyMoiiufQqHEA6HOFepPLLJWHmMC30BAGqvTqPRnJycI0eOnD59euXKlefOnRs/fnxwcPCvf/3r6Ojo/v37h4SEPPnkk6+++qqr1upUooukXP3X16uj2pSMWeX9i1X6wLbKK4O1TgEuzi0AXMXHH3/87rvvWr0KGxk7duxdd91VH49cpxAeP378wIEDU6ZMEUIcOHDAef297t27L1++/NVXX7148eK//vWvCRMmuGalQgghSnSR+P3rghcuOPLivRKX616qiMrQH+ytPh6hebvyBVQAuKJ169Z5e3snJCRYvRBbyMvLW7FihTuG8J577qny9oEDBw4cOLAuj3wlt32i39dDndz1u4munxRZN8vEXP2VwdqSw47wBXp6rHZLR7aGABpCnz59br31VqtXYQuVlZXZ2dn19OAedvj7WzfK1j89UcJPiqxRUhHiF93EymOOB9cbr+8V/4rVujQlhwCAa/OwQy5bV3W6oJ8UvlIIIW4KVnZMkCOD1YFZelqBUcbbBAEArsXDQnhNXqpIDVe3JcsD50X4An3xIY4pBQBcTWMLoVOIvzI3TpszTHtsszFqif7Zt+QQAFC1xhlCp+FByrZkOb6jOnSRnppvFFdavSAAEOJ8pTCu/Mu56RDnLn9fO9SvxhxC8f2kdOdEebZcXD9fn7ufC3MDsNjEFXrKaqPKFpoOcd9aY9wy3k+nQTXyEDoF+ylz47R347QXdpo35eiFTEoBWGdmf23pEfOOVYb+09/MHUI8sN743wHzmZhanhAdFBR04cIFFyxRCF3XV6xYUfe3JyouLm7Tps3VH6dXr16HDh2q4xPVhS1C6HRjkLItWSaEqjcs0lPzjQtMSgFYIaqNkjtWrjhqTln9QwsdQvxmnfHOfnPxaBkXVMuzv8rLy+ueLqeSkpJRo0YZRl0Pvm/SpMkTTzxx9QttunDZtWOjEAohpCpSw9VdE72YlAKw0GUtdEkFL+NwOObNm/fYY4+98cYbFRU/vOq4bNmyJ5988vXXXz937pzzlpUrV+7du9f58aFDhz7++GMhRGZmphDiP//5z+zZs0+ePFnNJy0uLn755Zf//Oc/5+TkOG9RFOXSuw8tXLjwyJEjb7311hNPPLF27VpX/C1dw14hdAryE3PjtPeGay/tMkdk63vOMikF0NB+3MJf57m4gkKI6dOnv/XWW927d58/f/7EiROdNz7xxBMzZswIDAzcunVrZGSks4Vz5sz59NNPnXfYuXPn888/X7tnLCkpGTp0aFFRUceOHR9++OGXX35ZCFFaWjp9+nTTNIUQzz777O23337gwIGAgICEhIR169a54O/pCh52ZRkXuiFQ2Zok//u5OTJHn9xV/WuM1szL6jUB8GSZX5mvFtZsztStmTL/oCmEiGytPLvdeHZ7Db723h7qbV2r3swUFRUtXLjw0KFDfn5+U6dO7dSp06ZNm8LCwl544YVdu3Z1795dCDF69OhZs2b94Q9/qPIRkpKShBD33ntv9d9/e968ecHBwa+99poQon///iNGjHjooYcuu098fPyjjz4qhDh58mRWVtaQIUOq+eD1yr4hFEJIVUzrqU7soj5dYPSar/8tRp0apnJlNgC1M7idGiBr8CPE4RAv7za9VIeqiBZNlEf6qDX6ARTe6or33rVrV58+ffz8/IQQ3t7ekZGRu3fvFkK0bNnSWUEhxLBhw5w31sixY8ecJfPx8ZkzZ86P/2jnzp2DBg1yfhwZGVleXn748OHWrVv/+D79+vVzfhASErJz586aPns9sXUInVp7i/RYbWqY44F1xhv7zFcGa32u/O0FAFfSzleMDKnuTw/n64JrvzaXjpHNvMSoJfobn4t5cZp0xQtWxcXFzrcDciotLfXz8/Pz8ystLb10Y0lJibOUmqZdOiimpKTk6o/crFkz56D15zvFixcvBgQEOD/Wdb2ystL5+D+mae74DkF2fI2wSjFtlPwEeU8P9eYlemq+wQmtAOrPZUfHVHkcaR1t3br1888/F0IcPHhw69atgwcPDgsL8/PzW7RokRDi4sWLCxYsiIuLE0J07ty5oKBACOFwOBYuXOj8cn9/fy8vr9OnT1/2sAEBAUlJSUlJSePHj//5k2ZlZTlb++GHH4aFhbVv3941f5l6Rgh/oCoiJUzdM8lLCNF7gT77M9PkMBoArlblMaIub2FERMSUKVMSEhIGDRr09NNPh4aGent7//e//50+ffrYsWPDw8MHDx48efJkIcR99923dOnS4cOHR0dHX9rnaZr2m9/8JiIiIiYmZs+ePdV80q5du8bGxo4bN27GjBmzZs1ywV+jQTAavVwrb5Eeq/0qzPHb9cZ/9pmvDtb6t2VSCsBl7ltr/O+AuXSMHBr4k58tUW2UJWPkzUv0O9eId+PqOkJs3779ggULtm/f3qlTp8DAQOeNo0ePLioq2rt3b2BgYEhIiPPGzp0779u3r7CwsFOnTq1atbo0U/3nP//5/PPPX7x4sWnTptV80ri4uLvuuquoqKhXr17Or2ratOmJEydUVRVCrFq1ytfX13nP+++/v+4nKboKIaxaVBtlXYJ8Z78Zv1y/OUR9eZDWxsfqNQFoFIa0V+4Mu7yCTgPaKsvHym2nXTOM8vb2/vl7pPv5+UVHR192Y0BAwIABA5wfe3n9cAB9kyZNmjRpUqMnbd269Y8PkFEUpV27ds6PmzVrdun2SycXugNCeEWKEClhanyomlZg9Jpf+ad+2oO9a3ZMFwD83F3XXe01qQFtlQG1nUJ17NjRufcKDw9v+P3WLbfccql5NRISElL9kzTqAyG8hpbeIj1Wu/M69bfrjXeKzFcGa4PaEUMA7mjHjh3ODyIjIyMjIxv42SdNmlS7L1yzZo1rV1JTHCxTLZGtlbx4+VBvNTlXT1ltnCq79pcAADwCIawu56R0761eLb1F7/mV6bvNq7yjGADAUxDCmmnRRKTHarm3yAVfmgOy9PyTxBAAPBshrI2IVsqa8TItSr39EyNltXGy9NpfAgBwTxwsU3vxoeqIYPWFnUafhZV/jNB+21vVOIwGsA1N0/71r399+OGHVi/EFs6cOTN06NB6enBCWCf+UqRFab/opj603nhrv/nvWK3Kc4MAND6PPPLIpbc3QgPo2LFjPT0yIXSBHs2VZWPlokPmL1cbNwQqLwzU2vtavSYA9axly5Y/PzMdnojXCF0mPlQtnCS7NhN9F1bO3GFWuOjKuQCAekUIXclPirQoLW+8XH3c7LtAzz3KMaUA4O4IoeuFNVeWjJEvDFSn5Rnxy/XDF8khALgvQlhf4kPVPRNldBul30I9rcAod5fLrAMAfoIQ1iPnpHRjotxyytF3ob7sCFtDAHA7hLDedW+mLB4tXxyo3r/OiF+uf1VMDgHAjRDCBnJpUhqVwaQUANwIIWw4vlKkRWnbJ8gD50X4Aj3nMFtDALAeIWxoHf2VuXHa60O032804pfrX14ghwBgJUJojZEhyo4JcmSwOiBLTyswypiUAoBFCKFlvFSRGq5uS/5uUrr4EFtDALAAIbRYiL8yN06bPVR7dJMRv1w/wKQUABoWIXQLI4KV7c5Jaaaemm8UV1q9IACwDULoLpyT0p0T5dlycf18fe5+LtoNAA2BELqXYD9lbpz2Tpz2wk5zZI6+91smpQBQvwihO4oLUrYly/hQddgiPTXfuMCkFADqDSF0U1IVqeHqroleTEoBoF4RQrcW5CfmxmnvDdde2mWOyNb3nGVSCgAuRgg9wA2BytYkeXs3dWSOnppvnGdSCgCuQwg9g1TFtJ7q7kleQohe8/W5+032hgDgEoTQk7T2FumxWuYo7dVCc3i2vusMNQSAuiKEniemjZKfIO++Tr15iZ6ab5yrsHpBAODJCKFHUhWREqbumeQlhOi9gEkpANQeIfRgrbxFeqz28Sjt9b3mwCx98ylqCAA1Rgg9XlQbZX2C/G0vNX65nrLa+KbM6gUBgEchhI2BIkRKmLp3kldLb9FrfmX6bpNRKQBUEyFsPFp6i/RYbdlY+eFBc0CWvvEkMQSAayOEjU1kayUvXj7UW03K1VNWG6eYlALAVRHCRui7SemtXi29Re/5lem7TYPNIQBcASFstFo0EemxWu4tcv6X5oAsPZ9JKQBUhRA2chGtlLXjZVqUevsnRspq42Sp1QsCADdDCG0hPlQtnCS7NhN9FjIpBYCfIIR24S9FWpS2ZrzMOWzGZOrrThBDABCCENpNj+bKsrHy6Wh1yiojZbVxgkkpANsjhHZ0aVLad2HlzB1mhWn1ggDAOoTQpvykSIvS8sbL1cfNiIV67lEmpQBsihDaWlhzZckY+fwAdVqeEb9cP3yRHAKwHUIIER+q7pkoo9so/RbqaQUGk1IAtkIIIcT3k9KNiXLzKUefBfqyI2wNAdgFIcQPujdTskfLFweq968z4pfrXxWTQwCNHyHE5S5NSqMy9LQCo9ywekEAUJ8IIargK0ValLY5SRaeFX0W6jmHv9saVlb18mGVNwKApyCEuKKuTZUPb9JeG6z9bqMRv1xf87V5U45+ruIn9zl4wTFssV6iW7REAKgzQohrGBmi7JwgRwark1YYof7KmKX6mfLv/uhQsWPSSuOVwZqftHSJAFAHhBDX5qWK1HB1W7JUFfFVsRi6SD9boRwqdiSvMGYN1WLaKFYvEABqj9/kUV0h/srcOG3ZEcfda/RemWqgr3h/hNavNRUE4NnYEaJmRndQVo2TUhVFF8TGU5xfAcDjEULUzKFixy9WGVkjzPSByu82GIm5xunya38VALgtQogauPS6YFQrx7QeYt5wufUbM3y+vuQwW0MAnooQorqOXnQk5RpvDPvh6JjETsprg7W2vmLaWmN6nsFJFAA8ESFEdbX1VT742dExCZ3U94ZruybJ4koxIEvffpqtIQAPQwhRXU1UcV3zKo4RDW+ptGgi5g3XnoxSxy7VZ+4wTWoIwHMQQrjMrV3UDYky57A5aol+hLc2BOAhCCFcqVOAsmqcTAhVozL0eUVchBSAByCEcDFVEanhas4Y+ddt5uSVxllOrgDg3ggh6kVMG2VbsgzyE5EZ+qfHGZMCcF+EEPXFV4r0WG3WUG3KaiM136hgUArALRFC1K/RHZRtyfLLYhGTqe86w9YQgNshhKh3bX1E1ijt933UETl6+m7OrQDgXgghGkhKmLopUX500ByzRD9WQg0BuAtCiIbTpamyepwcEaz2zzQyv+I1QwBugRCiQUlVPBqhLhypPbLRTFltFFdavSAAtkcIYYGB7ZSCZOkrRd+Fet7XjEkBWIkQwhpNvcSsodqLA9VJK/XHNhuVDEoBWIQQwkoTOqs7JnjtPuMYukjff46tIQALEEJYrL2vWDRa3tNDHbxIT9/NxhBAQyOEsJ4ixLSe6prx8u395i3L9K9LrV4QADshhHAX17dQNibKAW2VqIzKxYcYkwJoIIQQbsRLFWlR2vyb5MMbjOl5xkXd6gUBsAFCCLczuL1SkCxLdRGTqRd8w9YQQP0ihHBHzbzE3Djt6Wh17DI9rcAwqCGAekMI4b5u7aJuTpRrjjuGLdK/OE8MAdQLQgi3FhqgrBwnb+uqxn6sz/6MkysAuB4hhLtThEgNVz8ZJ18rNG9daZwpt3pBABoXQgjPEN5S2ZAouzUTfRboSw4zJgXgMoQQHsNHE8/1196+UZuWZ6TmG+WG1QsC0CgQQniYkSHKronyVJmIztS3n2ZrCKCuCCE8T4sm4r3h2pNR6til+swdpkkNAdQBIYSnurWLmp8gsw+bNy/Rj1wkhgBqiRDCg3VuqqwaJ0eFqFEZ+rwiTq4AUBuEEJ5NU8SjEWr2aPnXbebklcZZTq4AUEOEEI1B/7bKtmQZ5CeiMvU1XzMmBVADhBCNhK8U6bHa/w3R7lhlpOYbFQxKAVQPIUSjMrqDUpAkvywWQz7W951jawjg2gghGpt2viJrlPZgb3XoIj19N+dWALgGQojGKSVM3ZgoPzxojlmiHyuhhgCuiBCi0eraVPl0nIxtr0Rn6Flf8ZohgKoRQjRmUhVpUVrGKPn7jWbKaqO40uoFAXA/hBCN36B2SkGy9JWi70J93QnGpAB+ghDCFpp6iVlDtRcHqhNX6I9tNioZlAL4HiGEjUzorO6Y4LXrjGPoIn0/J1cAEEIQQthNe1+xeLS8p4c6eJGevpuNIQBCCPtRhJjWU/10vHx7vzlhhfFNmdULAmApQgib6tVC2Zgo+7YSfRdWLj7EmBSwL0II+/JSRVqUNm+4fGC9MT3PuKhbvSAAViCEsLvhQcrOCbJEF/0z9YJv2BoCtkMIAdG8iXgnTnsqWh27TE8rMAxqCNgJIQS+c2sXdXOi/PS4Y9gi/cAFYgjYBSEEfhAaoHwyTt7WVR2Upc/+jJMrAFsghMBPKEKkhqufjJOvFpqTVxpnyq1eEIB6RgiBKoS3VDYmyq7NRFSGvuo4Y1KgMSOEQNV8NPFcf23OMC1ltZGab5QbVi8IQP0ghMDVjApRdk2Up8pEdKa+4wxbQ6ARIoTANbRoIt4brv2hrzoyR5+5wzSpIdC4EEKgWlLC1M2JMvuwOXqpfvQiMQQaD0IIVFfnpsont8iRwWpkhv7eF5xcATQShBCoAamKRyPU7NHyqQJz8krj2wqrFwSgzgghUGP92yrbk2WQn4jM0Nd8zZgU8GyEEKgNXynSY7XXh2hTVhmPbTYqGJQCHosQArU3poOyNUkWnhVDPtb3nWNrCHgkQgjUSTtf8fHN2oO91aGL9PTdnFsBeB5CCLhASpi6drx8p8gcu1Q/VkINAU9CCAHX6NlC2ZAgB7VTojP0j7/iNUPAYxBCwGWkKtKitIxR8ncbzZTVRnGl1QsCUA2EEHCxQe2UrcnSV4qIhfq6E4xJAXdHCAHXa+YlZg3Vnh+gTlyhP7bZqGRQCrgxQgjUl4ld1O0TvHadcQxbrBedZ2sIuClCCNSjQF+xeLS8+zp16CJ99mdsDAF3RAiB+qUIMa2n+sk4+X97zYkrjNPlVi8IwE8RQqAh9GqhbEyUfVqJPgsqsw8zJgXcCCEEGoiXKtKitHfj5P15xvQ8o0S3ekEAhBCEEGhgI4KVXRNliS5iMvVtp9kaAtYjhEBDa95EvBOnPRWtjlmqpxUYBjUELEUIAWvc2kXdnChXH3fcsFg/cIEYApYhhIBlQgOUVePk5C7qoCz9nSJOrgCsQQgBKylCpIarK2+RL+40J680znByBdDgCCFgvT6tlI2JsmszEZWhrz7OmBRoUIQQcAs+mniuvzZnmDZ1tZGab5QbVi8IsA1CCLiRUSHKtmR5qFjEZOo7zrA1BBoCIQTcSxsfkTFKe6SvOjJHn7nDNKkhUM8IIeCOUsLUzYly8WFz9FL96EViCNQjQgi4qc5NlVW3yJHBamSG/v4XZtF5x8+vymY6xJ6zZBKoE2n1AgBckVTFoxHq8GBl6mqjlbfw1cTi0dLv+3+1pkPcvcbo2ULp3VKxdJmAZ2NHCLi7AW2V7clyQFtl91nH0EV6caUQ31ewazPlsQj+FQN1wj8hwAP4SpEeq82Nk4eKHb3nV54uV5wVfCKSf8JAXfGvCPAYYzoouyd5tfQR1y/y9tYEFQRcgtcIAU/Szkf0a6WqDv2tz4WPJp6N0QK8rF4T4OH4jRLwGJdeF/x0VMXrQ7SFB82eH+lz93O1bqBOCCHgGS5V0DkRvbuH+tIgrb2veH6neVOO/tm3nEQB1BKjUcAzvLTL7NZM+cuPXhec3FWtMEXe147rWyg3Zuv3X68+FqH5aBauEfBIDRTCoqKiL774onv37t26dWuYZwQamYfDVa+fTXB+2V29ravwUsWtXZXHNpl9FuivDNZGd+C0QqAGGmI0+tprrz388MOffvrp5MmT//GPfzTAMwKNz88r+OPbg/2UuXHaP2O1+9cZ8cv1w1yVDai2hgjhtGnTFi9e/Le//W3hwoWvvfZaAzwjYE/jOip7JsroNkq/hfrMHaZBDYFqaIgQSvndAHb79u09e/ZsgGcEbMtXirQobUOC/OSYGZOp558khsA1uPI1ws8//zw9Pf3Ht3Tr1m3GjBnOj4uKiv70pz9lZma68BkBVCmsubJsrFx0yLz9E+PGQOXlQVobH6vXBLgrF4Tw9OnT/v7+Pj4+gYGBKSkpP/6jpk2bOj/46quvJk2a9Pbbb3fv3r3uzwigOuJD1RsC1Se2Gr3mVz4To93XU+UoGuDnqhXCoqKiRx55pKCgoLS09OTJk5duP3LkSFJS0rFjx0pKSp544okZM2YMHDjw519+5MiR5OTkWbNmRUdHu2zhAKqheRORHqvdeZ36m3XGu0Xm60M03qoCuEy1XiPUNC0xMfGFF144e/bsj29/5JFH+vfvf+zYsS1btjz11FN79+6t8svnz59fWlr6wAMPxMTE3HDDDS5YNYCaiGytrIuX9/ZQR+boqfmG8/0rADgpDkd1X0vftWtXVFRUZeV3/4ZKSkpatmy5c+fOHj16CCGmTJnSpUuXZ555ptZLiY+P37t3b2hoqPNTTdPeeOONVq1aXen+xcXFAQEBtX461EVZWZmU8tJhUGhgtf7mP1mm/Hm7XHNS/XtkZXJHrs1WGw6Ho6SkxN/f3+qF2FRJSYm3t7emVffKET4+Ptf8SVX7H2RHjhzRdT0sLMz5aY8ePfbt21frRxNCdOzYMSQkZNKkSc5PfX19O3TooKpX3LM6HA5CaBX5PasXYlO1/uYPCBDvjRRrvnb8Zp3y4WHx71itc1MmpTXjcDgUReGHj1VUVa1RCKuj9j/Izp8/7+PjcylU/v7+3377bV2W4ufnFxQUNHLkyLo8CIBruiFQ2ZYsXys0B2Tpv+mlPh6heXNhNthY7c8jbNu2bUlJSUVFhfPTb7/9tl27di5aFYD65aWK1HC1IFnuOSv6LtRzj3K6Ieyr9iEMDg5u3br1li1bnJ9u3ry5T58+LloVgIbQwV/56CbtxYHqfWuNySuNE6VWLwiwQrVCWFFRsWLFig0bNjgcjhUrVqxdu1YI4eXldc899zz++OOfffbZ3LlzN2zYMHXq1HpeLQDXiw9VCyfJXi1F34WV6bu5MBtsp1qvEZaWls6ePVsIMWHChNmzZ7dq1WrYsGFCiKeffvovf/nL5MmTAwMDs7Oz27ZtW7+LBVA//KRIi9Lu6KY+sN6YW2S+PkQb0JaDaGAXNTh9or79/ve/DwoK+t3vflfN+1+4cOHSlWvQwDh9wlr1983vEOKd/eajm4yETurzA7TmTerjSTybw+G4ePEiR41apaanT1QH71AP4AeKEClh6t5bvXw00XuBPne/6S6/KQP1hhACuFyLJiI9Vssapb1SaI7I1vd+Sw3RmBFCAFWLbqNsSJC/6KYOW6Sn5hsXdasXBNQPQgjgilRFTOup7prodbZc9F2g5xxma4hGiBACuIYgPzE3TvvPDdrvNxrxy/VDxeQQjQohBFAtw4OUHRPk0PZqZIaeVmBUcMluNBaEEEB1eani0Qh1U6LcdMrRP1Nfd4KtIRoDQgigZro1U3JGy2di1CmrjJTVxqkyqxcE1A0hBFAb8aHqnkky2F/0nl+ZvtvkfEN4LkIIoJb8pXiuv5Z7i/zwoDkwS9/yDTGERyKEAOokopWSFy8f7K3GL9NT843zlVYvCKghQgigrpwXZtszyUsI0Xu+Pnc/R5TCkxBCAK7Rylukx2rvD9de3GXelKN/xoXZ4CEIIQBXGhqoFCTJhFD1xgY3UmkAABHuSURBVGw9rcAoM6xeEHAthBCAi0lVpIar25LlgfOizwJ92RG2hnBrhBBAvQj2U+bGaf+M1e5fZ8Qv1w9fJIdwU4QQQD0a11HZM1FGt1H6LdRn7jANagj3QwgB1C9fKdKitA0J8pNjZkymvuEkMYR7IYQAGkJYc2XZWPl0tHrbJ0bKauMbLswGt0EIATSc+FB15wTZ0lv0ml85+zOuywa3QAgBNKjmTUR6rLZsrHzzc/PGxfqes9QQFiOEACwQ2VpZFy/v7aGOzNFT841iLswG6xBCANZQFZESpm6f4HW2XFw/X//oIBdmgzUIIQArtfcVc+O0ecO1pwrM+OX6lxeYlKKhEUIA1rshUNmWLEcGqwOy9LQCo5wLs6EBEUIAbsFLFanhakGy3H1W9F2o5x5la4gGQggBuJEO/sr8m7QXB6r3rTUmrzROlFq9INgAIQTgduJD1cJJsldL0XdhZfpuLsyG+kUIAbgjPynSorS14+Xiw+aALH3TKWKI+kIIAbiv65ory8fK1N5q4nJ9ep5xrsLqBaExIoQA3JoiREqYuvdWLx9N9F6gz93PhdngYoQQgAdo0USkx2pZo7RXCs0R2freb6khXIYQAvAY0W2UDQnyF93UYYv01Hzjom71gtAoEEIAnkRVxLSe6q6JXmfLRcRCPecwW0PUFSEE4HmC/MTcOG3OMO33G4345fqhYnKI2iOEADzV8CBlxwQ5tL0amaGnFRgVXLUbtUIIAXgwL1U8GqFuSpSbTjn6Z+rrTrA1RI0RQgAer1szJWe0fCZGnbLKSFltnCqzekHwKIQQQCMRH6rumSSD/UXv+ZXpu03ON0Q1EUIAjYe/FM/113JvkR8eNAdm6Vu+IYa4NkIIoLGJaKXkxcsHe6vxy/TUfON8pdULgnsjhAAaIeeF2fZM8hJC9J6vz93PEaW4IkIIoNFq5S3SY7X3h2sv7jJvytE/48JsqAohBNDIDQ1UCpJkQqh6Y7aeVmCUGVYvCG6GEAJo/KQqUsPVbcnywHnRZ4G+7AhbQ/yAEAKwi2A/ZW6c9s9Y7f51Rvxy/fBFcgghCCEAuxnXUdkzUUa3Ufot1GfuMA1qaHuEEIDt+EqRFqVtSJCfHDNjMvUNJ4mhrRFCADYV1lxZNlY+Ha3e9omRstr4hguz2RUhBGBr8aHqzgmypbfoNb9y9mdcl82OCCEAu2veRKTHasvGyjc/N29crO85Sw3thRACgBBCRLZW1sXLe3uoI3P01HyjmAuz2QYhBIDvqIpICVO3T/A6Wy6un69/dJALs9kCIQSAn2jvK+bGafOGa08VmPHL9S8vOIQQp8uruOeFSlFBKz0fIQSAKtwQqGxLliOD1QFZ+pMFxqQV+uJDP3nt8HylGL9M33KKFxQ9HiEEgKp5qSI1XC1IlnvOiiMXxeNbjKyvvtsAnq8U8cv0/9dHHdxesXaRqDtCCABX08FfmX+T9vIg9Vy5uHet8fZ+84KuOCuY1IkfoY2BtHoBAOAB4kPVEcHqHzcb09aZ7by9Xh2iJlLBxoL/kQBQLf5S/DVG691cOV0uXtzl4I3vGw1CCADV4nxd8IkopTC+8sB5M+zDyo1cpLRRIIQAcG2Xjo5J6qS29XEUTvJq462MWaqn7+b8CY9HCAHgGhxCTMjVZ/zo6JjmTcS6BNm9uTJnn5mca5yt6ixDeApCCADXoAjx/gh52dExLZqIFWNlfoIMDRCRGbyXkwcjhABwbW19qrixeRPR1Eukx2r/jFUTluszd/DmFR6JEAJAXSV1UjcnycyvGJN6JEIIAC7QKUBZM172bCEiM/R8xqQehRACgGt4qeK5/lp6rJqUy5jUkxBCAHClxE7q5kSZ9ZWZlGucYUzqCQghALhYaICyZry8voWIzNDXn2Bn6O4IIQC4nlTFc/21f8eqySv0tAKDOak7I4QAUF8SOqlbkmTuUUdSrlHlW/vCHRBCAKhHHf2VT8fJqDYiKkNfx5jULRFCAKhfUhVpUdorg9UJjEndEiEEgIYQH6puSZIrjjpGL9VPlFq9GvwIIQSABtLRX1k9Tg5prwzI0vO+ZmPoLgghADQc55j0jWHabZ8YjEndBCEEgIY2MkTZmKitPOa4eYn+NWNSqxFCALBAB39l1S1yaKASlVG54igbQysRQgCwhnNM+k6cvHMNY1IrEUIAsNJNwcrGBO2TY45RjEktQggBwGIh/sqqcXJYoBKVUZnLmLTBEUIAsJ6miLQo7d04edcaI63AMKhhAyKEAOAuRgQrBUly/QnHqBz9eInVq7ENQggAbqSdr1gyRt4QpERlVC5nTNogCCEAuBfnmPS9EfKuT43HNjMmrXeEEADc0fAgZVuy3PaNY2SOfqyEGNYjQggAbqqdr1g6Vo7poEZn6MuO0ML6QggBwH0pQjwaob4/Qt6z1kjNNypNqxfUGBFCAHB3cUHKtmS575xj1BLGpK5HCAHAA7T1EUvGyLEd1OgMfclhWuhKhBAAPINzTPrBCHlfHmNSVyKEAOBJbgxStiXLz885RuboRy+yNXQBQggAHqatj8gZIyd0Vvtn6TmMSeuMEAKA51GESA1XPxghpzMmrTNCCACe6oZAZVuy3H/eMXSR/uUFtoa1RAgBwIO18RHZo+Ud3dTBi/RsxqS1QggBwLM5x6SZo+SD643UfKOCMWkNEUIAaAwGtFU2JcovzjuGLtIPMiatCUIIAI1EGx+xaLSc0k0dmKXPP8jGsLoIIQA0Hs4xafZo+egmkzFpNRFCAGhs+rdVCpLl0RIx5GP9AGPSayGEANAINW8iPrpJ+2V3dWCW/hFj0qsihADQODnHpEvGyMcYk14VIQSAxiymjVKQLI+XiMEf61+cZ0xaBUIIAI1c8ybifzdpU7urgz7WPzzAxvByhBAAGr9LY9I/bjGn5xnlhtULcieEEADswjkm/bZCDF6kFzEm/R4hBAAbaeYl/jdCm95Tjf1Y/x9jUiEEIQQAG5rWU106Rv5pi5my2ijVrV6N1QghANhRdBulIFlWmGKI7cekhBAAbKqZl/hghPZwuDpkkf7BF/YdkxJCALC1lDB1yWj5l61mymqjxJZjUkIIAHYX1UYpSJaVDjEgSy/81nZjUkIIABBNvcT7w7U/9FWHZ+vv2WxMSggBAN9JCVOXjpFpBfYakxJCAMAPIlsrW5Ok7hD9M/U9Z20xJiWEAICfaOol3huuPRqhjsjR5xU1/jEpIQQAVCElTF01Tv59h5my2rjYqMekhBAAULVeLZQNCdIUon+mvrvxjkkJIQDgigK8xLtx2mMRatxifc5njXNMSggBANeQEqaujZf/2tM4x6SEEABwbde3UDYlSm9NxGTqu840qjEpIQQAVIuvFHOGaY9HqMOz9fTdjWdMSggBADWQEqbmxcs3PjdTVhvFlVavxhUIIQCgZnq2UDYmSF8pYjL1nZ4/JiWEAIAa85Vi1lDtj/3Um3I8fkxKCAEAtZQSpq4dL//7uTnVk8ekhBAAUHs9WygbE2UrbxGTqe/wzDEpIQQA1ImPJtJjtT9HqiM9c0xKCAEALvDL7mpevHzzc3PSSuNchdWrqQlCCABwjR7NlQ2JMsRPDMjSt5/2mDEpIQQAuIxzTPq3GHXMUo8ZkxJCAICLTeyirhkv39pvTlxhfOv2Y1JCCABwveuaK/kJsoO/iMzQN5506zEpIQQA1AvnmPSlgWpirluPSQkhAKAeTeisbkiQ7x8wk3ONs+VWr6YqhBAAUL86N1U+HSdDA0RUpr7B/cakhBAAUO+8NZEeq708UE3K1dN3m24VQ0IIAGggyZ3VjYnyAzcbkxJCAEDD6RSgrB4nOwWIyAw93z3GpIQQANCgnGPSf8aqicv1mTusH5MSQgCABZI6qVuSZOZXZlKuccbSMSkhBABYIzRAWTteXt9CRGbo609YtjMkhAAAy0hVPNdf+3esmrxCTyswLJmTEkIAgMUSOqmbE2XuUYclY1JCCACwXmiA8uk4GdVGRGXo6xp2TEoIAQBuQaoiLUr792B1QsOOSQkhAMCNxIeqW5LkiqOO0Uv1E6UN8YyEEADgXjr6K6vHySHtlQFZet7X9b4xJIQAALfjHJO+OlibuFL/7Xr9mW2Xv4tTqS4eWG9UuOLNnQghAMBNjQ9VtibJHWfErM+MB/ONS7eXG+K2T4wBbZUmrogYIQQAuK8O/sqqW+TdPdS3Pjd/udoQQpQb4taVxsQuyq/CXJMw6ZJHAQCgnkhVPBWlDW2vJufqRd9qbXwdt3ZVXVVBwY4QAOARRoUoBUly9zn1fKXDhRUUhBAA4BHKDfH7TebLMUbvFspjm41rf0G1EUIAgLu79LrgL7uY/45VzlcIF7aQEAIA3NplR8coQrw6RHNhCwkhAMCtHS1x/KLbT44RVYR4ZbDW1kex9XmEJ0+efPPNN61ehX0tXbp0+/btVq/Cvl5++WVd161ehU0dPnx43rx5Vq/CXro2VX7R7btaZWVlFRYWCiFURfyuj2rr8wj379///vvvW70K+1q2bFl+fr7Vq7Cv9PT08+fPW70Km9qzZ8+CBQusXoV9ZWdnb9myxbWP6akhBADAJQghAMDWCCEAwNYUh6NB3wj4KsaNG7dp06ZmzZpV587l5eVnzpwJCgqq71WhSqdPn27SpEnTpk2tXohNffXVVx07dlRVfpG1QGlp6blz5wIDA61eiE2dOnXK19c3ICCgmve/4447/vrXv179Pm4UwuLi4uPHj2uaVs37l5eXe3t71+uScCW6riuKUv3/WXAtvvmtxX9/C1VWVmqaVv3fAoOCgnx9fa9+HzcKIQAADY/RCgDA1gghAMDWCCEAwNYIIQDA1jzvHeoNw9i7d+/27dsNw/jVr35l9XJs5+jRo9nZ2fv27Wvbtu2UKVM6duxo9YpsxOFwLFiwYMeOHRcvXuzVq9cdd9zh5+dn9aLsqKCgYMuWLbfffns1T/eCS+Tn5+/atevSp/fcc4+rDlz3vB1hdnb2mDFjXn311RkzZli9Fju688478/LyOnTocPDgwV69eu3Zs8fqFdmIaZrvvvuun59f586d33nnnbi4OMNw5duTojrOnz8/derU6dOnnzp1yuq12MtHH3309ttvH/ieC0958LzTJ0zTVFV1/fr18fHxp0+ftno5tlNWVubj4+P8eOLEid27d585c6a1S7KnkpKSZs2a7d69u2fPnlavxV6mTZvWr1+/Bx54oKioqFu3blYvx0ZmzJgREBDw9NNPu/yRPW9HyNU0rHWpgkKIsrKy6l/fAa61fv365s2bh4SEWL0Qe1m1atX+/fvvvvtuqxdiU9u2bZs5c+YHH3xQVlbmwoclKqilnJycjRs33nfffVYvxHYmTpwYGBg4YcKEDz/8kKvcNaSSkpLf/va3r7/+uqIoVq/FjoKCgkJCQs6fP5+enh4REXH27FlXPbLnHSwDd7Bx48Zf/epX7777LldcbHhvvfXW+fPn58+ff/vtt+/cuZMr7jaYxx9/fOrUqT179iwvL7d6LXb0yCOPOD8wTXPYsGGvvPLKX/7yF5c8MjtC1Ni2bdsSExPfeuutMWPGWL0WO2ratGlISEhqamrHjh1zc3OtXo6NzJs374MPPoiJiRk8eLAQIikpKTs72+pF2ZGqqrGxsQcOHHDVA7IjRM3s3LnzlltumTVr1rhx46xei+2UlZV5e3s753KnTp368ssvO3XqZPWibGTlypW6rgshKioqBg8e/MwzzwwaNMjqRdlIaWmp8/LZZWVlubm5KSkprnpkzztq9MiRI0lJScXFxV988UVERESXLl0++ugjqxdlIzExMUVFRd27d3d+OmrUqL///e/WLsk+cnNzf/3rX0dHRyuKsmrVqoSEhDlz5vB6VcMrLy/38fHhqNEG1qVLl969e7do0WLt2rVdunRZsmTJNd9Wopo8L4Tl5eW7d+++9KmPj0/v3r0tXI/dFBYWlpaWXvq0ZcuWXbt2tXA9drN3797CwkJFUcLDw6+77jqrl2NTDoejoKAgPDycN2NqSEeOHCkoKCgpKenevXtMTIwLH9nzQggAgAtxsAwAwNYIIQDA1gghAMDWCCEAwNYIIQDA1gghAMDWCCEAwNYIIQDA1gghAMDWCCEAwNYIIQDA1v4/dL4ZRwObiWMAAAAASUVORK5CYII=", "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.8.2" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.2", "language": "julia" } }, "nbformat": 4 }