{ "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 `AtomsBuilder`\n", "to build a bulk silicon lattice,\n", "(see AtomsBase integration for more details." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using AtomsBuilder\n", "using PseudoPotentialData\n", "\n", "pseudopotentials = PseudoFamily(\"dojo.nc.sr.lda.v0_4_1.standard.upf\")\n", "model = model_DFT(bulk(:Si); functionals=LDA(), pseudopotentials)\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 `ScfSaveCheckpoints`, which stores the state\n", "of an SCF at each iterations to allow resuming from a failed\n", "calculation at a later point.\n", "See Saving SCF results on disk and SCF checkpoints for details\n", "how to use checkpointing with DFTK." ], "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. This example is a bit artificial, since the norms\n", "of all density differences is available as `scfres.history_Δρ`\n", "after the SCF has finished and could be directly plotted, but\n", "the following nicely illustrates the use of callbacks in DFTK." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "To enable plotting we first define the empty 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 = 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`. The latter is the function\n", "responsible for printing the usual convergence table. Therefore if we simply did\n", "`callback=plot_callback` the SCF would go silent. The chaining of both callbacks\n", "(`plot_callback` for plotting and `ScfDefaultCallback()` for the convergence table)\n", "makes sure both features are enabled. We run the SCF with the chained callback …" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) α Diag Δtime\n", "--- --------------- --------- --------- ---- ---- ------\n", " 1 -8.457241783716 -0.89 0.80 5.0 104ms\n", " 2 -8.460040422389 -2.55 -1.72 0.80 1.0 188ms\n", " 3 -8.460186658041 -3.83 -2.86 0.80 1.2 16.5ms\n", " 4 -8.460214649769 -4.55 -2.90 0.80 3.2 21.6ms\n", " 5 -8.460214818404 -6.77 -2.98 0.80 1.0 15.5ms\n", " 6 -8.460215125208 -6.51 -4.79 0.80 1.0 15.6ms\n", " 7 -8.460215135187 -8.00 -4.38 0.80 3.5 23.0ms\n", " 8 -8.460215135738 -9.26 -5.64 0.80 1.0 15.9ms\n" ] } ], "cell_type": "code", "source": [ "scfres = self_consistent_field(basis; tol=1e-5, callback);" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "… and show the plot" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Plot{Plots.GRBackend() n=1}", "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nO3dd3wUdf7H8e+UzaZCEgglBCkBCSX00A8IBBBEqqJyomcDy1kO5ZDzLKAIWJA7Dw/xRxHbWQAhGpTQRaQISC9BghQVCAQSSJ3y+2MxIiCk7O7M7r6ef/AYhmTmk2zYd77f+RbJNE0BAECgkq0uAAAAKxGEAICARhACAAIaQQgACGgEIQAgoBGEAICARhACAAIaQQgACGgEIQAgoBGEAICAZqMgPHfu3LPPPmt1FbZQXFxsdQkoM141n2Oapq7rVleBstF13e0rg9ooCLOyst59912rq7CFgoICq0tAmfGq+Rxd14uKiqyuAmVTVFRkGIZ7r2mjIAQAwPsIQgBAQCMIAQABjSAEAAQ0HwvCnGKRdaURCZm5bC8MACgPHwvCD38weqZpl2Th2I36XasZAw0AKA8fC8JRCfLAOlLXz7Vf8i+cefo7fekxc2Ev1dK6AAC+yvfyY0IbRQi9xxfaihvVN3bpaUfMZf3UKk6rywIA+CbfC0IhxIQ2im7qifO1qsFi7U2kIACg/Hysa7SEZgjdNH/IMbeeYpgMAKD8fLJF6HoumDHM8dd1et8l2msdlEeb+mqiA/AbK1asmDFjhtVV+K0+ffrce++9nriy7wXh09/99lzww2Ql2inGbNB3nDZndFEUyeriAASwjRs3FhcXDx8+3OpC/NC6devS09MJQiGEmLHHuGR0zPROilMRc/YbR8+bH/ZQI4MsrQ9AYEtISLjlllusrsIPGYaxcOFCD13cx3oU72ggL79sjOjU9srmQWqLKlK7RdqeMzwyBACUgY8FYbhDRF9pjGj9CGlykvJcK7n7F1rqYTfv0AEA8GM+1jV6dX9uINevJN2yXN+dLca28LGMBwBYwt/SomM1acMAZf4h4/aVep5mdTUAANvztyAUQtQKk9b0V4Nk0SVVO3yOR4YA7C5fE1d/q+LXeo/ywyAUQgQr4p1uyl0N5c6p+voTZCEAWxu5Vh+1Vv+jt6pJ24yUNJLQg/wzCF0eaya//SdlULo2ez/DZwDY18S28vJj5n1rdOOyMHx1h/HyNn1aR6V8V65fv/6JEycqWp8QQgjDMJYtW6brFd3qp7i4uGrVqnl5eVf5mFatWmVkZFTwRqXnz0EohLghTvq6v/rqdmPUWr2YNARgS9eFSytvVFb9bN7/9e+y8NUdxsSt+ld91XYx5VwupLCw0DTd0yumaVqvXr3y8/Ov/aFXpSjKs88+GxR0tUnfRUVFhuG9t2w/D0IhRMPK0vqB6s95IiVNO3mlTX0BwHIlWThy7YUsrHgKXu7jjz9+6qmnZs6cWVhYWHJy+fLl48ePnz59+unTp11nVq1atWPHDtfxsWPHFixYIIT47LPPhBBz5syZOXPmzz//XMo75uXlTZs27emnn05NTS05GRwc7DpYvHhxZmbmvHnznnnmmZUrV1b46ysn/w9CIUQlh1jYS+lYXeq0WNuVzSNDAHbkysKVP5kj1+qvbHd/Cj722GNvvvlmgwYNUlNT+/fv7zr50ksvPfTQQzExMTt37mzVqtWpU6eEEO+8887y5ctdH7B3796XXnqpfHcsKCjo1q3b7t27r7vuuqeeemrixIlCCE3TRo0a5Uril19+efDgwfv27YuMjLzllltWrFjhhq+z7PxqHuFVKJKYnKS0iDZ6pGlvdVEG1QmI3wAAWGvZMXPK9rI9VKsZas47YBqG1CZGevq7sn3unQ3lEQ2u/OZ25MiRefPm/fjjj5UqVbrrrrvi4+PXrFnTunXrF198cdOmTU2bNhVCDBo0aPr06c8+++wVrzBo0CAhxN133x0eHl7Kej799NOIiIiZM2cKIbp27ZqUlPTEE0/I8u8qHDRo0PPPPy+EyM7O/uyzz3r06FHar9Z9AiUIXW6PlxMipcHp+voT5qQk1ugG4Fltqkpjm5dtnMuHB43vTxnOIBHjlB5rJpfpbapx5B/+044dOxISEipVqiSEcDgcbdu23bFjR2RkpNPpdKWgEKJLly4bN24sU7VCiJMnT44ePVoIIcvyO++8c/E/bd++vUOHDhdqa9xYVdXMzMz4+PiLP6Zly5aug7i4uNWrV5f17m4RWEEohGhVRfp2gDpkmXbbCn1OVyU04L4BALwnyilSapUhy17dYSzINFbcqNYIEclf6B8dNGZ2UcoWhn/g/PnzBQW/jZLIz88PDQ0NDQ11jaaRJKnkpBBCUZSS0aFXH94phAgLCxs6dKgQwnWRS25aMijGMIzCwkLX9S+mKOUcEOtGgdhDWDNUrLpRDVNFp8XaoVweGQKwhYtHx1z8vPDyORXls2PHjp07dwohjh49+s0333Tp0qVu3brVqlX79NNPhRD5+fkff/xx9+7dhRB16tTZsmWLEMI0zfnz57s+PSgoKDQ0tGRATYnQ0NBBgwYNGjRo4MCBl9908eLF58+fF0IsXLgwNja2du3a7vli3CoQg1AI4VTE7K7KQ03kjou1lT+ThQAsdvkY0ZIsvGRORbklJibed999AwYMSEpKGjduXMOGDVVVnT179uOPP963b9/ExMTExMQRI0YIIe65556vv/66W7dubdu2vXgawyOPPJKUlNS2bdutW7eW8qaNGjXq3Llz//79H3zwwbfeeuuSB4Q2EdA9gyMT5LoR0q3LtWdaKY+wxz0Ai0zZZkzZri/rp7ap+rvexevCpeX9lOQ0/YFv9JldKtqFWKVKlbS0tO+//z4uLi42NtZ1skePHgcOHNizZ09MTExJc61WrVp79+7dtWtXXFxcTExMSe/o5MmTX3jhhXPnzkVERJTypp06dXr44Yf379/fuHFj1xPKoKCg48ePh4WFCSGWLFlSMpXi7rvvdsWw9wV0EAoheteSvh2gDkzXd2ab/+mkOEhDAF4XX0mk9700BV3qRkgr+ylfHnVPx1VQUFC7du0uORkSEtK6detLToaGhiYlJbmOXQHm4nA4oqKiynTTqKio9u3bX3ymWrVqroOLA9XpdDqdV9pmz/N44xfxlaRvB6gn8kWPNO1ERddMAIAyu7mefMUUdKkbIT3QuJzv1XFxca7RKAkJCbfddls56yuvXr16derUqRyfGBsb63A43F7PHwn0FqFLhEMs6KW8vM1IWqQtSFGu8hMJAD5kw4YNrgPXI0Av39019bAc0tPT3VvJ1dEivEASYmwL+ZV2cr+vtPmZLEsKAIGCFuHvDKsvN4qUBqfrm7LMl9q6Z/oOAMDOaBFeqkW0tHGguuGEOWCpllNsdTUAAA8jCK+garBI76vWryR1SdUymXEPAH6NrtErU2Xx747KzL1Gp8Xae8lqz1g6SQFcg6qqc+fO9fJAjwCRnZ1dsmyp2xGEVzMyQW4cKd22Qn+0qTy2Ba1nAFczatSo5ORkq6vwW3FxcR66MkF4DX+qIa0fqAxO1w/kmNM7K0GkIYA/EBER0aZNG6urQJnxvn5ttcOkVTeq2UUi+QvtF2bcA4B/IQhLJdwhPumpDLhObveZtukkw2cAwH8QhKXlmnH/Zmel/1Lt3QPMuAcAP0EQlk3/66RlfdXnNxtPbXLbJmEAAAsRhGWWGC1tHKR+d9Lsv1Q7U2R1NQCAiiEIy6OKU3x5g9o8Wmq/SNt7hoYhAPgwgrCcVFlMTlKeTJT/9Ln2+WGyEAB8FUFYIfcnyIt6qQ98o0/ZxvAZAPBJBGFFdaoubRigzD9kDF+p52tWVwMAKCOC0A1qhUlr+quqJLp8rh05TzcpAPgSgtA9ghXxTndlWD2502J9wwmyEAB8BkHoNq4Z9zO7KAPTtbn7eWQIAL6BIHSzvrWlNf3VKduNUWt1jTQEANsjCN3v+srSupvUzFyz/1Itu9DqagAAV0UQekSUUyy5QW1ZRWq3SNvNjHsAsDGC0FMUSUxOUsa3kZO/0Bb9SCcpANgUG/N61vB4OaGyNGSZvvWU+VxrRbK6HgDAJWgRelzrqtK6AcqXR83bV+h5zLgHAJshCL0hNlRafaMaoorOqdqP53hkCAA2QhB6iVMRc7oqDzaWOy7WVv9MFgKAXRCEXjUyQZ7TVb11hTZrH8NnAMAWCEJv6xMnrb1JnbrDGLVWLyYNAcBqBKEFGlSS1g9Uj+eLnmnaiXyrqwGAwEYQWiPCIRb2UnrESp1StZ3ZPDIEAMsQhJaRhHi+tTKxrdwzTVtwiE5SALAGQWixW+vLX92gPrHB6LVEm595aRx+e8L8xybdksIAIEAQhNZrWUXaNFDNLhJ3rNIv3r/p2xPmgKVap+q8RgDgQbzJ2kLVYLHuJrVvbfn+tfq/d+ni1xSc01Xtfx3rsgGAB7HWqF0EyWJBivLsZvG39cbGXxxf/UwKAoA3EIT2MqGNUmSIV7aJW+NlUhAAvICuUXv59oQ5a5/xcIL+0UFjwlaGyQCAx9EitJGS54Ldogqujw567Fs9XBWjExWr6wIAf0YQ2sXFo2Nyc8Vfm8g5xebfNxohivRgExruAOApvMPaxa5s893uvxsd848WyqgEefz3xpkiC+sCAD9Hi9Au7mt0hV9KpndWHIo+OF378gbVSRcpAHgALUK7m9peiQmW7lqtG6xICgAeQBDanSyJ95KVU4XmI98yiBQA3I8g9AFBspifoq47br68nbW5AcDNCELfUMkhvuij/HePMWc/WQgA7sRgGZ8RGyql91W6fa7XCJH61mbRGQBwD1qEvqRBJWlhL+Wu1dq3Jxg5AwDuQRD6mHYx0pxu6pB0bd9ZshAA3IAg9D031pZeaa/0+1L/Jd/qUgDA9xGEPumOBvJfrpd7L9FYdAYAKogg9FXPtJKTa0pD0rVCphcCQAUQhD7s9Q5KFRadAYCKIQh9mCyJ95OVrAIWnQGA8iMIfVuQLBb0UtcdN19h0RkAKBeC0Oe5Fp15c48xl0VnAKDsCEJ/EBsqLb1Befo748ujPC0EgLIhCP1Ew8rSwl7KnatYdAYAyoYg9B8sOgMA5UAQ+pUba0svs+gMAJQFQehvRrDoDACUBUHoh1h0BgBKjyD0T65FZ/6yhkVnAOAaCEL/5Fp05mQ+i84AwDUQhH7LtejMN7+Yr+5goj0A/CGC0J9Vcoi0G5Tpu1l0BgD+EEHo52JDpS/6KE9t0ll0BgCuiCD0f00ipUW9VBadAYArIggDQvtq0uyuLDoDAFdAEAaK/tdJL7Rl0RkAuBRBGEDuayTfdb3cZ4l2lkVnAOBXBGFgebaV3K2mNJhFZwDgVwRhwJnWQYlm0RkA+BVBGHBkSbzXXTl23nyURWcAgCAMTMGKSO2trv3FfI1FZwAEPIIwQFUOEmk3KP/ZbbyTQRYCCGgEYeByLTozdqP+FYvOAAhgBGFAcy06M2KVtp5FZwAEKoIw0LkWnRmcru1n0RkAAYkghOh/nTShjdKXRWcABCSCEEIIcX+CfGdDFp0BEIgIQlzwXGsWnQEQiAhC/OZ1Fp0BEHgIQvxGYdEZAIGHIMTvsOgMgEBDEOJSlYPEF32UN3ax6AyAgEAQ4gpqhUlpN7DoDICAQBDiylyLzty1WtucRRYC8GcEIf5Q+2rS239S+n/FojMA/BlBiKu56Tp5Qhul31f6cRadAeCnCEJcw/0J8h0NpBu/0nKLrS4FADyAIMS1Pd9a6VhNGsSiMwD8EUGIUpnWUYlySnez6AwAv0MQolQUSbzbTTl63hy7kVYhAL9CEKK0QlSR2ltdesycyqIzAPwIQYgyqBwk0voo/95lzGPRGQD+giBE2bgWnfk7i84A8BcEIcqsZNGZLSw6A8D3EYQojwuLzixl0RkAPo8gRDnddJ08vjWLzgDweV4KQl3Xd+3atW/fPtOkAeE/ShadOceiMwB8luqFexQWFrZs2bJVq1a5ubl5eXlpaWlOp9ML94UXPN9aOVWgD0zX0vqoTsXqagCg7LzRInQ4HJs3b/7ggw9SU1ODgoLWrFnjhZvCa6Z1VCKDWHQGgK/yRhDKshwaGiqE0HX92LFjcXFxXrgpvEaRxHvdlSPnzac2segMAN/j5q7R6dOnnzx58uIzDzzwQI0aNVzHY8eOTUlJady4sXtvCsuFqOLz3mrXz7UaIcboREZgAfAlbg7Cpk2b5ubmXnzG1RYUQowfP/7kyZNz5sxx7x1hE65FZzqn6lWDxZ0NyUIAPqO0Qfj6668vXbo0IyNjwoQJw4cPLzn/xhtvTJgwoaCg4IYbbpgzZ0737t2v+Omvvfbazp07P/zwQ1nmLdJvuRad6fGFVj1E6hMnWV0OAJRKaWNJkqQ777yzcuXKOTk5JSe3bNny3HPPrVmz5uTJk+fPnx8/fvwVP/f06dMTJ07MyMjo0KFD27ZtU1NT3VA4bKlJpPQZi84A8ClSmSb2JScn33rrrQ888IDrr4888khRUdFbb70lhFi1atUtt9xyyQPCMtmzZ09SUlKfPn1KzgwZMmTw4MHlvqDvOnfuXHh4uNVVlF/aMemRjcpXPbUGlawuxYt8/VULQJqmFRcXh4SEWF0IyiA/Pz8oKEhRSjtbKygoSFWv0fdZoWeEGRkZ/fr1cx03a9YsKysrOzs7KiqqfFdzOp1hYWG33npryZkWLVoEBwdXpEIfVVxc7NNf+JB4cUozh6yRv75Rrh4wbzK+/qoFIE3TFEXhVfMtpmmWKQhL8zyuQkF45syZkl+BIyIihBCnT58udxC6ZlkMGzasIiX5B1mWff1h6qjG4qc8fUC6sfJGNdxhdTVe4QevWqCRf2V1ISgDT7xqFbpWlSpVzp496zo+c+aMECImJsYNRcEvjG+jdKgmDUzXiti7EICNVSgIGzduvH37dtfxtm3bYmNjK1UKpIdCuJYLi86sZtEZAPZV2iA8ePDg5s2bc3NzDx8+vHnzZlf775577lmwYIFr1OiLL754//33e7JU+B7XojNrjxt3rrp00Zkvjpj/3kVTEYD1ShuEb7/99qhRo4QQS5cuHTVq1NatW4UQTZo0mTFjxoMPPtiyZcvmzZuPGzfOg5XCN4WoYlpH9aNM46/rfsvCL4+a967R2scw1xCA9co2fcKjDh06lJycnJmZaXUh1svNzXUNPvIbH/xgjFilj20uv5SkfHnU/MtqbVEvtX01vwpC/3vV/B7TJ3xRWadPlIY3tmEChsfL5zXxwFr9ZKGU+qPuNymYfsz8Uw0p+Pf/JfecMR2yaFDJH75AIBAQhPCS+xvJ60+IWft0hyz6faVFBkmhqghRReUgEaZKoaqIcIgIhwhRRbgqVQoSIYoIc4jKQVKIIkJVERkkQlQRokhRdtrLcs5+Y+oOc2EvtSQL95wxU9L0qR1kghDwFQQhvOTLo+YXh/VpHZQXturvd1ebRYsCXeRrIrtQFOgiXzezC387/unMr/9aZFzyYacKRJEhghUR5RQhinThQBXBiogKki4+vnBw2YdVDZaC3DQH6d3uyl2r9YFLtUW9VSHEvrNm7yX6pCT51vpMTQN8BkEIb7j4uWBCpPTnldqnPdVuNS9uM5Wt/ZSviewisyQjLz525eXB3F/PXPZPWQWi2E1RqkjinW4XsvD5ZtKwtfrEtjKbbwC+hSCEx10yOqZ3Len9ZPXm5ZdnYRmEqCJELX+OntdEviZyis1zxSJfE7nFIrfYzNfFuWKRUyzyNXFeM/eeNfM1kaeJ7CIzXxP5ujhTKPI0M18XZ4tEqCpCFFE5SApziGBZHDwnkpc5k2LEjtPmpG1GVJCoEiyinVK0U0Q7RbRTigiM5XUAX0QQwuM2Z5mLe6vtLpos0buW9F53dVe2We4grKAwVYSpomqwe6J0zxnz4W+M2qHmmUIR7ZRzisxDueJ0oThdaJwudB2YhfqFRHT9GeW8+K+/HVQJliKD3Pu1ltOubLNp1KXfk8xcs3qIFMrbBvwLP9HwuKdbXqGrsE+cVNbssZWSKN131nxygzG1gzygev5ft4Su+tlY1FsNvtLQblc/7c954qe8Cw9Es4vM/WddB0bJmRP5wiGLKOeFflrXQWyYqBny219dBzHBksMzvbDFhrhthX5LffnZVr/dYP9Zs2ea/t/OSv/rfPiFAy5HEALlt++smZJ24blgbu5vzwuvmIWu7tzYUNHmWr8BuCKzJBpd2Xkw18zO+t35rAKhSJdG5sUHsaFSzVAR5RRVnJKzLNOuHLJY3k/tmaaZpniutSyEyDhr9kzTn24pk4LwPwQhUH6j1uovJckjGlxoNimSmNtV+fMq/b97jL81K39jrSQyhRBXbzdfEpklBwdzzZ+Pi5/yDNdfTxUI+aqReXF2uiaoVAu5kIVCiOHxUo80/emW8gONGQcEP0QQAuW3tK96yUwMVRbvd1ckb7WaSh+ZucXidKHpemZ5qsAseXh5+Jzr2DxdKLILxelCU5J+e2Z5Xbg0dYcxdYd4vaN8z/WkIPwTQQiU3xXnI6q2zIsIh4hwSHUu7B96tcjM036LzJ3ZxrrjIqfY3H7aLmsxAm5HEAL4nVBVhKpSXJjIOGuOWGVOSpKPnBP/2a1HBUmu54WAnyEIAVxBxlmz5LlggS7e/cF4J8MQQpCF8D/8TAO4VLEh+i/Vx7e5MDomWBFTkuRwVXx00FhwiF0k4W9oEQK4lEMWX/dXq120PdFt8fL03cbwBvLAOvz2DH/DzzSAK6j2+036JCFeba9M+t4o0P/gEwCfRRACKJUO1aTO1aXXdtA1Cn9DEAIorclJ8rSd+pHzTKWAXyEIAZRW3QhpZIL83GYahfArBCGAMni6pfLVUXNzFo1C+A+CEEAZRDjEM63kJzcwZgb+gyAEUDb3J8inCsWiH+kghZ8gCAGUjSKJ1zsoT2wwiohC+AWCEECZ9YyVGlYS/91NEsIfEIQAymNqB2Xi9/qpQqvrACqMIARQHo0jpZvryS9uZdQMfB5BCKCcxrdR3jtg7D/LVAr4NoIQQDnFBIsnE5WxG3lSCN9GEAIov78lyjuyzWXHaBTChxGEAMovSBYT28pjNuoGUQifRRACqJBb68vhDjEvgw5S+CqCEEBFvdZe+edm47xmdR1AuRCEACqqXYzUrYb0ynamUsAnEYQA3GBKO/k/u4zD53hUCN9DEAJwg7gwaVRj+Vm2KoQPIggBuMe4Fkr6MXPTSRqF8DEEIQD3CHeI51vLj6/XSUL4FoIQgNvc20gu0MXCQ3SQwpcQhADcRpbEK+2UMRuMQgaQwncQhADcqUes1DhSTGerQvgOghCAm73aXpm8Tc8qsLoOoHQIQgBulhApDasvT2CrQvgIghCA+01oo3x00NiVzQBS+ACCEID7RTvFmObKuE08KYQPIAgBeMSjTeW9Z810tiqE7RGEADwiSBYvtZXHbNCZYA+bIwgBeMrN9eRKQWLufjpIYWsEIQAP+ldH5ZnNek6x1XUAf4wgBOBBrapIKbEyWxXCzghCAJ41KUn+727jR7YqhF0RhAA8q1aY9FAT+Z/f8aQQNkUQAvC4vzdXVv5sbmSrQtgSQQjA48IdYkIb+fFvmUkBOyIIAXjDXxrKRYb4NJMOUtgOQQjAG2RJvNpeGbuRrQphOwQhAC/pXlNKjJb+vYtGIeyFIATgPVPbyy9v14/nW10HcBGCEID3xFeShsfLL7BVIeyEIATgVc+3Vj7JNHayVSFsgyAE4FVRTvFUC2XsRhqFsAuCEIC3/bWJfCBHfHWURiFsgSAE4G0OWUxOkkev1zUGkMIGCEIAFhhcV64ZKmazVSFsgCAEYI1X2yvPsVUhbIAgBGCNllWkPnHy5O8ZNQOLEYQALPNSkjxzr3Eol1EzsBJBCMAysaHSX5vK/2CrQliKIARgpbHNlW+Om98cp1EIyxCEAKwUoorxreUnN7BVISxDEAKw2J0NZc0QHx+kgxTWIAgBWMy1VeGYDUaeZnUpCEgEIQDrdasptY1hq0JYgyAEYAuvtpdf26H/wlaF8DqCEIAt1I+QRjSQx29hfj28jSAEYBfPtFIWHjJ2nGYAKbyKIARgF1FO8Y+Wyt/W0yiEVxGEAGzkocbyz3liyREahfAeghCAjaiymJQkP7GBrQrhPQQhAHsZUEeOCxNv7yMJ4SUEIQDbeaW9MmGLfrbI6joQGAhCALbTIlrqV1uetI1RM/AGghCAHb2UpMzaZxzIYdQMPI4gBGBH1UPEo02Vp9mqEJ5HEAKwqScT5Q0nzLW/0CiEZxGEAGwqRBUvtJWfYKtCeBhBCMC+7mggm0J8+AMdpPAgghCAfUlC/Kuj8tRGtiqEBxGEAGytYzWpfTXp9Z00CuEpBCEAu3u5nfz6Dv3nPKvrgJ8iCAHYXb0I6e7r5efYqhCeQRAC8AFPt1IW/WhsyWIAKdyPIATgAyKDxDOtlDEbaRTC/QhCAL7hgQT5eL74gq0K4W4EIQDfoMpiSjvlifV6MQNI4VYEIQCfcWNtqU64eGsvSQh3IggB+JLXOyovbNVPF1pdB/wIQQjAlzSJlAbWkSd9z6gZuA1BCMDHvNhGmZthZJxl1AzcgyAE4GOqhYi/NVPGsVUh3IQgBOB7RifKm7PMNWxVCHcgCAH4nmBFTGwrP/6tbhCFqDCCEIBPuj1eDlXF+2xViAojCAH4JEmIV9sr/9hknGerQlQMQQjAV3WoJnWqLk3dQaMQFUIQAvBhU5Lkf+/Sf8rjUSHKjyAE4MPqRkj3NpKfYSqFEEKIAzlX+IXgiidxMYIQgG97uqWy5KixOeC3KszMNTst1pYd+933If2Y2Wmxdig30L85V0cQAvBtEQ7xTCvlyQ2BvuhavQhpYS/1z6t+y8L0Y+bwldrHPdW6EZK1tdkcQQjA541MkLMKxOIfA72DtHN1aUHKhSx0peAnPdXuNUnBa1CtLgAAKkqRxCvtlce/1fvWlh2B/et95+rSB8nqgK80RRZLblA7VycFry2wf2QA+Isb4qS6EeK/ewK6Ubj3jDl6vT50mVZkivPF4iDDZEqHIATgJ17voLy4VT8VeFsVFhnik0yj1xIt+QvtpzwhSyK9r/p4omLHdE8AABK/SURBVHzP1/riH8nCa6NrFICfaBwpDa0nT9yqT+2gWF2LlxzIMf9vnzFnv9GwkvRYMzncId25SluQonavKXWvqRzIEcNWaKm91V616CC9GlqEAPzHhDbKuweM/f6+VaFuitTDRq8lWvtFWnahWNFPXXuT2raqdNcq7bNev42O+bSn0qiyNDid6RPXQIsQgP+ICRZPJCpPbTIWpPhno/DYefO9A+b03Ua1EDEyQV7cSw759V28XoS0ZbAaF/Zb488hi6/6qm0WFu/MFnUjrCnYJ9AiBOBXRifK20+bl8wr93WGKZYdM4ct11ss0A7mml/0Ub4bpI5M+C0FXS5OQZcaIWJ+inr3Gm33Gb/6hrgXLUIAfiVIFhPbymM26psHqbLvPxr7Kc98N8N8c48REyxGJshzuzlCy/i23aGa9Fp7ZUi6vmGgWjnIM1X6OFqEAPzNrfXlcIeYl+HDUylKmoCJ87WDuebi3heagGVNQZc7G8o9a0l3rmYf4ysjCAH4oVfbKf/c7JNbFf6cJ6ZsM+I/1p7apKfUkg7f7niri9IiuqJt22kdlJwic/yWQF+I7ooIQgB+qH01qWsN6ZXtPvO+X9IEbDa/+GCuuSDlQhMwzE3Prxyy+DRFfe+A+UmmDzeUPYRnhAD808vt5JYLtHsbybUvG0JiK7/ki3f2GzP3GlFOMTJBnt3VEe7wyI2qOMX8FKXPl1rjSKlZlK2/J15GixCAf4oLk0Y1tu9Whaa40ARs9HHxrmzzk54XmoAeSkGXllWkqe2VAUv1rAIP3sXn0CIE4LfGtVAafaJ9l2W2rWqjBlB2ofgk0/jXTsMU4q6G8ltdHFFO7939zw3kLafM4Su1JTeoio2+K1aiRQjAb4U7xHOtZftsVbg5yxy1Vo//uHjZMfNfHZXdN6tjW8jeTEGXl9spiiSe/s4u3xbLEYQA/Nm9jeTThWLhISs7SM8UiZl7jcT52h2r9PoRUsYwx8c9lRTr1v9UJPF+sjo/0/zwB5v2G3sZXaMA/JkiiWkdlPu/1vvVlp1eX3Ztc5Y5c6/xSaaREiu/3kHpWUuySWdktFPMT1F6pmnXV5ba2Knf2BK0CAH4uR6xUkKkeNOLWxWeLRIz9xotFmjDV+r1I6T9t1xoAtoqcJpHSzO7KDcv108G/MAZ7wWhaZqfffbZoUOHvHZHAHB5rb0y6XtvDJV0PQWs91HxsmPma+2VvbeoY1vIVYM9ft/yGVxXvr2+NGSZVhzYXaTeC8LZs2c/8MADq1ev9todAcAlIVIaVl9+YaunhofkFIuZe41WC7XbVuj1I6R9tmwCXtGLbZXIIDFmY0APnPFSEP7000+LFy+++eabvXM7ALjEhDbK/w4abt+E4UIT8H/Fy46Zr7RT9g9Tx7aQY+zaBLycLIkPktX0o+asfYHbKvRSEP7tb3+bNGmSLPNIEoA1op3iyURl3Cb3vN0X6GJehtHmM21Quh7lFNuGqL7SBLxchEMs6KU8/Z2+8WSArsntzlGjW7Zsefzxxy8+k5iYOH369Pfeey8xMbFJkyZuvBcAlNVjzeSZe7X0Y2avCkxd2H3GnJdhzNpntKwiPdVCHlJX9oNp6Y0qSzO7KEOW6RsHKrGhvv/1lFEZgjA7O/vMmTP16tW7+GRubu7mzZtjYmKaNm3asmXLL7/88uJ/VRRFCDFnzpysrKzPPvvsyJEjqamp0dHRN910k1uqB4DSC5LFpCR5zAZ98+AyL6pSoIvUw8bMvcaeM+KOBtLmQep14X4VGAPqyFtPiZuX6StvVL0/z8RapQrCbdu2DR069NChQ5IkFRcXl5zfsmVLv379mjdvnpGRkZycPHv27NDQ0Ms/ffny5a6DRx99tE2bNqQgAKvcXE/+9y5j7n7j3kalfVKz94w5N8OYvc9oUUUamSAPriOrfvqQ59nW8q5s8/H1+n87B1YSlur1rFmz5kcffbRhw4ZLzo8bN+7RRx9dunTpli1blixZsm7duqtfp2XLlnXr1i1foQDgFv/qqDyzWc8pvsaHFerik0yj1xKtZ5ouhNg4UE3vq95Sz29TUAghCTG7q7L2F3Pm3sAaOCOZZmmfju7YsaN169YlLcKzZ89GRUUdPnw4Li5OCDFy5Mjw8PCpU6eWu5Tdu3c3b948IiLiQmWSNHbs2IcffrjcF/Rd586dCw8Pt7oKlA2vmk/IKZZu/1qNChIJlY2nGhcWFxeHhIQU6uLubx1jmmitok0hxIFcad5B+b1MpVmkeXe8flMtw4/D73I/npd6pjve6ax1jrFjHObn5wcFBbmeu5VGcHCwql6j77P8g2WOHTumKEqtWrVcf61Tp862bdvKfTUhRGhoaK1atb7//vuSM+Hh4Q6HJ7cksSvTNHlL9Tm8aj4hXIiB9Yz/7DbWnJBHNlZrhmuKM2T4Mt3pEG1inV8cMWbuNXZlm3c2lDcMlOtF+NVTwFJqGi7eTTbvXiOvH6DE2W8rR0VRyhSEpVH+ICwoKHA4HNKvo4WdTmd+fn4Fq5FlOSoqqoIXAYCrGJ0oCyEmbjPGbDJndxIjlunFhkiIFPU/Kk6MlkYmyIPqyI5AagJerlct6a9N5IHp+tr+akgArEhd/i+xevXq+fn5eXl5rgEyWVlZNWrUcF9hAOApoxPlQl08u0Vv+7mao5mSabaoIn9zkxpfyXYNIKuMbSF/f8p8aJ0+p6v/D5wp/689NWvWrFOnTsmSaatXr+7QoYObqgIAzxqdKNcIEYfPi2dbSYeHOyYnKaTgxSQhZnVVtmaZb+yy45NC9ypVizAvL++NN9745ZdfDMOYMmVKRETEQw89JMvy6NGjH3300ZdeemnDhg1Hjx697bbbPF0uAFRckSFuWa63j5HaROn/2iUNqmPWst/DMMuFqWJRb6XjYq1ZtJRc05+/P6UKQsMwsrOznU7nmDFjsrOzdf3C8qyPPPJIdHT0okWLqlevvnbt2rCwME+WCgBuUGSIm5fpTkW8+ydh6oYzyJGcpq/sp5CFl6sTLr3XXb19hbZugFrff4cOlWH6hKcdOnQoOTk5MzPT6kKsl5ubWzKNBL6CV80nlKTgh8mKMDTX9ImpO4wZew2y8I9M22nM3W98M0ANs8HAmbJOnyiNwB4aBSDAGKboUkP6MFm5eGrg6ET5yUTZLm0C+3m8mdymqnTXat1fv0UEIYAAEqyIvze/wuowIxNkG86Zs4/pnZUj58xXt/vnwBmCEABwDcGK+KyX+u9dRtoRP2wWEoQAgGurGSo+6qHcs0Y7kONvWUgQAgBKpVN16bnWyk1Lr71kuW8hCAEApfVgY7lrDemuVX41cIYgBACUwX86KacLzYlb/WfgDEEIACgDhyw+7qm+vc9IPewnWUgQAgDKpnqIWJii3Pe1vivbH7pICUIAQJm1riq92l4Zskw/U2R1KRVGEAIAymNEA7lPnHTbCs3XR84QhACAcpraXinSxfNbdKsLqRCCEABQTqosPklRPzhgfnzQhwfOEIQAgPKr4hTzU5RHvtV3nPbVHlKCEABQIS2rSK93UAak61kFVpdSLgQhAKCihsfLN9eVbl+paT7YRUoQAgDcYEo7xSGLf3znewNnCEIAgBvIkni/u7rwkDl3v4+1CglCAIB7RDnF4t7K2E36d1m+NHCGIAQAuE3jSOmtLsrNy/QT+VaXUmoEIQDAnQbVke9oIA1drhX5SBcpQQgAcLMJbZTIIPHkBt8YOEMQAgDcTJbEB8nq8mPm/+3zgVYhQQgAcL8Ih1jQS/nnd/qGE3YfOEMQAgA8olFl6Z1u6tDl+rHzts5CghAA4Cl94qRRCfIty/VCGz8uJAgBAB70z1Zy7XBp1Fr7JiFBCADwIEmIWX9StpwyZ+yx6cAZghAA4FnhDrG4lzJhq776Zzs+LCQIAQAeVzdCmtdNHb5SP2K/gTMEIQDAG1JqSY83kwcu1fM0q0v5PYIQAOAlY5rLCZG2GzhDEAIAvGfWn5S9Z8xpO200cIYgBAB4T4gq5qcoU7bpXx21y8NCghAA4FXXhUuf9FT/slr7IccWWUgQAgC8rUsNaVxLZcgy/bwNBs4QhAAACzzaVE6Kke5cpVveKiQIAQDW+E8n5Vie+fI2iwfOEIQAAGsEK2Jhijp9t/HFESubhQQhAMAyNUPFRz2Ve9doGWcty0KCEABgpY7VpMlJyoB0/WyRNQUQhAAAi/3lerlbDemu1bphRbOQIAQAWO+NTsqZIvOFrRYMnCEIAQDWc8jiox7qrH3Gp5nezkKCEABgC9VDxKLeysPr9F3ZXu0hJQgBAHbRqor0antlyDL9jBcHzhCEAAAbGdFA7ltbunW55rUlZwhCAIC9vNpOKTbEs5u9tG0hQQgAsBdVFp+mqB8dND866I2BMwQhAMB2op1ifory13X6liyP95AShAAAO2oRLb3VRRm6XM8q8OyNCEIAgE0NqSsPqyfdtkLTPNlFShACAOxrUpLiVMTYTR4cOEMQAgDsS5bE35sr/7fXmLXvd63CrAJx83I9p9gdt3DDNQAA8JgO1aTWVaWH1+nrj18YOHOyQPRM0+pHiEoON1yfIAQA2JpTEV/eoCZGi55L9J/zpJMFIiVN6xMnvdxOccv1VbdcBQAAz3EqYm1/R/MFxW3TlFqhRt/asrtSUNAiBAD4BKcilvV1FGgiX5cmJbktBQVBCADwCScLRP+l2v3XGwmVzdtX6G6cUEEQAgDsruS54OTWxsfJcr5uDl/ptiwkCAEAtnbJ6BinIj7tqboxCwlCAICtnSk0h8f/bnSMUxGf9FTrVxKFBKG/mjp1qqZpVleBMsjJyZkxY4bVVaBsdu7c+fnnn1tdBa6tYWVpbIsLaZWamrpr1y4hRLAiJicpYe6Y+kAQ2tG0adNyc3OtrgJlcOzYsVmzZlldBcpm48aNaWlpVleBsklLS9u4caN7r0kQAgACGkEIAAhoBCEAIKBJpunxzX9LKSMjo1mzZnFxcVYXYr0ff/yxdu3assyvKT6juLj4+PHj/PT6lnPnzhUUFFStWtXqQlAGWVlZwcHB4eHhpfz44cOHv/DCC1f/GBsFoRDi4MGDVpdgC4WFhU6n0+oqUDa8aj7HMAxd1x0Od+xfAG8pLi5WFKX07YSaNWuGhIRc/WPsFYQAAHgZnW8AgIBGEAIAAhpBCAAIaAQhACCgsUO9jZimuWnTpmXLlp0+fbpZs2bDhw8PCgqyuiiUVkZGxsqVK/v37x8bG2t1Lbg2TdP+97//bd68OTo6eujQoU2aNLG6IlyDaZqpqanffPNNcHDwgAED2rRp464r0yK0kSNHjgwbNuzMmTO1a9eeMWNGSkqKrutWF4VS0TRtxIgRjz322L59+6yuBdemaVq/fv3efPPNuLg4TdO+/fZbqyvCtU2YMGH06NHx8fGhoaHJycnLli1z15VpEdpIbGzsgQMHVFUVQtx9990xMTE7d+5s0aKF1XXh2iZPntynT59Dhw5ZXQhKZfbs2b/88suWLVtc/93gEz799NMXX3zxtttuE0IcOXJk/vz5KSkpbrkyLUIbUVW15L+lpmm6rpd+9QRYaO/evZ988sm4ceOsLgSllZaWNmLEiCVLlrz++us0B31FkyZNtm7dKoQoLi7euXNn06ZN3XVlgtCmHnvssSFDhsTHx1tdCK7BMIz7779/+vTpwcHBVteC0srMzHz77bfnz5+fm5s7dOjQN954w+qKcG1vvfXW8uXL69evHxsbW69evYcffthdV6ZbwI6eeeaZbdu2rVq1yupCcG1Tp05NTEzs0qWL1YWgDGRZbtSo0dy5c4UQSUlJd9xxxyOPPGJ1UbiGxx9/vHbt2nPnzs3JyRk5cuSsWbPuu+8+t1yZJdZsZ+LEiR988MHKlSurVatmdS24ttatW+fn54eFhQkhtm/fXq9evTFjxrjr/yc85KabbmrSpMmUKVOEEIcPH65Tp05OTk5ERITVdeEPFRQUhIaG7t279/rrrxdCzJ49e8aMGe7aoZcWob1MnTp13rx5q1atIgV9xfvvv5+Xl+c67t279+jRo2+88UZrS8I1DRo06N133zVNU5Kk9evX165dmxS0OafTGR4e/sMPP7iC8IcffqhSpYq7Lk6L0Eb27duXkJBQr1696Oho15nXXnutW7du1laF0qtRo8aHH36YnJxsdSG4hvz8/J49ezqdzvj4+MWLF7/99tsDBw60uihcwxtvvPH8888PGTIkOzt75cqVixcv7ty5s1uuTBDaSH5+/u7duy8+Ex8fHxkZaVU9KCtX1yhtC59QXFy8cuXKc+fOdejQgTUQfEVmZuaWLVuCg4M7deoUFRXlrssShACAgMb0CQBAQCMIAQABjSAEAAQ0ghAAENAIQgBAQCMIAQABjSAEAAQ0ghAAENAIQgBAQCMIAQABjSAEAAS0/wdvvZm3IYL7cwAAAABJRU5ErkJggg==", "text/html": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ], "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ] }, "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", "> **Debugging with callbacks**\n", ">\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.11.4" }, "kernelspec": { "name": "julia-1.11", "display_name": "Julia 1.11.4", "language": "julia" } }, "nbformat": 4 }