{ "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 - 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 Eₙ-Eₙ₋₁ ρout-ρin Diag\n", "--- --------------- --------- -------- ----\n", " 1 -7.844874202097 NaN 1.97e-01 4.0 \n", " 2 -7.850310197326 -5.44e-03 2.93e-02 1.0 \n", " 3 -7.850646849666 -3.37e-04 3.04e-03 3.2 \n", " 4 -7.850647502290 -6.53e-07 4.48e-04 2.2 \n", " 5 -7.850647510934 -8.64e-09 1.41e-05 1.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+gvaeTAAAgAElEQVR4nO3deXwU9f3H8e8cmzvk4AinEAwQbhEiZyCEECDcSEBFRSzg8bMe9Cfgr1aRVoX+tIpWW/y1Wm0plUsFDEcSzghKSMIpl0AQOQqBEJIQQmZmf3+sxci5STY7uzuv5x8+JuNk5pO47jvfne/nO5LdbhcAAFiVbHYBAACYiSAEAFgaQQgAsDSCEABgaQQhAMDSCEIAgKURhAAASyMIAQCWRhACACzNfUF48eLFS5cuue1yAAA4w01B+MILL4wYMaJ///7vvvvuzY65cuXK7t27nT+nruuuKA3VZBiG2SVYF798c/HmYyLDMFy+Mqg7gvDAgQNr1qxZt27dxo0b582bd/78+RsedvLkyREjRjh/WsaX5uL3b6JLly6xSrCJePGbqLy83OV/CLojCLds2ZKYmCjLckBAQFxc3Pbt291wUQAAnOGOIDx//nydOnUc2+Hh4efOnXPDRQEAcIY7gjAyMrKoqMixfeHChXr16rnhogAAOKNGQZibm/vss88mJCRMnjy58v7Dhw8nJyfXr18/ISFh7969vXv3Xr9+vWEYly9fzs7O7tq1a81qBgDAZWoUhCdOnAgPD2/Tps2BAwcq77/vvvu6det2+PDhlJSUMWPGxMTEDBkypH///n379n3uueciIyOrfcXiCnH5RtO1zl6u9ikBAJYm1Xzu2f/93/998sknmzdvdnyZl5fXt2/fs2fPBgQEGIbRpEmTBQsWJCYmlpSUKIoSGBh4s/Pk5+d369Zt1KhRlXfWq1dv1qxZV79ccVx8+J28sK8RoIji4uLQ0FAhxIeHpA3/lj7pw4RytyopKQkJCTG7CosqKSkJDg6WJMnsQizq6psP3K+srMzPz09RFCeP9/Pzk+XbDPnUGld1rX379rVt2zYgIEAIIctyp06d9u3bl5iY6Mybps1mu+aD0zp16lT+gUe1EEWauG+TsihBKIqiKMpHh8TKE+LTBOH87wUu4fj9m12FRTl++QShWXjxm0j5DyePd+Z/E9cH4fnz5yv/rRQWFlZQUODk9wYEBDzxxBO3PmZyW6EoxviNxsfdbZ8cUb44biwbqAbwmnQ7m81ms9nMrsKiHL98gtAsvPhNpGmazWZz7R8irg/CiIiIkpKSq18WFRXVrVvXtZeY1Fq2C5G0zj8iwFg7WPUnBQEA1eX69omYmJiDBw9WVFQIIex2+759+1q1auXyq+iGuKyJr8/Yn9yi/1DKEhsAgGqq0Yjw0qVLp0+fLigouHz58pEjR4KDg6Oiou65557GjRu/884706ZN++ijj4QQAwYMcFG1P/rrAePzY8ZXg8qXnAp+e4/eeakxtqX8265Kg5tOxAEAF8vIyPjggw/MrsJCBg0a9Itf/KI2zlyjIPz666+nTJni2B44cODgwYPfe+89SZIWLlw4adKkF198MSYmZvHixarqyg9g/3rAWJZvLE1SKy6JJ9rKAYr4x3dGqE20W1IxOVb+9V1KKB/dA6h927Zt0zTt/vvvN7sQS9iyZUt6eronBmFiYuLhw4ev39+pU6ecnJyanPlmPj9mfH7MWJak+iuiQgjxn/uFa3+w54xWX9thtFpU8VwH5dkOMjcOAdS2tm3bpqamml2FJRiG8dlnn9XSyb3swbxJjeUlA66dHfNoa/m93krzEGl+HyUjRc0psLdZrH2w39C5dQgAuB0vC8IQm7jhUK+u/48bHSKkRQOUf/ZX/vGd0XmZtvgoXfYAgFvxsiB0Uq8oadMw9e0eyu/yjN4rtE2nGRsCAG7MN4PQIamJlDdafbaDPGmjPnCVtvM8cQgAuJYvB6EQQpZEarS8L1VNjZaHrNbGZepHiolDAMBPfDwIHfxkMTVWPphq61pP6v6F9liW/u8ys2sCAHgGSwShQ4hNzOgs7xtri/AX7ZdUzMzWL1aYXRMAi0k7br/trIVvL9g/OcREP/exUBA61AsQc+KU3NFqYblovahi7k7jhg84BIDa0K2e9Pw3+vpTN83C/RfsE9br3eqzorr7WC4IHe4Ikeb3UdYNpekQgFs1CBQrktWZ226chfsv2O9fry/or7QLr04QJiUlXbhwocY1CiGE3W4vLCys+QNrr1y5Eh8ff+XKlVscM2bMmO+//76GF6oJiwahQ7twadEA5dNEZeFho+NSbfFRgzQEUNtuloU1TEEhxJ49ezRNc0WNory8PDIysqyspvMpZFnu0KHDrR+Nu2/fvvLy8hpeqCYsHYQOPRpI64eq7/RUXt9h9Fqubbj5RxYA4BLXZ2HNU/B6O3bs+Pjjj7OysirvPHLkyKeffrpu3bqrkbl3797jx487ts+dO7dt2zYhRHZ2thBi/fr1GRkZzo8y7XZ7enr6xx9/fODAAcceRVEeeOABx+MDs7OzCwoKtm3btmDBgvz8fBf8hC5CEP4oqYm0fZQ6raM8ebM+cJWWd444BFCLKmdhbaTgW2+9NX78+J07dz7xxBNXH3j+97//vXfv3llZWbNnz+7Xr9/ly5eFEHPnzr26jGdOTs5TTz0lhMjMzBRCLFu2bPHixU4+XF3X9bFjx7766qt5eXkJCQn//Oc/hRDl5eV9+/Z1DPiefvrpRx99dO7cuRs2bOjSpcvu3btd9cPWkOsfzOu9HE2Ho5rLHx00hq3Re0dJr8XJMXW4ZQ3AWU9t0SuqMt+zTZg0NlMzDDGoqTRvT9Vmir7R/aYP2zl79uwrr7yyb9++Ro0aFRUVRUdHP/HEE61bt37mmWe++OKL+Ph4wzDi4+M/+OCDp59++oZnmDlz5iuvvPLuu+8GBQU5Wc+KFSvy8/Ozs7NlWR4zZsy4cePGjRt3zTEtW7Z8++23hRBBQUF///vff//73zv709YmgvBaNllMjZUnxMh/3Gv0XK6NaSHPultp5OwrAYCldakrVWnm3akye+ApSSj2ZiFSqyr+2a3e/PCcnJzY2NhGjRoJIcLCwnr27OnIJ7vdHh8fL4SQZXnIkCHbt2+v0hWFEBcuXHj//feFEIqizJgxo/K/2rJlS3JysuN2YN++fYuLi7///vvGjRtXPiYhIcGxERMT880331T16rWEILyxYFXM6CxPjpX/d5fecWnF5Fh5Zmcl3M/ssgB4tl+0qcL9pv0X7PevN9YOUer5S8PXainN5P6NXPMRVFFRUeUvJUlyzP80DOP6nVc3hBC3nt559Ruv/rOyixcvhoeHV95z/aRTP78f30YdqezEj+IO3CO8lbr+Yk6csmOMWlgu2iyumLvTKHPNhCwAVlf5vuCteyqqZ9euXf/+97+FEMXFxVu3bo2Li7vzzjsVRdm6dasQwm63r169ulu3bkKIRo0afffdd47v2rRpk2PDz8/PZrOVlpZec9rw8PAXXnjhhRdemD59+vUXzczMdGRtVlZWSEhI8+bNXfXj1CpGhLfXNFia30d5roP8Uo7RZrH2Yhf50dayyp8QAKrr+tkxjiwcvlabc4/iknFhs2bNRo4c2a9fvzVr1owePbpz585CiLfeeuvee+8dP378rl27Kioqpk6dKoR48MEHHa1+Fy5cOHPmjOPbZVkeNmzY0KFDW7RoMXv27NjYWCevO3jw4M6dOy9YsOCNN95QVdVV7Ry1iiB0Vmy4tGiAsu2s/YVs/a09xuyu8thomYk0AKrqZnNEXZuFzZs3/9vf/rZ27drBgwdfvTM3ceLE3r17f/311ykpKQkJCTabTQjRoUOH3NzcjRs3xsTEtG3b9vDhw46DlyxZsm/fvlOnTjVs2NDJi44ePXro0KF5eXlr167t0KGDEMLf33/Dhg3+/v5CiHfeeSc6OtpxpCOka/gzugpBWDX31JcyU9SME/aZ2fobu4zX45TExqQhAGfdulPCtVnYuHHjRx555JqdMTExMTEx1+yMjo6+GlF169Z1bMiy3L59+/bt21fpop07d3aMPh0kSboaeHFxcVf3N23atGnTplU6c+0hCKsjqYmU3URdctR4/Cu9eYiYE6d0rUccAri9oivi00SlddhN3zEaBIrlyWoNn58aFRXVpUuXmpyhGmJjY6Oiotx8UZcgCKtJEj81HY5M13s1kF7tJre6+YsbAIQQ3Rvc/l0iKlAkN6nOm8nRo0cDAwOFEAkJCVc/DnWbZ599tnrfmJeX5/js1CxM+agR249POlS71pN6r9Aey9JPXvKUCcEArMaRgl4nICDg+mYMdyIIXSBIFTM6ywdSbRH+otNSbWa2Xmjm+rEAgCogCF0mwl/MiVN2/qfpcFauXsyDfwHA4xGELtYkWJrfR8karh65KFotqpi707jCg6YBwIMRhLWidZj0SYKSPuSnB//yqEMA8EwEYS3qGCktGqD8I0H5+3dG52Xa4qOMDQHA49A+Uet6R0mbh6kZJ+zTvtHn7TFej1PiG9JlAfiCrVu3zp071+wqLGHXrl21d3KC0E2Smkg7RqtL842JG/U764g3uyudIolDwItNmzattLS0sLDQ7EIsoVmzZo4HSNUGgtB9HA/+Hdlc/ttBY9AqLb6hPPceOTqUOAS8UkBAwKuvvmp2FXAB7hG6m58spsbKh8bZutaT7vlCeyxLP1Nmdk0AYGEEoTlCbGJGZ3nfWFuEv2i3pGJmtn6RpkMAMANBaKZ6AWJOnJIzWi0sF60XVczdaZTrZtcEABZDEJqveYg0v4+SmfJT06FO0yEAuAtB6CnaR0iLBigLE5UFNB0CgBsRhJ6lZwNp4zD17R7Kq3lGr+XaxlOMDQGgdhGEniipiZQ7Wn2uo/zoJn3gKm3HOeIQAGoLQeihHE2H+1LV1Gg5ZY02LlM/fJE4BADXIwg9WuWmwx7Ltcey9NM0HQKASxGEXiBYFTM6y/tTbRH+osOSipnZetEVs2sCAF9BEHqNuv5iTpySN0YtLBetF1fM3WlcpukQAGqMIPQyzYKl+X2UDUPVnAJ760U0HQJATRGEXqltuLRogLJ4gPKvw0bHpdriozz3FwCqiSD0Yt0bSOuGqu/0VObsNHou19bTdAgAVcdjmLxeUhNpexN1yVFj6ma9RaiYG6fcXY9HOwGAsxgR+gJJiNRo+duxamq0PHytPi5T/46mQwBwDkHoO2w/Nh2qXetJvZZrj2Xppy6ZXRMAeDyC0NcEVWo67Li0Yma2foGmQwC4OYLQN0X6izlxyo4xamG5aLO4Yu5Oo0wzuyYA8EgEoS9rGizN76NsGqbmFNhbL9Y+2G9oPNwJAH6OIPR9bcKkRQOUZUnKoiNGB5oOAeDnCEKriKsvZaSof+ylzN1pdP9CyzxJGgKAEPQRWk1SEym7ibrkqPHkV/odIeL1OKUbTYcArI0RoeU4mg733qumRsuj0vXha7Xd5xkdArAugtCiVFlMjZWPjFeH3yEPXKWNy9Tzi4lDAFZEEFpa5Qf/dvtceyxLP8ODfwFYDEEIEWoTMzrL+1JtEf6i3ZKKmdl6cYXZNQGAuxCE+FH9ADEnTskZrRaWi1aLKubuNMp58C8ACyAI8TPNQ6T5fZSMFDWnwN5msfbBfoOuQwC+jSDEDXSIkBYNUP7ZX/nHd0bnZdrioz8uSLO38AapeKLUXlju3voAwHUIQtxUryhp0zD1rR7K7/KM3iu0jaftL+ca7+792SptJ0rtI9P1o8w4BeC1aKjHbSQ1kfJGq0vzjUkb9ZZ1xKofhF0YT7eXhRA/lNpHpet/7MWjgAF4MYIQtydLIjVaHtlc/ttBY1auvrdQOlduf6CpmLBF/2MvpUcDUhCAF+OjUTjL0XS4P9X2YIz0ep7Rc7XfG91JQQBejyBE1dSxiSfbym0jJH9FujdDTzvO3UEA3o0gRNX8UGofma7P76PsH1HepZ706CZ9+Frth1LiEIC3cl8QFhQU/PDDD267HGrD1dkxPRpIiiTWDFb6NZI0u+j6uTZvDx2HALySm4Jw5MiRAwYMSElJcc/lUBvsQkzapL9XaXaMIol/JChBivS/9yhfHDP6rtS+vUAYAvAybgrCzz//fMWKFe65FmqJJETaILX7z2fH2GTxaaLyUCs5c6g6NVZO/FKbma2zNhsAL+KmIJQk5hb6AtuNXi+qLCQhJCEebiXvGGM7WSo6LNUyTzI0BOAdXNxH+PLLLxcXF1fe89JLL4WHh7v2KvBYDQPFJwnKyu/tv9ik920o/aGHUi/A7JoA4JZcEIRXrlzRNC0oKEgIMWTIkCtXrlT+t4GBgTW/BLzLsDukxMbq7Dy98zLt9Tj54VZMTgbguZwNwunTp2/YsOHw4cMffvjhyJEjHTvtdvv06dM/+OADWZYTEhIWLFjQo0ePWisV3iRIFXPilPEt7VM36/86YrzfS2kRysfjADyRs3+qN2/e/He/+13dunUrD/jS0tIWLVr03XffnTlz5uLFi2+++ebNvn3mzJmjRo26fPlyt27dmDVjHV3qSltHqIOayPd8oc3K1SuM238LALiZZLdXYVJDx44dX3rppdTUVMeX48ePj42NfeWVV4QQy5cv/9WvfnXo0KFql5Kfn3/33XffddddlXc2atRo/vz5Nzy+pKQkJCSk2pdDDZWWlgYHBzt58LFS6dnt6r/LpHfv0bpGkoc1VVpaGhQUxBw0s/DmY6KysjI/Pz9FUZw8PigoSJZvM+Sr0T3Cw4cPjxo1yrEdGxubn5+v67rz9V0vKCjof/7nfyrvCQwMvNkLzm6381o0l/O///YhIn2oWHzUuG+zNKK5/EZ3JdRWq6X5vuDgYILQLLz5mEhRlCoFoTNqFIRFRUVXxwQhISGappWUlISFhVX7hDabLSkpqSYlwZOlRssDm8gztultl2jv9JTHtGASDQDz1eidqH79+hcuXHBsX7hwwd/fv06dOq6oCj4r3E/M76Ms7K+8uN1gkVIAnqBGQdi+ffvc3FzHdk5OTvv27fmsBs6IbyjtHKP2iZJZpBSA6Zz9aDQ3N/f8+fMlJSW7d++OiIjo1q1beHj41KlTk5OTR44cGRUV9dprr02bNq1Wa4UvscliRmd5TAvpsSx98VFjfh+lfQR/RQEwgbNB+Omnn+bm5sbExGzdunXr1q1/+MMfwsPD4+Li3n///RkzZpSVlT300EOTJ0+u1Vrhe1qFSZlD1b8fMgakaY+0lmfdrQS48hY4ANxe1donalV+fn7//v2PHj3q5PHFxcWhoaG1WhJuwbUzyE+Xienf6FvP2P/UW0lqwtDwNkpKSpg1aiLefExU1fYJZzBtDx7BsUjp2z2VKVn6wxv0gstmFwTAMghCeJChzaS996qNg0W7JRUf7KfvHoA7EITwLI5FStcOUf9vv9H/S+1gkad8dA/AVxGE8ER31ZW2jlBHNZf7rNBm5epXGBwCqDUEITyUKotnOsjfjFS/OWOP+1z75gxDQwC1giCER4sOlVYNVl/sIo9K1x7L0osrzC4IgM8hCOEFUqPlfam2AEW0W6ItPcrnpABciSCEdwj3E/N6Kv/sr/wmh0VKAbgSQQhvwiKlAFyOIISXcSxSmjVMXf69Eb9S21tIGAKoEYIQXqlVmJSRoj4WKyelaTOz9cu62QUB8FoEIbyVJMTDreS8MbaTpaLDUi3jBENDANVBEMK7ORYpfaenMnmzPi6TRUoBVBlBCF+Q0kz6dqzass6Pi5QyNgTgPIIQPuLqIqV/OWAkfqkdYJFSAM4hCOFT7qorbRmujmoux7NIKQDnEITwNY5FSreNVLedtcd9rn3NIqUAbokghG9qESqlDVJ/100el6k/lqVfZJFSADdBEMKXDb9D3nWvGqCI9ku0JSxSCuBGCEL4OMcipQv7Ky/nGMPXasdZpBTAzxGEsIQ+DaUdY9Q+UXK3z7V5ewydNATwHwQhrMKxSOlXw9Xl3xt9V2p7WKQUgBCCIITVxNT5cZHSgWnaM1v1Us3sggCYjSCE5TgWKd0xxlZYLjov09JZpBSwNoIQFhX1n0VKp2bp4zL1syxSClgVQQhLS2km7b1XbVlHtGeRUsCqCEJYnWOR0vQU9a8HjP4sUgpYD0EICCFE50hp6wj1gTtZpBSwHIIQ+JEsiamxct5oddd50Y1FSgHLIAiBn2kSLC1LUl5lkVLAMghC4AauLlLajkVKAV9HEAI35lik9F8sUgr4OoIQuJWri5TetUybu5NFSgEfRBACt+FYpPSbkWr6CeOeL7TcAsIQ8CkEIeCUmDpSeor6THt56BoWKQV8CkEIOOuaRUrXskgp4BMIQqBqHIuUvttTeYxFSgGfQBAC1TGkmfTtvWq7CNFxKYuUAt6NIASqKVAVs+5W1gz5cZHS/RdIQ8ArEYRAjVxdpLTvShYpBbwSQQjUlGOR0h1j1N3nRbfPta0sUgp4FYIQcI3GQdLSJOXVbvJ4FikFvApBCLjS8Dvk3f9ZpPSTQ3xOCngBghBwsTA/Ma+n8mmi8vtdxvC12vclfFIKeDSCEKgVvaOkvNFqnyi5y2csUgp4NIIQqC1XFynNOGnc84WWwyKlgEciCIHaFVNHWjtEfaa9PIxFSgGPRBACta7yIqWdlmprfmBoCHgQghBwE8cipX/pqzy9lUVKAQ9CEAJu1b+RtGO02i5CdGCRUsAzEISAuzkWKU0fon540EhYySKlgMkIQsAcnSKlLcPVCTFyvy+1Wbl6uW52QYBVEYSAaRyLlOaNVvcUskgpYBqCEDBZ4yBpyQDltTj5vnX6wxv08+VmFwRYDEEIeIThd8i7xqgR/qLTMhYpBdyKIAQ8hWOR0kUsUgq4F0EIeJZeLFIKuBdBCHgcxyKl20aqGSeNuM+17SxSCtQmghDwUHfWkdKHqC/cJQ9foz2zVU//wf7m7mvvHZ4rF1OzGDQCNUIQAh4tNVreO9Z2WRePf6WnHbe/tuOnLCwsF6PTtXHRsmRifYD3IwgBTxfpL+b3Uf7aV/mh1P7BfuPX2YYQorBcjEzXXuqiJDUhB4EaUc0uAIBTEhpJeaPVl3P1P+zRdxTYinWdFARcgiAEvEaQKv73HmVIUzlljRZmE3wkCrgEH40C3qSwXMzK1f/Vp6JNuBi/Th+4igffAzVFEAJe46f7go2MdSnKoCZSiCqNTNfHZeqHiohDoJoIQsA7XDM7RpHExwlKiE1MaSN3rSf1XqE9lqWfvEQcAlXmpiA8cuTIJ598snDhwqKiIvdcEfAxJy7Zr5kdo0jio75K3QAxvbN8INUW4S86LdVmZusXrphYJuB93BGEZWVlEydOPHPmzP79+7t27VpYWOiGiwI+pkOEdP0cUVUWT7WTJSEi/MWcOGXnGLWwXLRZXDF3p1GmmVIm4H3cMWs0MDBw8+bNju09e/Z8/fXXQ4YMccN1AatpEizN76NM6yj/ZrvRerH2my7yo61llRsgwC259X+Rixcv7t69u2PHju68KGA1bcKkRQOUZUnKoiNGh6Xa4qMGdw6BW3DxiPCpp576/vvvK+95++23W7ZsKYS4cuXKAw888Pzzzzdt2tS1FwVwvbj6UkaKmnHCPjNb/99dxutxyoDGNB4CN1C1ILxw4UJgYKC/v3/lnfn5+RcvXmzfvr2iKC+//HJFRUXlf9ugQQMhhKZpDzzwQGJi4pQpU2peNAAnJTWRspuoS44aT3ylNw8Rc+KUrvWIQ+BnnP1odMqUKQ0aNIiIiPjoo4+u7rTb7ZMmTerdu/ekSZPatm177Nix+vXrN/45VVUNw3jooYfuvPPOSZMmFRYWXrnCnDbAfSTHyt33qqnRMk2HwPWcDcLU1NSsrKy+fftW3rl27dp169Z9++23OTk5SUlJL7/88g2/t7i4uKCgIDc3d9y4cePGjduwYUMNiwZQVTZZTI2VD6aqNB0C13D2o9Hk5GQhhCz/LDg//fTT1NTUsLAwIcTkyZPj4+M//PDDa44RQoSFhaWnpztzlfPnz0dHR1fe06xZs7S0tBseXFpaKkl8yGOa0tJSs0uwrtLSUrvdXr3X/3/dKR5oJr21T+64RJl4p/6rdkaYjUSsGt58TFRWVubn56coipPHBwUFXZ9K16jRZJnvv/++a9euju2WLVteunSpoKDAcVOwesLCwjIzMyvvsdlsISEhNzzYbrff7F/BPfj9myg4OLja78UhQvyht/jVXfbZecrdaca0DsrT7eVAVuB3Gm8+JlIUpUpB6IwavfYvXboUEBDg2HZslJSU1CQIFUVxTDEFUNtoOgQcavSqj4qKOn/+vGP73LlzQoiGDRu6oCgA7kLTIVCjIOzateuWLVsc21u2bGnXrl1QUJArqgLgVo6mwz/2UubuNLp/oWWeJA1hIc5+NJqZmXn48OGTJ09mZWXJspycnNyiRYtHH330zTfffPvtt2NjY2fOnDlz5sxarRVAraLpENbk7IjwyJEjOTk5CQkJwcHBOTk5Fy5cEEI0btx43bp1ubm577///q9//Wua5QFvR9MhLEiy2z3lVZ6fn9+/f/+jR486eXxxcXFoaGitloRbKCkpYeKcWUpKSmoya9RJlzTx7l7jzd366Bbyy3fLjYMYHf6INx8TVbV9whlMEQNwY0GqmNFZ3s+TDuHrCEIAtxLJkw7h6whCALfnaDrcNEzNKbC3Xqx9sN/QDLNrAlyEIATgLJoO4ZMIQgBVQ9MhfAzLCwKoDpoO4TMYEQKoJpoO4RsIQgA1wpMO4e0IQgAuQNMhvBdBCMBlaDqENyIIAbgYTYfwLgQhgFpB0yG8BUEIoBbRdAjPRx8hgFpH0yE8GSNCAO5A0yE8FkEIwH1oOoQHIggBuBtNh/AoBCEAc9B0CA9BEAIwE02HMB1BCMB8NB3CRAQhAE9B0yFMQR8hAM9C0yHcjBEhAI9D0yHciSAE4KFoOoR7EIQAPBpNh6htBCEAL0DTIWoPQQjAa9B0iNpAEALwMjQdwrUIQgBeiaZDuAp9hAC8GE2HqDlGhAC8G02HqCGCEIAvoOkQ1UYQAvAdNB2iGghCAL6GpkNUCUEIwDfRdAgnEYQAfBlNh7gtghCA76PpELdAHyEAq39eKNUAABTCSURBVKDpEDfEiBCAhdB0iOsRhAAsh6ZDVEYQArAomg7hQBACsDSaDkEQAgBNh5ZGEALAj2g6tCaCEAB+5pqmw3Un7UKI6dv0H0qvjcWlR42lRxk5ej2CEABuIKmJlD1Kfb6T/PhX+sBVWudIacRaPb/4pyxcetR491sjuSnvol6P/4QAcGOVmw5nZBt1/aWha37MQkcKrkhWQ21mV4kaIwgB4FauNh0mNZFOXbLHfa794VuFFPQlBCEA3J6j6fC78ba2EdJv99p+2U4mBX0GQQgAzlp/0pAlMaujft96nWkyPoNFtwHAKVfvC4rLZVGhfvet1/+miQmtGE54PYIQAG6v8uyY4stiSqwsS+KRTXqZbp8cq5hdHWqEIASA2ztVJq6ZHfOLNrIqi2lf65EB0pgWjAu9GEEIALf3VLsbRN3EVnKnSClltXZJEw/GkIXeiiAEgOrrUldaN1QdtEqvMMSk1mShVyIIAaBG2oZL6SnKwDS9VLvxwBEejiAEgJpqEyZtGqYkrdIrDPFcB7LQy/AfDABcoEWotH6o8qd9xuw8+gu9DEEIAK7RLFjaPExdctSYma2bXQuqgCAEAJeJChSZKerq4/bp28hCr0EQAoAr1Q8QG4apm0/bn/hK58G+XoEgBAAXC/cT6UPUAxfsj2WRhV6AIAQA1wuxiZWD1PwS+4QNusbsGc9GEAJArQhSxcpk9ZIm7s3Uy7lj6MEIQgCoLf6KWDxAUSUxJkO7TBZ6KoIQAGqRnywWDVDq+ktDVmslFWZXgxtxx8oyhmF8/PHHhw4dCg4OnjBhQosWLdxwUQDwEIokPuqnTNmsp6zRvhyk8mh7T+OOEaGu60VFRQkJCY0bN+7Xr9/FixfdcFEA8ByKJP7aV7mrrpT4pXa+3Oxq8HPuCEKbzfbss88mJydPmjSpfv36p0+fdsNFAcCjSELM66nEN5SS0rSCy2ZXg0rcdI+woqJi+vTpI0eOTEpKat26tXsuCgAeRRLiDz2U4XdIfVdqpy6ZXQ3+w5X3CK9cudKrV69rdm7fvl0IoShKamrqgQMH5syZ89xzz0VFRbnwugDgRV7pqgSpRv8vtYwUpWmwZHY5qEoQnjp16uTJk7GxscHBwVd3lpaWrlu3TpblxMTEwMBAR+xdT5bluLi4uLi4devWbd68eezYsTUtHAC81ozOsiKL+JV6ZorSMpQsNJlTQVhRUREdHX3hwoWysrItW7Z0797dsf/06dO9evVq06aNruu/+tWvvvrqq7p1617/7ceOHTt06FCbNm2OHDmycePGF1980ZU/AQB4of/uKIeoIvFLPX2I0iqMLDSTU/cIFUVZvXp1UVFRSEhI5f3z5s3r1q3bqlWr1q5d27p16/fff/+G3+7n57dixYpf/vKXCxYsWLhwYcuWLV1QOAB4ucfbyi92kRO+1PcUsiCpmZwaEcqy3KFDh+v3L1++fPbs2Y7t8ePHv/XWW7/5zW+uP6xRo0bz5s1z5kIFBQWS9LO/jKKjo3ft2nXDg0tLS685GO5UWlpqdgnWVVpaarfbef2bxYVvPvc1EaomJ6fZlva70jGcOLy9srIyPz8/RVGcPD4oKEiWbzPkq9FkmR9++KFZs2aO7aZNm544caImZxNC1KtXr7i42MmD7Xb7NSNUuBm/fxMFBwcThGZx7ZvPI+1F3VBj9Ea/Lwaq3Rvw3/Q2FEWpUhA6o0btE5qmXa3GZrNduXLFFSUBgLUMv0P+sK86Ml3b8m8GhSaoURA2bNjw7Nmzju0zZ840btzYFSUBgOWkNJP+laiOydAyTpCF7lajIIyPj8/IyHBsZ2RkxMfHu6IkALCihEbS4gHqA+u15cd4gKFbOXuP8I033igoKCgvL3/vvfc+++yzGTNmREREPPfcc/Hx8ZGRkbqu/+Mf//j6669rtVYA8G3xDaW0werwNZpuF6Nb8HQgN3E2COvUqaPr+iuvvOL40jEJp3PnzllZWQsWLJBlecuWLbGxsbVVJgBYQ7d6UtpgNWW1VqqJB2PIQndwNginTp16w/2dOnXq1KmT6+oBAKvrUldaN1RNXqVrhnikNVlY69zxPEIAQJW0DZcyUpSBaXqJJp5qRxbWLoIQADxRmzBp0zAlaZWuGeLZDmRhLeKXCwAeqkWotC5FeX+f8ds85pHWIoIQADzXHSHSpmHqoiPGzGzd7Fp8FkEIAB6tYaDITFFXH7dP30YW1gqCEAA8XYNAsWGYuum0/YmvdIOVZ1yNIAQALxDuJ9YMVnedtz+WRRa6GEEIAN4hzE+sHaLml9gf3KBrzJ5xHYIQALxGsCpWJqulmnhgvV5BFroIQQgA3sRfEYsHKLpdjE7XLjN7xhUIQgDwMn6yWDRAifCXUlZrJRVmV+P9CEIA8D6KJP7WT2kRKqWs0YrJwpohCAHAKymS+GtfpXOkNCBNO19udjXejCAEAG8lCfFOL6V3lDRwlVZw2exqvBZBCABeTBLirR7KsGZSv5XaqUtmV+OdCEIA8HqvdFUebiX3/1I7UUqzfZURhADgC2Z0ln/RRo5fqR8pJgurhucRAoCPeL6THGoTiV/qGSlKTB3J7HK8BkEIAL7j8bayKot+K/U1Q5QOEWShUwhCAPApk9vIIaoYtEpPG6x0jiQLb497hADga+67U/5Tb3nwKu2bM9wvvD1GhADgg0Y0l1VZGpmuLUtSe0UxLrwVRoQA4JtSmkn/SlTHZGiZJxkX3gpBCAA+K6GRtHiAOmG9tuYHsvCmCEIA8GXxDaWVg9RHNmqf5fMAwxvjHiEA+Lhu9aQvB6kpa7RSTTwYw/jnWgQhAPi+u+tJ64eqyat0zRCPtCYLf4YgBABLaBsuZaQoA9P0Ek081Y4s/AlBCABW0SZM2jRMSVqla4Z4tgNZ+CN+EQBgIS1CpXUpyvv7jN/mMXfmRwQhAFjLHSHSpmHqoiPGzGzd7Fo8AkEIAJbTMFBkpqirj9unbyMLCUIAsKQGgWLDMHXTafuTX+mGtbvtCUIAsKhwP7FmsLrzvP1xa2chQQgA1hXmJ9YOUY8W2x/coGtWnT1DEAKApQWrYmWyWlIhHlivV1gyCwlCALA6f0UsSVI0uxidrl223uwZghAAIPxksXiAEuEvjUrXyjSzq3EvghAAIIQQiiT+1k9pFCQNXq0VV5hdjRsRhACAHymS+LCv0ilSGpCmnS83uxp3IQgBAD+RhHinl9I7Shq4Siu4bHY1bkEQAgB+RhLirR7K0GZSv5XaqUtmV1P7CEIAwA3M7qo83Eru/6V2otTHm+0JQgDAjc3oLP+ijRy/Uj9S7MtZyPMIAQA39XwnOdQmBqTp6UOUmDqS2eXUCoIQAHArj7eVVVn0W6mvGaJ0iPDBLCQIAQC3MbmNHKyKQav0VYOVTpG+loUEIQDg9u6/U1ZlMXCVtnyg2r2BT2Uhk2UAAE5JjZY/7KuOTNe2/Nun5s4QhAAAZw1tJi1MVMdkaJknfScLCUIAQBX0byQtHqBOWK+tPeEjWUgQAgCqJr6htCJZnbhB+yzfFx5gyGQZAECVxdWXvhykpqzRSjXxYIx3j6kIQgBAddxdT1o3VB20StftYmIrL85CghAAUE3twqWMFGVgml5SIf6rnbdmIUEIAKi+NmHSpmFK0iq9whDPdvDKLPTKogEAnqNFqLQuRXnvW+O3eV45d4YgBADU1B0h0ubh6qIjxsxs3exaqowgBAC4QMNAkZmirj5un7HNy7KQIAQAuEaDQLF+qLrxtP3Jr3QvarYnCAEALhPhL9YMVneetz+WpRteEoYEIQDAlcL8xNoh6pGL9gc36Jo3zJ4hCAEALhasiuXJ6vly+4QNeoXHZyFBCABwvSBVLE9WKwwxOl277NmzZwhCAECt8JPFokQlwl8ana6VaWZXc3MEIQCgtqiy+Fs/pWGQNHi1VlxhdjU34b4gvHz58tNPP71161a3XREAYDpFEh/2VTpFSgPStPPlZldzI+4LwtmzZ2/YsOHgwYNuuyIAwBNIQrzTS+kdJQ1cpRVcNrua67gpCHfs2HH27NmEhAT3XA4A4FEkId7qoQxtJvVbqZ26ZHY1P+eOINQ0bfr06a+99pobrgUA8FizuyoPt5IT07QTpR7UbO/KxzBlZ2c//vjjlffcddddf/3rX1977bVJkybVr1/fhdcCAHijGZ1lWRLxK/XMFCU6VDK7HCGcD8KzZ8/m5uYWFxePHTu28v6cnJyMjIx69erdf//9cXFxOTk513/vqVOnli9f/uabbx4/fnzFihWRkZHDhw93Qe0AAC/0fCc5xCYS0/T0IUpMHfOz0KmPRjdt2nTHHXdMmzbt/vvvr7x/2bJlgwYNKi0tXbZsWb9+/TTtxn0if/rTn7Zv3759+/bx48fPmjWLFAQAi3uirfw/neV+K/W9heZ/RurUiLB79+7FxcX79u27++67K++fPXv2W2+99dBDD2ma1rFjxxUrVowePfoW5xkxYkTDhg1vcYCu60eOHKm8x2azNWvWzJkiAQBeZEqsHGITyav0VYOVTpFmjgudCkJ/f//rd545c2bnzp3Dhg0TQqiqOnjw4PT09FsHYVJS0q0vVFRUNGDAgMp7mjVrlpaWdsODS0tLJcn8MbVllZaWml2CdZWWltrtdl7/ZuHNx1WGRwmti5yUpi7uq3WNdGpN0rKyMj8/P0VRnLxEUFCQLN/ms8/qT5Y5efKkv79/RESE48uGDRt+/fXX1T6bQ2Rk5NGjR5082G63h4SE1PCKqAl+/yYKDg7mvdgsvPm40ENtRWSIffwm6bOBas8Gt39JK4pSpSB0RvXbJ2RZNoyfAtwwDP63BABU1dBm0sf91BFrtcyT5twvrH4QNmrUqKKi4ty5c44vT58+3ahRIxdVBQCwkEFNpaVJ6oT12toTJmRh9YOwfv36Xbt2Xb58uRCioqJi1apVgwcPdl1hAAAL6dtQWpGsTtygfX7M3Q8wdOoeYVFR0ZQpU4qKinRdHzduXGRk5J///GchxKxZsyZOnLhnz568vLyoqKiUlJRarhYA4LPi6ktfDlKHrdU0Q4yNdt9S2E4FYWBg4NSpU4UQzz//vKg0iXTYsGFZWVnp6ek9evQYMWKEa+9eAgCs5u56UkaKOmiVXqqJia3clIVOBaGfn9/NOh/atm3btm1bl5YEALCuduFSRoqSlKaXVIj/aueOLHTlWqMAANRcmzBp8zAlaZVeYYhnO9R6FvKEegCAx2kRKq1LUd771vhdXq3PnSEIAQCe6I4QafNw9dMjxsxsvVYvRBACADxUw0CRmaJ+fswev+Lahzpc1sV/bdHLXRGRBCEAwHM1CBRbhquHiu1dPtOuNttfMcS4TL1jhOTvimYFghAA4NEi/cX+sbbTZfa7lmmGXVwxxNgMPaWZ9Hhb10QYQQgA8HThfuLQONv5cnvXL5XUTLsLU1AQhAAArxCiih1jbAXlkhB2F6agIAgBAF7hiiEmbdRf6mTcESpN3+bKeaQEIQDA0129LzillfFOD6lUEy7MQoIQAODRrpkdIwnxx16KC7OQIAQAeLQTpfbxLX82O8aRhVGBkkv6CFlrFADg0aJDpehQ6ZqdkhC/6mj59olnnnnmyJEjZldhUbm5ub/5zW/MrsK6Zs2alZ2dbXYVFnXs2LGnnnrK7Cqsa968eWvWrHHtOb04CLOzs0tKSsyuwqLOnTu3a9cus6uwrl27dp07d87sKiyqpKSEv0JMtH///pMnT7r2nF4chAAA1BxBCACwNMlut9/+KLc4d+7csGHDzpw54+TxhYWFoaGhqsp8HxNUVFRcunQpLCzM7EIsqqioKDAw0M/Pz+xCrEjX9YsXL0ZERJhdiEUVFxfbbLaAgAAnj1+5cmXbtm1vfYwHBaEQ4uzZs8XFxWZXAQDwEU2bNr3tn4yeFYQAALgZ9wgBAJZGEAIALI0gBABYGkEIALA0rwzC119/fdCgQTExMStXrjS7FsvJysrq379/aGhoWFjYvffee/r0abMrspYnn3yyYcOGsiy3aNHi7bffNrsci5o3b96dd965bds2swuxltdff/3OSsrKylx1Zq8MQiHEo48+KssyS6y538WLF59++ukTJ07k5+drmvbEE0+YXZG1TJgwYe/evYZhLFy48KWXXtq0aZPZFVnOwYMHP/7448LCQhe+EcMZ58+fHzRoUPp/ON9KeFte2Y3+wgsvCCFef/11swuxopSUlKvbkydPfvrpp00sxoJ69+7t2OjZs2ebNm2OHTtmbj1WYxjGlClT3nnnnXvvvdfsWqwoIiKiZcuWLj+tt44I4QnS0tJ69OhhdhWWs2PHjk8//fS///u/7Xb7iBEjzC7HWt59990OHTr06dPH7EIs6i9/+UujRo169OixZMkSF57WK0eE8ARLly5dtmxZTk6O2YVYTl5eXlpa2s6dOxMSElhlzZ3y8/Pff//9b775xuxCLCo1NXXSpEn16tVbt27dxIkTIyMjExMTXXJmRoSojrS0tCeffHL16tVNmzY1uxbLmTRp0uLFi7/99tvdu3e/9957ZpdjIc8888yECRMOHz6ck5OjadrBgwdPnTpldlEWcs8997Rr165Bgwb33XffQw895MJBISNCVFlGRsakSZOWL1/epUsXs2uxLlVVO3XqdPz4cbMLsZCGDRtu3rx58+bNQoiSkpL58+cHBwc/8MADZtdlRZLkyvVBvTIIDx06dPHixUuXLh09ejQnJ6d169ahoaFmF2UVWVlZo0aNeu2111RVzcnJkSTp7rvvNrsoq7h06dK//vWvpKSkwMDArVu3Llq0aOHChWYXZSHz58+/uh0VFfXmm2/269fPxHqs5sMPP+zXr1+dOnU2bNjwySeffPHFF646s1cuuv3LX/5y69atV7/885//3K1bNxPrsZQ///nPf/nLX65+6efnt2XLFhPrsZSysrKHH354+/btZWVl0dHRzzzzzH333Wd2URaVnJw8Z84c/gp0p0ceeWTTpk2lpaUtW7acNm1aamqqq87slUEIAICrMFkGAGBpBCEAwNIIQgCApRGEAABLIwgBAJZGEAIALI0gBABYGkEIALA0ghAAYGkEIQDA0ghCAIClEYQAAEv7f6yzZp7pJSQFAAAAAElFTkSuQmCC", "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" ], "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" ] }, "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.6.1" }, "kernelspec": { "name": "julia-1.6", "display_name": "Julia 1.6.1", "language": "julia" } }, "nbformat": 4 }