{ "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.846687520186 -0.70 0.80 4.2\n", " 2 -7.852281958897 -2.25 -1.54 0.80 1.0\n", " 3 -7.852621379818 -3.47 -2.53 0.80 3.2\n", " 4 -7.852621945757 -6.25 -3.36 0.80 2.2\n", " 5 -7.852621955393 -8.02 -4.26 0.80 1.8\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+gvaeTAAAgAElEQVR4nO3deWAU5f348WdmnpyQcJOECCFALEIgkgMIBElA5ExCuL5+tcYT0Naf6U/rF/RrlZ+VCra1pmotqC1FaBXJgQlnuOQKh4Qb5LZcCogokHtm9vfHWlRECGGzs7vzfv2VrJvdT0qad57ZZ2cUh8MhAACwK9XqAQAAsBIhBADYGiEEANgaIQQA2BohBADYGiEEANgaIQQA2BohBADYGiEEANiam0JYXFz81FNPFRQUuOfpAACoI3eEsKCg4PXXXx83btysWbPmzp37U3c7ceLE6dOn6/6whmG4YjrUh2manJzPQvzwW8s0TatHsK+G+OXjjhDOnDnz+eef79Wr14svvjhz5syfuttrr702Z86cuj9sRUWFK6ZDfdTU1PC72EL88FvI4XDwv7+FqqqqXP6HiDtCePjw4ZiYGCFETEzM4cOH3fCMAADUkVs3yyiKwiEFAIBHudkQnjhxYuXKlVe8tldZWVlQUDB37tyzZ88KIaKioj777DMhxNGjR9u3b3+TzwgAgAvdVAi7d+/epUuXYcOGrV69+vKNly5d6tWr15tvvrlkyZIuXbrs27fv/vvvnzZt2sGDB6dOnfrAAw/c5MQAALiQvJkvXrBgQVRUVHx8/Pdv/Mc//hEaGrps2TJVVX/961+//PLLs2fPNk1z6tSpaWlpjzzyyDUesLKy8vz5886P/f39GzVqdMUdzlaJVoFXflW5LhQhgm/qWwEA2NRN1SM6OvrHNxYXF48ePVpVVSHEmDFjhg0bJoTIzs7Ozs6+9qPV1tb+7ne/++Mf/+j8dOjQoTNmzLjiPvev9R8TZYxtZwghLl26JISoMJT/Wuv3ZGc9LZxXH92nqqpKSiklf31Yw/nDD0s4d43y9iGrVFZW1tTUaJpWx/sHBwdf986u/0V28uTJyMhI58eRkZHnz58vLy//8drux/z8/H77298+9dRT17jP/LtEVonuF6De10kVQmhBIT8v0cffpmZ0CnDJ8KgjPz8/QmitkJAQq0ewKYfDoapq48aNrR7EpjRNCwgIqHsI68L1v8gMw3AuB4UQzlld+IazYCkKBsmsEl0IMbilcm+J/kCMem8nThQHAKgn1yckIiLizJkzzo9Pnz7duHHj0NBQFz6+s4V/P2CmLPN76FYqCAC4Ka6vSGpq6tKlS50fL1u2LDU11eVP4RDCdIizVerkLcamMxypBwDU300dGs3Nzd27d++JEyfefvvtlStXPvXUU7feeuuECRPeeOONxx9/PCIi4ve//31xcbGrZnWq0MXIEn38z9RBLSuGrAoYtlQf3lb9U7LWglcJAbjRRx99dENnhcRNGjp06IMPPtgQj3xTIYyJiQkKCkpISHB+6nz1uHXr1lu3bn3vvfcqKytXr159++23u2DM/6jQReZ/Xhe8eNGxLt0vfal+usrRPU9/OUnNjuEwKQA3Wb9+fUBAQEZGhtWD2MK6deuWL1/uiSF0vjXixyIjIydPnnwzj/xT/mulPv5n6rgO3wYvWIqiwTKzRP/f29U/7zFnHTDfStF+1kRpiKcGgCt069Zt7NixVk9hC7W1tQsXLmygB/ey7e+z+ssrDoEGS7FgkFSEmNBZfXOv2a9I/0UX9dnbNX8WhwCAOvCyXFz1hcBgKYKkkKrIiVU3ZcrNZx3d8vRVn7OJBgBwfV4WwuuKDlEWDZZ/6KU+8LGRvdo4V231QAAAz+ZrIXRKb6fuHCWbBYjuefrsgyZrQwDAT/HNEAohmviL3GSt6C7t9T3mgIX6p19TQwAe4UKtMH76F5LpEN/UuHEa+HAIneJbKqUZcmSUekexPqXMqHbZud4AoJ5GL9ezVxtXbaHpEOPXGsOX6m4fytZ8PIRCfLuJZluW3H1edMvXV5xiaQjAStOTtCUnzHtWGfoPL5njEOKXG4wPjpgvJdbzjNIREREXL150wYhC6Lq+fPnym7/IxqVLl1q2bHntx+nSpcuxY8du8oluhu+H0CmykTJ/oPbHXurDa4zs1caXVVYPBMCu4lsqJUPl8pPmvau/a6FDiF+sN947aBYPlqkR9XwzdHV1tauuD1VRUTFo0KCbv2SCv7//888/ryjX+o5cOHb92CWETunt1D1jZJtGosv82pmfsokGgDWuaKFLKngFh8Mxd+7cyZMnv/vuuzU1373quHTp0hdeeOGtt9765ptvnLesWLFi3759zo+PHTv20UcfCSEKCwuFEO+8887MmTMvX0fhui5duvTqq68+99xzixYtct6iKEpg4LeXU8/Pzz9x4sSsWbOef/75tWvXuuK7dA17hVAI0UiKaUna0qHy3f1marG+j000AKzw/RY+us7FFRRCTJw4cdasWZ06dZo/f/7o0aOdNz7//PNPPvlkeHj41q1be/To4Wzh22+//fHHHzvvsHPnzldeeaV+z1hRUZGSknLo0KG2bdv+6le/evXVV4UQlZWVEydONE1TCDF16tS77777yJEjjRs3zsjIWL9+vQu+T1fwsjPLuEqPFkpphnxnv9mvSH+ks/r/4rUAV17lEYAd6aYYeoP7XDqGKvOPmkKIHi2UqduNqdtv4GsLB8lGP/Er/NChQ/n5+ceOHQsODr7vvvuioqI2b94cExPz+9//fteuXZ06dRJCDB48eMaMGf/zP/9z1UcYOXKkEOKRRx6p+/W3586d26ZNm7/85S9CiKSkpAEDBjzxxBNX3Cc9PX3SpElCiDNnzixYsKBv3751fPAGZdMQCiFURUzorI5op/yq1IzN09/qq90ZyUlKAdSfqohJ3W/gb2qHQ7y62/RTHaoimvorT3dT1Rv5JRTw00f0du3a1a1bt+DgYCFEQEBAjx49du/eLYRo1qyZs4JCiH79+jlvvCGnTp1yliwwMPDtt9/+/n/auXNn7969nR/36NGjurr6+PHjLVq0+P59Ll+GITIycufOnTf67A3EviF0ahOszBuoFR0zH1lr9GylvNlXaxVo9UwAvJOqiLr/Pe18XXDtF+aSITLUTwxarL97QMxN1aQrXrC6dOlSVdV3ewIrKyuDg4ODg4MrKysv31hRUeEspaZplzfFVFRUXPuRQ0NDnQdaf7xSLC8vd16DSAih63ptba3z8b9P0zzx4JvtXiO8qvR26t4xskOo6MomGgAN74rdMVfdR3qTtm7deuDAASHE0aNHt27d2qdPn5iYmODg4KKiIiFEeXl5Xl6e88Lp7du3LysrE0I4HI78/Hznlzdq1MjPz+/cuXNXPGzjxo1Hjhw5cuTIESNG/PhJFyxY4GztvHnzYmJiwsLCXPPNNDBC+K1gKaYlaSXD5N8OmP2L9b1sogHQMK66R9TlLYyLi7v33nszMjJ69+794osvtmvXLiAg4G9/+9vEiROHDh0aGxvbp0+fcePGCSHGjx+/ZMmStLS0hISEy+s8TdN+8YtfxMXFJSYm7tmzp45P2qFDh+Tk5OHDhz/55JMzZsxwwbfhFnY/NHqFuObKhnT5zn7zjiL9kc7qlHgt0BPX8QC82Pi1xgdHzCVDZEr4D46jxrdUFg+Rdy3WH1gj5qTe7K+esLCwvLy87du3R0VFhYeHO28cPHjwoUOH9u3bFx4eHhkZ6byxffv2+/fv37t3b1RUVPPmzS8fU33ttddeeeWV8vLykJCQOj5pamrqgw8+eOjQoS5duji/KiQk5PTp06qqCiFWrVoVFBTkvOdjjz12829SdBVCeCXnJpr0duqkzUZsnv6XvtpdbKIB4Dp9w5QHYq6soFPPVsqyoXLbOdcckQoICOjVq9cVNwYHByckJFxxY+PGjXv27On82M/P7/Lt/v7+/v7+N/SkLVq0+P4GGUVRWrdu7fw4NDT08u2X31zoCQjh1UUEi9mp2sLjjonrjKSWyht9tNZBVs8EwCc8eOu1XpPq2Urp2aqef3y3bdvWufaKjY11/3pr2LBhl5t3QyIjI+v+Jo2GQAivZXhbJW20fGWn0S2/9tk47f90vbHNzQDgTjt27HB+0KNHjx49erj52ceMGVO/L1yzZo1rJ7lRbJa5jmAppsRrJUPl+0fM/gv1PefZRAMAPoUQ1kn35sqGDDn+Z+rARXpOqVHONVIAwFcQwrpShMiOUbeP8jtfLbrn6UtPsDQEAF/Aa4Q3JjxIzE7VVn3ueHSdcWsT8de+WmQjXjYEAC9GCOsjLULZniWn7zTiC3U20QD21KFDh9/+9rfz5s2zehBb+Oqrr1JSUhrowQlhPQVJMSVeG91efXS9Me+oOSNFi21GDAEbGT9+fGJiotVT2Ejbtm0b6JEJ4U3p1lxZly7fO2gOWqSP66BOTdQa+13/qwD4AFVVf/zOdHgjNsvcrO9voonL1xcfZxMNAHgTQugaYUFidqr27h3a/91opC/TT5STQwDwDoTQlVIjlG1ZMqGlklCo5+42DWoIAB6PELqYcxPNuhHyo2NmzwX6J18SQwDwaISwQcQ0UZYPkzld1fSlek6pcbHW6oEAAD+BEDYU5yaaPWP8qgxx23w976iLLjsNAHApQtiwmgeIGSna3FTtua1m+jL9OJtoAMDDEEJ36B+h7BwlU8LU2/P16TvYRAMAHoQQuomfKibFqZsy5fJTZlKhvuUsMQQAj0AI3apTqLJsqPxVrJq+TJ+4jk00AGA9Quhuzk00+8b4CSFum6/PZxMNAFiKEFqjWYCYkaL9K017YauZvkw/dokjpQBgDUJopX7hyvZRMiVM7VHAJhoAsAYhtJhzE83mTLnylJlQoG86QwwBwK0IoUfoGKosHSr/t4eaWaJPXGdcYBMNALgLIfQgY6PVfWP8AjXRZb4++yCbaADAHQihZ2kWIHKTtffTtFd2miOW6v9mEw0ANDBC6IlSwpVtWXJQpJpUqE8pM2pYHAJAgyGEHspPFTmx6sZMufGMI6lQ38gmGgBoGITQo3UIUZYMkS8lquNWGNmrja+qrR4IAHwOIfQC6e3UnaNlswDRPZ9NNADgYoTQOzT1F7nJ2gcDtN/vNIcv1T+7yJFSAHANQuhN+oYp27LkXZFqzwVsogEA1yCEXkaqIidW3ZQpN591JBbqG06zNASAm0IIvVJ0iLJosJyaqN690shebZxjEw0A1Bch9GLp7dRdzk00eWyiAYB6IoTerYm/yE3Wiu7SXt9jpi3U93/DkVIAuDGE0BfEt1RKM+TIKLVfkT6lzKg2rB4IALwHIfQRzk0027Lkrq9E93x95SmWhgBQJ4TQp0Q2UvLu1P7QS31wjZG92viyyuqBAMDjEUIflN5O3TtGtmkkusyvnfmpydoQAK6BEPqmRlJMS9KWDpXv7DfTFuqffk0NAeDqCKEv69FC2ZAus6LUvkX65C1sogGAqyCEPs65iWbXaHn4guiWry8/ydIQAH6AENpCm2Dlw4HaH3upj6w1slcbZ9lEAwD/QQht5PImmq5sogGA/yCE9hIsxbQkbdlQ+e5+M7VY38smGgC2Rwjt6PYWSmmGvLeTekeRPnmLUcUmGgA2RghtSlXEhM7qrtF+p8pFtzy9hE00AOyKENpaRLCYnar9qbc2fq0xboVxptLqgQDA7QghxIh2yt4xsksz0S2/Nne3yS4aALZCCCGEEMFSTInXlg2V7x8x+y/U95wnhgDsghDiO3HNlfXp8r5Oav9iPafUKNetHggAGh4hxA84N9HsHuN3vlrE5etLT7A0BODjCCGuIjxIzE7V/pysPbreGLfCOM0mGgC+ixDiJw1rq+wdLbs0E93za3N3m8cvOWYfNK+4j0OI13abFRxEBeC1pNUDwKMFSTElXhvdXp24zph3VARqyoVa8/Eu3/795BDiV6WGpohgfo4AeC1+geH6ujVX1mfI9w6aT282/n1RqTbEL2O+q+CrvTWrBwSA+iOEqBNFiOwY9a5b1F9vMl7Yah44r/prDj9NoYIAvB0hxA0IDxJzUrWi444xJY5WQeLAWCoIwOuxWQY3xiHE8hPmg53MAEV0mKfv/Ir3VwDwboQQN+Dy64Kv9XTsGqXcEiT6Fum5u6/cSgoAXoQQoq6u2B3jr4oNmbJ3K/UPu8wxK4yva6yeDwDqhRCirj676Aj1/8EeUX9VLByiDW+nhAWKHgX6xjMcJgXgfdgsg7qKDlF+m3Dl7hh/Vfy1ryaEKPjMzCrRJ8dpT8SqihXjAUD9sCKEa2S1Vzdmyg+OmFklxlfVVk8DAHVGCOEyUY2VNSPk7S1EfIG+/jSHSQF4B0IIV5KqmBKvvd5HHb1cn1JmcI1fAJ6PEML10tupW0bKFaccg5foX3DlCgCejRCiQbRtpKwaJvuGKfEFtSUnWRgC8FyEEA3FeZh0Tqp8aI0xpcwwqCEAj0QI0bAGtFG2jpQbzzjuXKSfqiCGADwOIUSDax0kFg+RI6PUxEJ98XFaCMCzEEK4gyJETqz6rzQ5YZ2RU2rUcnZSAB6DEMJ9+kco27Lk4QuOlCL96EWWhgA8AiGEW7UMFEWD5T0d1V4L9A+PsjAEYD1CCHdzHiZdOFg+s8XMKTVqqCEASxFCWCOplVKWJb+oFH0+0g9d4DApAMsQQlgm1E98MEB7oquaUqT/6zALQwDW4DJMsFh2jBrbTLl7lbH4uOOvKVowP5IA3IsVIawX31LZOlKaQiQV6rvPc5gUgFsRQniEED8xJ1WbFKemLdRnfsphUgDuQwjhQbJj1DUj5Bt7zezVxqVaq6cBYA+EEJ7ltqbKpgzZLEAkFuo7vuIwKYAGRwjhcYKkyE3WXkpUBy/Wc3dzmBRAwyKE8FBjotU1I+Ssg+bo5cbXNVZPA8B3EUJ4rlubKBsz5C2NRI8CfeMZDpMCaBCEEB4tQBO5ydqfequZJfr0HSYxBOByhBBeYGSUuiVTfnTMHFlifFVt9TQAfAshhHdo11j5eLjs0ULEF+jrvmBlCMBlCCG8hlTFlHjt7X7a3auMKWUGx0kBuAQhhJcZFKlszNBWnnLctVj/otLqaQB4P0II73NLI2XVcJkSrsQX1C47ycIQwE0hhPBKmiKmxGtz0+RDa4zJWwyDGgKoL0IIL5YWoWzLkjvOOQYu0k+WE0MA9UEI4d1aBYpFQ2RWlJq0QF90nBYCuGGEEF5PESInVi24Uz6+wcgpNWo4OymAG0EI4SN6tVY2Z8rDFxwpRfqRiywNAdQVIYTvaBkoigbLezuqvRfoHx5lYQigTgghfIrzMOmiIfKZLebEdUa1YfVAADweIYQPSmyplGXJb2pEnyL90AUOkwK4FkII3xTqJ94foOV0VVOK9H8e5jApgJ9ECOHLsmPUFcPk77ab2auNCt3qaQB4JEIIH9e1mVKaIU0hEgv1XV9xmBTAlQghfF+In5iTqk2OUwcs0nN3c5gUwA8QQthFdoy6doT82wEze7VxqdbqaQB4DEIIG+ncVNmUKZsFiMRCffs5DpMCEIIQwm4CNZGbrE1NVIcs4TApACEIIexpdLS6MUP+64g5arlxvtrqaQBYihDCptqHKB8Pl20biR4FeukZDpMC9kUIYV8BmshN1nKT1ZEl+pQyw6SGgC0RQthdZpS6JVMuP+kYWWKc4zApYD+EEBDtGiurh8v4liK+QF/7BQtDwF4IISCEEFIVU+K1d/tp/73K4DApYCuEEPjOnZHK1pFyw2nHoMX65xVWTwPALQgh8ANhQWLxENkvXIkvqF16goUh4PsIIXAlTRFT4rV/DpAPrzVySo1a3nYP+DRCCFxdWoSyLUsevOC4c5F+spylIeCzCCHwk1oFioWD5aj2atICfeFxWgj4JkIIXIsiRE6sWjhIPrHByCk1ajhMCvgcQghcX89WSlmWPFUh+n6kH7nI0hDwKYQQqJMm/mLeQO3nndTeC/QPjrAwBHwHIQTqynmYdPEQ+dwnZvZqo1K3eiAArkAIgRuT0FLZmiVrHaJvkX7wGw6TAl6PEAI3LNRP/CtN+1Ws2q9Yn3uIw6SAdyOEQD1lx6grh8tpO8zs1UY5h0kBr0UIgfrr0lTZnCkDpUgs1Hd9xWFSwCsRQuCmBEkxM0V7Jk4dsEjP3c1hUsD7EELABbJj1LUj5N8PmGNXGN/UWD0NgBtBCAHX6NxU2Zgp2wSLngv0bec4TAp4DUIIuEygJnKTtZeT1GFL9NzdXNwX8A6EEHCxUe3VjZny/SPmqBLjfLXV0wC4HkIIuF5UY2XNCPmzpqJHgb7hNCtDwKMRQqBB+KliWpL252Q1a7k+pczgOCngsQgh0IAyotRPRsrlJx2Dl+inK62eBsDVEEKgYbVtpKweLvuGKT0X6Gu+YGEIeBxCCDQ4qYop8drf7tDuWWVMKTMMagh4EkIIuMnANkrZSFl62jFokX6qghgCnoIQAu7TOkgsHiIH36ImFOhLTtBCwCMQQsCtVEVMilP/NUCOX2vklBq1nJ0UsBohBCyQGqGUZclDFxz9ivXPLrI0BKxECAFrtAoUxYPlf3dQ+xTpxcdoIWAZQghYRhEiJ1ZdMEjmlBo5pUYNh0kBKxBCwGJJrZSyLPl5hej7kX74AktDwN0IIWC9Jv5i3kBtfGc1+SP9/cMsDAG3IoSAp5jQWV0yRP5mq5m92qjQrZ4GsA1CCHiQ+JZKWZbUHaLnAn3PeQ6TAu5ACAHPEuIn/pmm/U93deAifc4hDpMCDY4QAp4oO0ZdNVxO32FmrzbKOUwKNCRCCHio25oqmzNlswCRWKjv/IrDpEBDIYSA5wqSIjdZeyZOHbhIzywx3t1/5ZHSZScdOaWGJbMBPoMQAp4uO0Zdly4/u+j4TZn56q7vWrjspOOFrcbz8ZqFswE+gBACXuBnTZRNmXJUlPKbreakzYb4TwWLB8sWAVYPB3g5afUAAOokUBNv9NH6hJkPrTHWfxFgCCoIuAYhBLzJPR3VCl08vsG4pZE4XeloEaBYPRHg9Tg0CniTZScd7+43d42oDvUXyQv0nFLjUq3VMwFejhACXuPy64LhgY6NGTI5TN10xhGXry8+zpsrgPojhIB3uGJ3jL8qPrpLCwtWMqKUJzcZ6cv04+XkEKgPQgh4h2b+4ordMf6qmDdAy4xSt2XJhJbK7fn69B2mQQ2BG+SOEG7evHnUqFEJCQnjxo07fPiwG54R8D1JrZQf7xEN0ERqhBKoiSnx2qZMueKUmVSobz5LDIEb4I4Qmqb54osvfvLJJxkZGQ8//LAbnhGwoU6hyrKh8pnb1Yxl+sR1xgU20QB1444Q9u7dOzY2VlGUvn37njlzxg3PCNjW2Gh13xi/QE10ma/PPsjFK4Drc+X7CKuqqk6fPv39WwIDA8PCwpwfOxyOZ5555oknnnDhMwL4sWYBIjdZGxvteGy9Me+I+UYfrX0IbzcEfpIrQ3jw4MEXXnjh+7fcdtttU6dOdX7861//umXLlo8++qgLnxHAT0kJV7ZlyTf3mj0X6L/ooj57u+bP3jjgauoUwpqamk8++WTbtm0VFRVPP/309//TvHnzli5d2rp168cff7xbt275+flXfYTnnnvu/Pnz77zzjgtGBlA3UhU5sWpGlPL4BqNbnv5WX21AG5aGwJXq9CdiaWnpI488Ulxc/Oyzz37/9rfeemvSpEmpqamXLl1KSUmprKy86pe/++67b7zxRsuWLZ999tkrlowAGlp0iLJwsPxDL/WhNUb2auPLKqsHAjyM4nDUdaf1rl274uPja2u/3YtmmmZMTMxrr72Wnp4uhOjZs+cvf/nL+++//8dfePjw4aNHjzo/1jQtLS3tqo//1FNPlZeXjxkzxvlp69atu3fvfo15Ll68GBISUsfh4VpVVVVSSik5V6016vfDX66L324z/rbffClRG99ZZW1YPw6Ho7y8vHHjxlYPYlMVFRUBAQGa5sqrj9X/F9kXX3xx5MiRAQMGOD9NS0tbv379VUPYsWPHjh07XvcBDcNYvnz5gQMHnJ/26dOnQ4cO17h/eXm5ovD/ZWsQQmvV+4f/udtERriS84mcvV+8lqR3DuUdhzfM4XBUVFRYPYV9VVRU1NbW1j2EwcHBqnqdY583FcLAwMBGjRo5P23VqtXevXvr/WhCCCnlY4899tRTT9Xx/g6Hgz/KrCL/w+pBbOpmfvj7NBabbhHv7DeHrFAf6axOidcCubLvjXA4HIqi8MvHKqqqunxFWP9tZIGBgbquXz6yWl1dHRQU5KKpADQgVRETOqu7RvudKhexefrSE6wLYWv1D2GbNm0Mw/j888+dnx4/fjwyMtJFUwFocBHBYnaqNjNFe6LUSF+mn+Sc3bCr+oewadOmaWlpc+fOFUJcuHChuLg4KyvLdYMBcIcBbZTtWTKhpRJfqOfu5pzdsKM6vcZz7ty5nj171tTUGIbRsWPHsLCwDRs2CCF+97vfpaenr127dv/+/f3797/jjjsaeFoArhckxZR47Z6O6i/WG7MPmTNStMSWbEODjdQphE2bNi0pKfnua/6zRaJXr1779+/ftGlTWFhYjx49GmRAAG5xaxOlZJh876CZsUwfG62+lKiF+Fk9E+AWdQqhpmk/9U6GZs2aDRkyxKUjAbCGIkR2jDqinfrMFuO2+fprvdUx0ZyWDb6Pn3IAP9A8QMxI0f6Vpr2w1Uxfpv/7Ei8bwscRQgBX0S9c2T5K3tlGTSrUp5QZNVzQCb6LEAK4Oj9V5MSqGzPlpjOOpEK99AxLQ/gmQgjgWjqEKIuHyJcS1btXGtmrjXPVVg8EuBohBHB96e3UnaNkswDRPY8L38PXEEIAddLEX+Qma0V3aW/sNdMW6p9+zZFS+AhCCOAGxLdUNqTLkVFq/4X6lDKj2rB6IOCmEUIAN8Z54fttWfLIBRGbp5ecZGkI70YIAdRHm2Bldqr2p97a+LXGuBXGmUqrBwLqixACqL8R7ZS9Y2SXZqJbfm3ubtNkcQgvRAgB3JRgKabEayVD5QdHzDuK9d3niSG8DCEE4ALdmyvrM+SEzuqgRXpOqXGp1uqBgDojhABcw3nO7u2j/M5Xi9vm600gXyMAABGySURBVPmf8XZDeAdCCMCVwoLE7FRtTqr23Cdm+jL9OBe+h8cjhABcr3+EsmOUTAlTb8/Xp+/gwvfwaIQQQIPwU8WkOHVTplx5ykws1Ddxzm54KkIIoAF1ClWWDpXP3q6OLNEnrjMusIkGnocQAmhwY6PVfWP9AjXRZT7n7IbHIYQA3KGpv8hN1j4YoP1hlzlwkX7gG46UwlMQQgDu0zdMKRspM9qpKUWcsxueghACcKvL5+zefV50z9dXnGJpCIsRQgAWiGykzB+o/aGX+vAaI3u1cbbK6oFgY4QQgGXS26l7x8g2jUTX+bUzP+WU3bAGIQRgpWAppiVpJcPk3w+Y/Yv1vVz4Hm5HCAFYL665sj5d/ryT2r9Yzyk1ynWrB4KdEEIAHkFVxITO6q7RfuerRVy+vuQES0O4ibR6AAD4TniQmJ2qrfrc8dg6I6aJeKuvdksjxeqh4ONYEQLwOGkRyrYsmdBSSSjUc3dzzm40LEIIwBMFSTElXls3QhYfN3su0LecJYZoKIQQgOeKaaIsGypzuqqZJXpOqXGRc3ajARBCAB7NeeH7PaP9hOCc3WgQhBCAF2gWIHKTtX+maa/sNNOX6f++xJFSuAwhBOA1+oUr27LknW3UpEJ9SplRw+IQrkAIAXgTP1XkxKobM+WmM46kQr2UC9/jphFCAN6nQ4iyeIh8KVG9e6WRvdo4V231QPBmhBCAt0pvp+4cJZsFiO55+uyDnLMb9UQIAXixJv4iN1krHqy9sddMW6jv45zduHGEEIDX69FC2Zgh7+mo9ivSJ28xqrjwPW4EIQTgC5zn7N45Wp4qF93y9GUnWRqirgghAN/RJliZnar9qbc2cZ0xboVxptLqgeANCCEAXzOinbJntOzSTHTLr83dbbKLBtdGCAH4oGAppsRry4fKeUfNO4r1XV8RQ/wkQgjAZ3VrrqxLlxM6q3ct1nNKjUucsxtXQwgB+LJvz9k9xq/KELfN1/OOclo2XIkQAvB9zQPEjBRtTqr2m61m+jL9GOfsxvcQQgB20T9C2TFKpoSpPQr06Tu48D2+RQgB2IifKibFqZsz5cpTZkKBvpFzdoMQArChjqHK0qHyf3uoWSV69mrjK87ZbW+EEIBNjY1W9431axYguudz4XtbI4QA7Kupv8hN1hYM0l7fYw5YqO//hiOldkQIAdhdQkulNENmRqn9ivQpZUY15+y2GUIIAEKqIidW3ZYl95wX3fL15Zyz204IIQB8K7KR8uFA7Y+91EfWGtmrjbNVVg8EtyCEAPAD6e3UvWNkm0ai63zO2W0LhBAArhQsxbQkrWSYfP+ImbpQ333e8fCaqywQp+8w135BJ70eIQSAq4trrmzIkI/8TL1zkX6hVgxZrH/xvQsc/mGXue2cI7m1Yt2AcA1CCAA/yXnO7u2j/II08UWluKPo2xb+YZf5yVnHnFRN8kvU+0mrBwAATxceJGanaktPOB5ea8Tm6fdE+52toYK+g39GAKiTwbcoB8fKqMbKW5+qKeEKFfQZ/EsCQF29uc+8tYnyepI+abMxcLF+spydMr6AEAJAnThfF3yvv3pfR6P4LnngvLi9UJ/5KScp9XqEEACu74rdMWltlH8O0MIClD/vNceu4K333o0QAsB1OIRoLMXctB/sjukXrrzeV5vdX+sYKuLya/OOsjT0VuwaBYDrUIR49LarLBvSIhQhRHwLbWSU+uAa44MjjrdStBYBbp8PN4cVIQDcrN6tlW1ZskOo6J6nF/6bpaGXIYQA4AKBmpiWpH04UJu02Ry3gqveexNCCAAu0ydM2Z4lO4SKbnl60TGWht6BEAKAKwVJMS1J+2CA9uRGM3u1cbHW6oFwPYQQAFwvJVzZmiWDpOier684xfvuPRohBIAGEeonZqRoM1K0h9YYE9cZl1gaeipCCAAN6K5IZecoKYTonq+v+pyloScihADQsJr4ixkp2pt9tPtXGxPXGeW61QPhhwghALjD0LbKztFSCBGXr6/huvaehBACgJs09RczUrTcZO2eVUZOqVFtWD0QhBCEEADcbHhbZVuWPFUh4gv0LWdZGlqPEAKAu7UKFB8O1KYkqBnL9MlbWBpajBACgDXGRqvbR/kd+EYkFuplX7I0tAwhBADLhAWJ/Du15+PVYUv1yVuMGk7KZgVCCAAWcy4N930tkgr17edYGrobIQQA64UHiQWDtOd6qIOX6FPKDIMauhEhBABPMTZa3ZIp15929C3SP/2aGLoJIQQAD9KusbJsqHzoVvWOYn36DpOloRsQQgDwLIoQEzqrmzPl0hNmvyL9wDfEsGERQgDwRO1DlBXD5QO3qn2L9Ok7TJMaNhhCCAAeyrk03JQpF58w7yjWD10ghg2CEAKAR+sQoqwcJrNj1OSP9NzdrAxdjxACgKdTFTGhs/rxCDn3sDl4sX68nBq6EiEEAO/QpamyIV0ObKMmFOgzP+UkNC5DCAHAa0hVTIpTVw2XMz81hy7RT7A0dAVCCABepmszZWOGTI1QEwpZGroAIQQA7+NcGpYMlX/dZw5fqp+qYGlYf4QQALxV9+bKpkx5R7jao0B/7xBLw3oihADgxfxUMSlOLbpLvrzdHLfC+LLK6oG8ECEEAK/Xs5VSliU7hIru+bX5n7E0vDGEEAB8QaAmpiVpeXfKZ7eY41YY56qtHsh7EEIA8B3JrZVtzqVhnr7g3ywN64QQAoBPCZJiWpI2b6D29GZz3ArjPEvD6yGEAOCD+oYp27Nkh1DRLV8vPsabK66FEAKAbwqWYlqS9vc7tMc3GBPXGRdrrR7IUxFCAPBlgyKVnaOlEKJ7vr7yFEvDqyCEAODjQv3EjBTtr321Bz42Jq4zLrE0/CFCCAC2MPgWZddoKYSIy9c//pyl4XcIIQDYRRN/MSNFe72P9vPVxsR1RoVu9UCegRACgL0Ma/vt0rB7vr72C5aGhBAA7Kepv5iRor2cpI5boU/eYlQbVg9kKUIIADY1NlrdMcrv4DcioVD/5Ev7Lg0JIQDYV+sgkXen9kK8OmKpPnmLUWPLk7IRQgCwO+fScP/XIrFQ33bOdktDQggAEGFBomCQ9pse6pAl+uQtRq2dloaEEADwrbHR6vYsv73nRVKhvuMruywNCSEA4DsRweKju7Qnu6mDFunTd5iGDWpICAEAV8qOUbeMlCUnzZQiff83Ph5DQggAuIqoxkrJMPngrWq/In36DtP03RoSQgDA1SlCTOisbsqUS06Y/Yr1gz66NCSEAIBriQ5RVg6X98eofXx0aUgIAQDX4VwabsyQi46bqQv1wxd8KoaEEABQJx1DlVXD5ej2avJH+sxPfWdlSAgBAHWlKiInVl01XL79qTl0iX683BdqSAgBADemazOlNEOmRagJBfrMT73+JDSEEABww6QqJsWpK4fLGZ+aw5bqJ715aUgIAQD1FNtM2ZQh+4er8YVevDQkhACA+nMuDUuGyrf2menL9M8rrB7oxhFCAMDN6t5c2ZwpE1oqCYW1eUe9bGlICAEALuCniinxWsGd8rmt5rgVxpdVVg9UZ4QQAOAyvVor27Jkh1DRPb+28N/esTQkhAAAVwrUxLQkbf5AOWmzOW6F8VW11QNdDyEEALhenzBle5bsECq65elFxzx6aUgIAQANIkiKaUnaBwO0Jzea41YY5z11aUgIAQANKCVc2TFKRgSL+EJ9xSlPfN89IQQANKxgKXKTtRkp2kNrjInrjEu1Vg/0Q4QQAOAOd0UqO0dJIUT3fH3V5x60NCSEAAA3aeIvZqRof+mr3b/amLjOKNetHkgIQQgBAG425BZl52gphOiep6/5wvqlISEEALhbU38xI0X7cx/tnlXGxHVGhaVLQ0IIALDG8LbKrlGyQhdJhfqWs5YtDQkhAMAyzQLEe6nalAQ1Y5k+eYtRbVgwAyEEAFhsbLS6Y5TfgW9EYqG+9Ut3Lw0JIQDAeq2DRP6d2vPx6vCl+uQtRo0bT8pGCAEAnmJstLp9lN++r0VSob79nJuWhoQQAOBBwoPEgkHacz3UwUv0KWWG0fA1JIQAAI8zNlrdkinXn3b0LdIXHjenbr/yUGmlLn65wTVHUAkhAMATtWusLBsqH7pVffBjo/iYOWnLdztKqw3xXyuNnq0Uf1dEjBACADyUIsSEzurmTBmkiX8cMCesM4QQ1YYYu8IYHa3cH+OahEmXPAoAAA2kfYiyYric+an5f0uNvV/JpgGOsR1UV1VQsCIEAHg+RYiJndUNGXLPN4pUHC6soCCEAACvUG2I57ea0+ONiGBl8hZXnoGGEAIAPN3l1wV/Hm2+nqxcqBEubCEhBAB4tCt2xyhCvNlXc2ELCSEAwKOdrHD8d8cf7BFVhHijj9YqULH1+wirqqo+//xzq6ewr3Pnzl24cMHqKezr2LFjpunGUzHieyoqKk6fPm31FPbSIUT5747f1urLL7+8ePGiEEJVxFPdVFu/j3Dr1q3333+/1VPY10svvTR37lyrp7Cvfv36ff3111ZPYVNr1qx57LHHrJ7Cvp599tmCggLXPqa3hhAAAJcghAAAWyOEAABbUxwOd18L+KfMnTv36aefDgoKqsuda2try8vLmzZt2tBT4aouXrwopazjPxZc7ty5c82bN1cUxepB7KimpqaysrJJkyZWD2JTFy9e9PPzCwwMrOP9i4uLb7vttmvfx4NCaJrmZ599ZvUUAADfccstt/j7+1/7Ph4UQgAA3I/XCAEAtkYIAQC2RggBALZGCAEAtuZ9V6g3DGPfvn3bt283DIOzrLnfyZMnFy5cuH///latWt17771t27a1eiIbcTgceXl5O3bsKC8v79Klyz333BMcHGz1UHZUVlb2ySef3H333aGhoVbPYiOlpaW7du26/OnDDz+saZpLHtn7VoQLFy4cMmTIm2+++eSTT1o9ix098MAD69atu+WWW44ePdqlS5c9e/ZYPZGNmKY5Z86c4ODg9u3bv/fee6mpqYbhysuToi4uXLhw3333TZw48ezZs1bPYi8ffvjhP/7xjyP/4cK3PHjf2ydM01RVdcOGDenp6efOnbN6HNupqqq6/FbW0aNHd+rUafr06daOZE8VFRWhoaG7d+/u3Lmz1bPYy4QJE26//fZf/vKXhw4d6tixo9Xj2MiTTz7ZuHHjF1980eWP7H0rQlX1vpl9yfdP6FBVVdW4cWMLh7GzDRs2NGnSJDIy0upB7GXVqlUHDx586KGHrB7EprZt2zZ9+vT333+/qqrKhQ9LVFBPixYt2rRp0/jx460exHZGjx4dHh4+atSoefPmhYSEWD2OjVRUVDz++ONvvfUWJ7ezRERERGRk5IULF3Jzc+Pi4s6fP++qR/a+zTLwBJs2bbr//vvnzJkTHh5u9Sy2M2vWrAsXLsyfP//uu+/euXNnRESE1RPZxTPPPHPfffd17ty5urra6lns6Omnn3Z+YJpmv3793njjjd/85jcueWRWhLhh27Zty8zMnDVr1pAhQ6yexY5CQkIiIyNzcnLatm1bUlJi9Tg2Mnfu3Pfffz8xMbFPnz5CiJEjRy5cuNDqoexIVdXk5OQjR4646gFZEeLG7Ny5c9iwYTNmzBg+fLjVs9hOVVVVQECA87jc2bNnP/vss6ioKKuHspEVK1boui6EqKmp6dOnz0svvdS7d2+rh7KRyspK5xVvqqqqSkpKsrOzXfXI3rdr9MSJEyNHjrx06dLhw4fj4uKio6M//PBDq4eykcTExEOHDnXq1Mn56aBBg15++WVrR7KPkpKSRx99NCEhQVGUVatWZWRkvP3227xe5X7V1dWBgYHsGnWz6Ojorl27Nm3adO3atdHR0YsXL3bVleC8L4TV1dW7d+++/GlgYGDXrl0tnMdu9u7dW1lZefnTZs2adejQwcJ57Gbfvn179+5VFCU2NvbWW2+1ehybcjgcZWVlsbGxAQEBVs9iIydOnCgrK6uoqOjUqVNiYqILH9n7QggAgAuxWQYAYGuEEABga4QQAGBrhBAAYGuEEABga4QQAGBrhBAAYGuEEABga4QQAGBrhBAAYGuEEABga/8f71jkXFjIFDAAAAAASUVORK5CYII=", "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.1" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.1", "language": "julia" } }, "nbformat": 4 }